updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_vtox_transforms / da_transform_vvtovp_adj.inc
blobc615c01c151579c4a92c3da13c5157684e36c08f
1 subroutine da_transform_vvtovp_adj(grid, evec, eval, vertical_wgt, vp, vv, mz, levels)
3    !---------------------------------------------------------------------------
4    ! Purpose: Adjoint of da_transform_vvtovp.
5    !---------------------------------------------------------------------------
7    implicit none
8    
9    type (domain), intent(in) :: grid
10    integer, intent(in) :: mz                         ! # vertical modes.
11    integer, intent(in) :: levels                     ! no. of vertical levels
13    real*8, intent(in)  :: evec(jds:jde,kds:kde,1:mz) ! Eigenvectors.
14    real*8, 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.
19    integer :: i, j, m, k, ij             ! Loop counters.
20    real    :: temp
22    if (trace_use_dull) call da_trace_entry("da_transform_vvtovp_adj")
24    !-------------------------------------------------------------------
25    ! [1.0] Apply inner-product weighting if vertical_ip /= vertical_ip_0:
26    !------------------------------------------------------------------- 
28    if (vertical_ip /= 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
33    !-------------------------------------------------------------------
34    ! [2.0] Perform vv(i,j,m) = L^{1/2} E^T vp(i,j,k) transform:
35    !------------------------------------------------------------------- 
37    !$OMP PARALLEL DO &
38    !$OMP PRIVATE ( ij, m, k, j, i, temp )
39    do ij = 1 , grid%num_tiles
40       vv(:,grid%j_start(ij):grid%j_end(ij),:) = 0.0
41       do m = 1, mz
42          do k = kts, levels
43             do j = grid%j_start(ij), grid%j_end(ij)
44                temp = evec(j,k,m) * eval(j,m)
45    
46                do i = its, ite
47                   vv(i,j,m) = vv(i,j,m) + temp*vp(i,j,k)
48                end do
49             end do
50          end do
51       end do
52    end do
53    !$OMP END PARALLEL DO
55    if (trace_use_dull) call da_trace_exit("da_transform_vvtovp_adj")
57 end subroutine da_transform_vvtovp_adj