updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_vtox_transforms / da_transform_vtox.inc
blob77bb57f4b2424afca8dd3acac88b99e9db2a152c
1 subroutine da_transform_vtox(grid, cv_size, xbx, be, ep, cv, vv, vp &
2 #if (WRF_CHEM == 1)
3                         , vchem &
4 #endif
5                         )
7    !--------------------------------------------------------------------------
8    ! Purpose: Control variable transform x' = Uv. 
9    !--------------------------------------------------------------------------
11    implicit none
13    type(domain),   intent(inout) :: grid
14    integer,        intent(in)    :: cv_size ! Size of cv array.
15    type(xbx_type), intent(in)    :: xbx  ! For header & non-grid arrays.
16    type(be_type),  intent(in)    :: be   ! background errors.
17    type(ep_type),  intent(in)    :: ep   ! Ensemble perturbations.
18    real,           intent(in)    :: cv(1:cv_size)   ! control variables.
19    type(vp_type),  intent(out)   :: vv   ! grdipt/eof cv (local).
20    type(vp_type),  intent(inout) :: vp   ! grdipt/level cv (local).
22 #if (WRF_CHEM == 1)
23    type(xchem_type), optional, intent(out) :: vchem   ! Grid point/EOF equivalent.
24 #endif
26    if (trace_use) call da_trace_entry("da_transform_vtox")
28    call da_zero_x (grid%xa)
29    
30 #if (WRF_CHEM == 1)
31    call da_zero_xchem_type (grid%xachem)
32 #endif
34    if (.not. use_background_errors) then
35       if (trace_use) call da_trace_exit("da_transform_vtox")
36       return
37    end if
39    !----------------------------------------------------------------------
40    ! [1.0]: Perform vv = u_h cv transform:
41    !----------------------------------------------------------------------
43 if ( cv_size > 0 ) then
44    if (global) then
45       call da_transform_vtovv_global(cv_size, xbx, be, cv, vv)
46    else if ( (fg_format == fg_format_wrf_arw_regional .or.   &
47               fg_format == fg_format_wrf_nmm_regional) .and. &
48               (.not. cv_options == 3) )then
49 #if (WRF_CHEM == 1)
50       if (present(vchem)) then
51          call da_transform_vtovv(grid, cv_size, be, cv, vv, vchem=vchem)
52       else
53 #endif
54          call da_transform_vtovv(grid, cv_size, be, cv, vv)
55 #if (WRF_CHEM == 1)
56       end if
57 #endif
58    end if
59 end if
61    !----------------------------------------------------------------------
62    ! [2.0]: Perform vp = u_v vv transform:
63    !----------------------------------------------------------------------
65 #if (WRF_CHEM == 1)
66    if (chem_cv_options >= 10 .and. present(vchem)) then
67       call da_transform_vchemtox(grid, vchem, be)
68    end if
69 #endif
71    if (  cv_options == 3 ) then
73       call da_apply_be( be, cv, vp, grid)
74       call da_transform_bal( vp, be, grid)
75    
76    else
78    if (vert_corr == vert_corr_2) then      
79       call da_vertical_transform(grid, 'u', be, grid%xb % vertical_inner_product, vv, vp)
80    else
81       vp % v1(its:ite,jts:jte,kts:kte) = vv % v1(its:ite,jts:jte,kts:kte)
82       vp % v2(its:ite,jts:jte,kts:kte) = vv % v2(its:ite,jts:jte,kts:kte)
83       vp % v3(its:ite,jts:jte,kts:kte) = vv % v3(its:ite,jts:jte,kts:kte)
84       vp % v4(its:ite,jts:jte,kts:kte) = vv % v4(its:ite,jts:jte,kts:kte)
85       vp % v5(its:ite,jts:jte,kts:kte) = vv % v5(its:ite,jts:jte,kts:kte)
87       if ( cloud_cv_options >= 2 ) then
88          vp % v6(its:ite,jts:jte,kts:kte)  = vv % v6(its:ite,jts:jte,kts:kte)
89          vp % v7(its:ite,jts:jte,kts:kte)  = vv % v7(its:ite,jts:jte,kts:kte)
90          vp % v8(its:ite,jts:jte,kts:kte)  = vv % v8(its:ite,jts:jte,kts:kte)
91          vp % v9(its:ite,jts:jte,kts:kte)  = vv % v9(its:ite,jts:jte,kts:kte)
92          vp % v10(its:ite,jts:jte,kts:kte) = vv % v10(its:ite,jts:jte,kts:kte)
93       end if
94       if ( use_cv_w ) then
95          vp % v11(its:ite,jts:jte,kts:kte) = vv % v11(its:ite,jts:jte,kts:kte)
96       end if
98       if (be % ne > 0) then
99 !        vp % alpha(its:ite,jts:jte,kts:kte,1:be%ne) = vv%alpha(its:ite,jts:jte,kts:kte,1:be%ne)
100          vp % alpha(its_int:ite_int,jts_int:jte_int,kts_int:kte_int,1:be%ne) =  &
101              vv%alpha(its_int:ite_int,jts_int:jte_int,kts_int:kte_int,1:be%ne)
102       end if
103    end if
105    !----------------------------------------------------------------------  
106    ! [3.0]: Perform x = u_p vp transform::
107    !----------------------------------------------------------------------
109    call da_transform_vptox(grid, vp, be, ep)
111    end if
112    
113    if (trace_use) call da_trace_exit("da_transform_vtox")
115 end subroutine da_transform_vtox