updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_vtox_transforms / da_transform_vvtovp_inv.inc
blobfa620d7f615e056a7aa3d2f490d82a1f85bb2ec3
1 subroutine da_transform_vvtovp_inv(grid, evec, eval, vertical_wgt, vp, vv, mz, levels)
3    !---------------------------------------------------------------------------
4    ! Purpose: Inverse of da_transform_vvtovp.
5    !          based on da_transform_vvtovp_adj
6    !
7    ! Author: Zhiquan (Jake) Liu, NCAR/MMM, 2015-09
8    !---------------------------------------------------------------------------
10    implicit none
11    
12    type (domain), intent(in) :: grid
13    integer, intent(in) :: mz                         ! # vertical modes.
14    integer, intent(in) :: levels                     ! no. of vertical levels
16    real*8, intent(in)  :: evec(jds:jde,kds:kde,1:mz) ! Eigenvectors.
17    real*8, intent(in)  :: eval(jds:jde,1:mz)         ! Eigenvalues.
18    real, intent(in)    :: vertical_wgt(ims:ime,jms:jme,kms:kme) ! Weighting.
19    real, intent(inout) :: vp(ims:ime,jms:jme,kms:kme)! CV in level space.
20    real, intent(out)   :: vv(ims:ime,jms:jme,kms:kme)! CV in EOF space.
22    integer :: i, j, m, k, ij             ! Loop counters.
23    real    :: temp
25    if (trace_use_dull) call da_trace_entry("da_transform_vvtovp_inv")
27    !-------------------------------------------------------------------
28    ! [1.0] Apply inner-product weighting if vertical_ip /= vertical_ip_0:
29    !------------------------------------------------------------------- 
31    if (vertical_ip /= vertical_ip_0) then
32       vp(its:ite,jts:jte,kts:levels) = vp(its:ite,jts:jte,kts:levels) * &
33          vertical_wgt(its:ite,jts:jte,kts:levels)
34    end if
36    !-------------------------------------------------------------------
37    ! [2.0] Perform vv(i,j,m) = L^{-1/2} E^T vp(i,j,k) transform:
38    !------------------------------------------------------------------- 
40    !$OMP PARALLEL DO &
41    !$OMP PRIVATE ( ij, m, k, j, i, temp )
42    do ij = 1 , grid%num_tiles
43       vv(:,grid%j_start(ij):grid%j_end(ij),:) = 0.0
44       do m = 1, mz
45          do k = kts, levels
46             do j = grid%j_start(ij), grid%j_end(ij)
47                temp = evec(j,k,m) / eval(j,m)
48    
49                do i = its, ite
50                   vv(i,j,m) = vv(i,j,m) + temp * vp(i,j,k)
51                end do
52             end do
53          end do
54       end do
55    end do
56    !$OMP END PARALLEL DO
58    if (trace_use_dull) call da_trace_exit("da_transform_vvtovp_inv")
60 end subroutine da_transform_vvtovp_inv