updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_vtox_transforms / da_transform_vtox_inv.inc
blob56c56c2433c6e1b3bb3fad79a7bbfedf1ef2694e
1 subroutine da_transform_vtox_inv(grid, cv_size, xbx, be, ep, cv, vv, vp)
3    !--------------------------------------------------------------------------
4    ! Purpose: Inverse control variable transform v = U^{-1} x'. 
5    !--------------------------------------------------------------------------
7    implicit none
9    type(domain),   intent(inout) :: grid
10    integer,        intent(in)    :: cv_size ! Size of cv array.
11    type(xbx_type), intent(in)    :: xbx  ! For header & non-grid arrays.
12    type(be_type),  intent(in)    :: be   ! background errors.
13    type(ep_type),  intent(in)    :: ep   ! Ensemble perturbations.
14    real,           intent(out)   :: cv(1:cv_size)   ! control variables.
15    type(vp_type),  intent(out)   :: vv   ! grdipt/eof cv (local).
16    type(vp_type),  intent(inout) :: vp   ! grdipt/level cv (local).
18    if (trace_use) call da_trace_entry("da_transform_vtox_inv")
20    call da_zero_x (grid%xa)
21    
22    if (.not. use_background_errors) then
23       if (trace_use) call da_trace_exit("da_transform_vtox_inv")
24       return
25    end if
27    !----------------------------------------------------------------------
28    ! [1.0]: Perform inverse of balance tranform: vp = u_p^{-1} dx
29    !----------------------------------------------------------------------
31    if ( cv_options /= 7 ) call da_transform_vptox_inv(grid, vp, be, ep)
33    !----------------------------------------------------------------------
34    ! [2.0]: Perform inverse of vertical transform: vv = L^{-1/2} E^T vp
35    !----------------------------------------------------------------------
37    !if (  cv_options == 3 ) then
38    !
39    !   call da_apply_be( be, cv, vp, grid)
40    !   call da_transform_bal( vp, be, grid)
41    !
42    !else
44    if (vert_corr == vert_corr_2) then      
45       call da_vertical_transform(grid, 'u_inv', be, grid%xb % vertical_inner_product, vv, vp)
46       !call da_write_vp(grid,vv,'vv_afterUvTransf')
47    else
48       vv % v1(its:ite,jts:jte,kts:kte) = vp % v1(its:ite,jts:jte,kts:kte)
49       vv % v2(its:ite,jts:jte,kts:kte) = vp % v2(its:ite,jts:jte,kts:kte)
50       vv % v3(its:ite,jts:jte,kts:kte) = vp % v3(its:ite,jts:jte,kts:kte)
51       vv % v4(its:ite,jts:jte,kts:kte) = vp % v4(its:ite,jts:jte,kts:kte)
52       vv % v5(its:ite,jts:jte,kts:kte) = vp % v5(its:ite,jts:jte,kts:kte)
53       if ( cloud_cv_options >= 2 ) then 
54         vv % v6(its:ite,jts:jte,kts:kte) = vp % v6(its:ite,jts:jte,kts:kte)
55         vv % v7(its:ite,jts:jte,kts:kte) = vp % v7(its:ite,jts:jte,kts:kte)
56         vv % v8(its:ite,jts:jte,kts:kte) = vp % v8(its:ite,jts:jte,kts:kte)
57         vv % v9(its:ite,jts:jte,kts:kte) = vp % v9(its:ite,jts:jte,kts:kte)
58         vv % v10(its:ite,jts:jte,kts:kte) = vp % v10(its:ite,jts:jte,kts:kte)
59       end if
60       if ( use_cv_w ) vv % v11(its:ite,jts:jte,kts:kte) = vp % v11(its:ite,jts:jte,kts:kte)
61       if (be % ne > 0) then
62 !        vv % alpha(its:ite,jts:jte,kts:kte,1:be%ne) = vp%alpha(its:ite,jts:jte,kts:kte,1:be%ne)
63          vv % alpha(its_int:ite_int,jts_int:jte_int,kts_int:kte_int,1:be%ne) =  &
64              vp%alpha(its_int:ite_int,jts_int:jte_int,kts_int:kte_int,1:be%ne)
65       end if
66    end if
68    !----------------------------------------------------------------------
69    ! [3.0]: Perform inverse of recursive filter: cv = u_h^{-1} vv
70    !----------------------------------------------------------------------
72    !if (global) then
73    !   call da_transform_vtovv_global(cv_size, xbx, be, cv, vv)
74    !else if ( (fg_format == fg_format_wrf_arw_regional .or.   &
75    !           fg_format == fg_format_wrf_nmm_regional) .and. &
76    !           (.not. cv_options == 3) )then
78        call da_transform_vtovv_inv(grid, cv_size, be, cv, vv)
80    !end if
82    !end if
83    
84    if (trace_use) call da_trace_exit("da_transform_vtox_inv")
86 end subroutine da_transform_vtox_inv