da_transform_vptovv.inc
References to this file elsewhere.
1 subroutine da_transform_vptovv (evec, eval, vertical_wgt, vp, vv, mz, &
2 kds,kde, ims,ime, jms,jme, kms,kme, its,ite, jts,jte, kts,kte)
3
4 !---------------------------------------------------------------------------
5 ! Purpose: Transform from fields on vertical levels to fields on vertical
6 ! EOFS.
7 !
8 ! Method: Perform vv(i,j,m) = L^{-1/2} E^{T} vp(i,j,k) transform.
9 !---------------------------------------------------------------------------
10
11 implicit none
12
13 integer, intent(in) :: mz ! # vertical modes.
14 integer, intent(in) :: kds,kde ! domain dims.
15 integer, intent(in) :: ims,ime, jms,jme, kms,kme ! memory dims.
16 integer, intent(in) :: its,ite, jts,jte, kts,kte ! tile dims
17 real, intent(in) :: evec(jts:jte,kds:kde,1:mz) ! Eigenvectors.
18 real, intent(in) :: eval(jts:jte,1:mz) ! Eigenvalues.
19 real, intent(in) :: vertical_wgt(ims:ime,jms:jme,kms:kme) ! Weighting.
20 real, intent(inout) :: vp(ims:ime,jms:jme,kms:kme)! CV in level space.
21 real, intent(out) :: vv(ims:ime,jms:jme,1:mz) ! CV in EOF space.
22
23 integer :: i, j, m ! Loop counters.
24 real :: ETVp ! E(k,m)^{T}*vp(i,j,k)
25
26 ! Do not add trace, as routine used by gen_be
27 ! if (trace_use) call da_trace_entry("da_transform_vptovv")
28
29 !-------------------------------------------------------------------
30 ! [1.0] Apply inner-product weighting if vertical_ip /= vertical_ip_0
31 !-------------------------------------------------------------------
32
33 if (vertical_ip /= vertical_ip_0) then
34 vp(its:ite,jts:jte,kts:kte) = vp(its:ite,jts:jte,kts:kte) * &
35 vertical_wgt(its:ite,jts:jte,kts:kte)
36 end if
37
38 !-------------------------------------------------------------------
39 ! [2.0] Perform vv(i,j,m) = L^{-1/2} E^T vp(i,j,k) transform:
40 !-------------------------------------------------------------------
41
42 do m = 1, mz
43 do j = jts, jte
44 do i = its, ite
45 ETVp = sum(evec(j,kts:kte,m) * vp(i,j,kts:kte))
46 vv(i,j,m) = ETVp / eval(j,m)
47 end do
48 end do
49 end do
50
51 end subroutine da_transform_vptovv
52
53