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