da_transform_vvtovp.inc

References to this file elsewhere.
1 subroutine da_transform_vvtovp(evec, eval, vertical_wgt, &
2                                 vv, vp, mz, levels)
3 
4    !---------------------------------------------------------------------------
5    ! Purpose: Transform from fields on vertical EOFS to fields on vertical 
6    ! levels.
7    !
8    ! Method:  Perform vp(i,j,k) = P E L^{1/2} vv(i,j,m) transform.
9    !---------------------------------------------------------------------------
10 
11    implicit none
12    
13    integer, intent(in)   :: mz                         ! # vertical modes.
14    integer, intent(in)   :: levels                     ! # no. of levels  
15 
16    real, intent(in)      :: evec(jds:jde,kds:kde,1:mz) ! Eigenvectors.
17    real, intent(in)      :: eval(jds:jde,1:mz)         ! Eigenvalues.
18    real, intent(in)      :: vertical_wgt(ims:ime,jms:jme,kms:kme) ! Weighting.
19    real, intent(in)      :: vv(ims:ime,jms:jme,kms:kme)   ! CV in EOF space.
20    real, intent(out)     :: vp(ims:ime,jms:jme,kms:kme)! CV in level space.
21    
22    integer               :: i, j, k, m                 ! Loop counters.
23    real                  :: temp
24 
25    if (trace_use) call da_trace_entry("da_transform_vvtovp")
26 
27    !-------------------------------------------------------------------
28    ! [1.0] Perform vp(i,j,k) = E L^{1/2} vv(i,j,m) transform:
29    !------------------------------------------------------------------- 
30 
31    vp = 0.0
32    do k = kts, levels
33       do m = 1, mz
34          do j = jts, jte
35             temp = evec(j,k,m) * eval(j,m)
36 
37             do i = its, ite
38                vp(i,j,k) = vp(i,j,k) + temp*vv(i,j,m)
39             end do
40          end do
41       end do
42    end do
43    
44    !-------------------------------------------------------------------
45    ! [2.0] Apply inner-product weighting if vertical_ip > 0:
46    !------------------------------------------------------------------- 
47 
48    if (vertical_ip /= 0) then
49       vp(its:ite,jts:jte,kts:levels) = vp(its:ite,jts:jte,kts:levels) / &
50                                     vertical_wgt(its:ite,jts:jte,kts:levels)                          
51    end if
52 
53    if (trace_use) call da_trace_exit("da_transform_vvtovp")
54 
55 end subroutine da_transform_vvtovp
56 
57