1 subroutine da_transform_vvtovp(grid, evec, eval, vertical_wgt, vv, vp, mz, levels)
3 !---------------------------------------------------------------------------
4 ! Purpose: Transform from fields on vertical EOFS to fields on vertical
7 ! Method: Perform vp(i,j,k) = E L^{1/2} vv(i,j,m) transform.
9 ! Zhiquan (Jake) liu's note: 2015-09
10 !-------------------------------------------------------------------------
11 ! 1. evec/eval assumed to vary in y direction (jds:jde) though it may not
12 ! be true in BE file (e.g., likely domain-averaged BE with bin_type=5).
13 ! 2. evec/eval truncated to number of EOF mode mz<=levels
14 ! 3. eval here is in fact square root of eigen values (see da_allocate_background_errors)
15 ! 4. by default, vertical weight not calculated/used
16 !---------------------------------------------------------------------------
20 type (domain), intent(in) :: grid
21 integer, intent(in) :: mz ! # vertical modes.
22 integer, intent(in) :: levels ! # no. of levels
24 real*8, intent(in) :: evec(jds:jde,kds:kde,1:mz) ! Eigenvectors.
25 real*8, intent(in) :: eval(jds:jde,1:mz) ! Eigenvalues.
26 real, intent(in) :: vertical_wgt(ims:ime,jms:jme,kms:kme) ! Weighting.
27 real, intent(in) :: vv(ims:ime,jms:jme,kms:kme) ! CV in EOF space.
28 real, intent(out) :: vp(ims:ime,jms:jme,kms:kme)! CV in level space.
30 integer :: i, j, k, m, ij ! Loop counters.
33 if (trace_use_dull) call da_trace_entry("da_transform_vvtovp")
35 !-------------------------------------------------------------------
36 ! [1.0] Perform vp(i,j,k) = E L^{1/2} vv(i,j,m) transform:
37 !-------------------------------------------------------------------
40 !$OMP PRIVATE ( ij, k, m, j, i, temp )
41 do ij = 1 , grid%num_tiles
42 vp(:,grid%j_start(ij):grid%j_end(ij),:) = 0.0
45 do j = grid%j_start(ij), grid%j_end(ij)
46 temp = evec(j,k,m) * eval(j,m)
49 vp(i,j,k) = vp(i,j,k) + temp*vv(i,j,m)
57 !-------------------------------------------------------------------
58 ! [2.0] Apply inner-product weighting if vertical_ip /= vertical_ip_0:
59 !-------------------------------------------------------------------
61 if (vertical_ip /= vertical_ip_0) then
62 vp(its:ite,jts:jte,kts:levels) = vp(its:ite,jts:jte,kts:levels) / &
63 vertical_wgt(its:ite,jts:jte,kts:levels)
66 if (trace_use_dull) call da_trace_exit("da_transform_vvtovp")
68 end subroutine da_transform_vvtovp