Merge remote-tracking branch 'origin/release-v4.6.1'
[WRF.git] / var / da / da_vtox_transforms / da_transform_vtovv_adj.inc
blob3a667ee0b4e722f2e8be2746e6fd7238423b8020
1 subroutine da_transform_vtovv_adj(grid, cv_size, be, cv, vv &
2 #if (WRF_CHEM == 1)
3                         , vchem &
4 #endif
5                         )
7    !-----------------------------------------------------------------------
8    ! Purpose: TBD
9    !-----------------------------------------------------------------------
11    implicit none
13    type(domain),  intent(inout) :: grid
14    integer,       intent(in)    :: cv_size ! Size of cv array.
15    type(be_type), intent(in)    :: be   ! Background error structure.
16    real,          intent(inout) :: cv(cv_size)   ! control variables.
17    type(vp_type), intent(inout) :: vv   ! Grid point/EOF control var.
19 #if (WRF_CHEM == 1)
20    type(xchem_type), optional, intent(inout) :: vchem   ! Grid point/EOF equivalent.
21    integer :: ic
22 #endif
24    integer :: s(4)   ! Index bounds into arrays.
25    integer :: n      ! Loop counter.
26    integer :: mz     ! Vertical truncation.
27    integer :: ne     ! Ensemble size.
29    logical :: scaling
31    if (trace_use) call da_trace_entry("da_transform_vtovv_adj")
33    if( .not. use_rf .or. do_normalize ) s(1:2)=1
35    !-------------------------------------------------------------------------
36    ! [2.0] Perform VToVV Transform:
37    !-------------------------------------------------------------------------
39    ! [2.1] Transform 1st control variable:
40    mz = be % v1 % mz
41    s(3)=s(1)+mz-1
42    if( do_normalize )call da_transform_rescale(mz,be%sd(:,:,s(1):s(3)),vv%v1)
43    if( use_rf .and. mz > 0 .and. len_scaling1(1) /= 0.0) then
44       call da_transform_through_rf_adj(grid, mz, be % v1 % rf_alpha, be % v1 % val, vv % v1)
45    elseif( mz > 0 ) then
46       s(4)=s(2)+nij(0,0,2)*nij(0,1,2)*mz-1
47       call da_transform_through_wavelet_adj(grid,mz,be%wsd(:,:,s(1):s(3)),cv(s(2):s(4)),vv%v1)
48       s(2)=s(4)+1
49    else
50 !      print'(a,": be%v1%mz=",I0)',__FILE__,mz
51    endif
52    s(1)=s(3)+1
54    ! [2.2] Transform 2nd control variable:
56    mz = be % v2 % mz
57    s(3)=s(1)+mz-1
58    if( do_normalize )call da_transform_rescale(mz,be%sd(:,:,s(1):s(3)),vv%v2)
59    if( use_rf .and. mz > 0 .and. len_scaling2(1) /= 0.0) then
60       call da_transform_through_rf_adj(grid, mz, be % v2 % rf_alpha, be % v2 % val, vv % v2)
61    elseif( mz > 0 ) then
62       s(4)=s(2)+nij(0,0,2)*nij(0,1,2)*mz-1
63       call da_transform_through_wavelet_adj(grid,mz,be%wsd(:,:,s(1):s(3)),cv(s(2):s(4)),vv%v2)
64       s(2)=s(4)+1
65    else
66 !      print'(a,": be%v2%mz=",I0)',__FILE__,mz
67    endif
68    s(1)=s(3)+1
70    ! [2.3] Transform 3rd control variable
72    mz = be % v3 % mz
73    s(3)=s(1)+mz-1
74    if( do_normalize )call da_transform_rescale(mz,be%sd(:,:,s(1):s(3)),vv%v3)
75    if( use_rf .and. mz > 0 .and. len_scaling3(1) /= 0.0) then
76       call da_transform_through_rf_adj(grid, mz, be % v3 % rf_alpha, be % v3 % val, vv % v3)
77    elseif( mz > 0 ) then
78       s(4)=s(2)+nij(0,0,2)*nij(0,1,2)*mz-1
79       call da_transform_through_wavelet_adj(grid,mz,be%wsd(:,:,s(1):s(3)),cv(s(2):s(4)),vv%v3)
80       s(2)=s(4)+1
81    else
82 !      print'(a,": be%v3%mz=",I0)',__FILE__,mz
83    endif
84    s(1)=s(3)+1
85    
86    ! [2.4] Transform 4th control variable
87       
88    mz = be % v4 % mz
89    s(3)=s(1)+mz-1
90    if( do_normalize )call da_transform_rescale(mz,be%sd(:,:,s(1):s(3)),vv%v4)
91    if( use_rf .and. mz > 0 .and. len_scaling4(1) /= 0.0) then
92       call da_transform_through_rf_adj(grid, mz, be % v4 % rf_alpha, be % v4 % val, vv % v4)
93    elseif( mz > 0 ) then
94       s(4)=s(2)+nij(0,0,2)*nij(0,1,2)*mz-1
95       call da_transform_through_wavelet_adj(grid,mz,be%wsd(:,:,s(1):s(3)),cv(s(2):s(4)),vv%v4)
96       s(2)=s(4)+1
97    else
98 !      print'(a,": be%v4%mz=",I0)',__FILE__,mz
99    endif
100    s(1)=s(3)+1
102    ! [2.5] Transform 5th control variable
104    mz = be % v5 % mz
105    s(3)=s(1)+mz-1
106    if( do_normalize )call da_transform_rescale(mz,be%sd(:,:,s(1):s(3)),vv%v5)
107    if( use_rf .and. mz > 0 .and. len_scaling5(1) /= 0.0) then
108       call da_transform_through_rf_adj(grid, mz, be % v5 % rf_alpha, be % v5 % val, vv % v5)
109    elseif( mz > 0 ) then
110       s(4)=s(2)+nij(0,0,2)*nij(0,1,2)*mz-1
111       call da_transform_through_wavelet_adj(grid,mz,be%wsd(:,:,s(1):s(3)),cv(s(2):s(4)),vv%v5)
112       s(2)=s(4)+1
113    else
114 !      print'(a,": be%v5%mz=",I0)',__FILE__,mz
115    endif
116    s(1)=s(3)+1
118    if ( .not. use_rf .and. cloud_cv_options > 0 ) then
119       call da_error(__FILE__,__LINE__,(/"no da_transform_through_wavelet_adj for v6-v11"/))
120    end if
122    if ( use_rf .and. cloud_cv_options <= 1 ) then
123       vv % v6  = 0.0
124       vv % v7  = 0.0
125       vv % v8  = 0.0
126       vv % v9  = 0.0
127       vv % v10 = 0.0
128       vv % v11 = 0.0
129    end if
131    if ( use_rf .and. cloud_cv_options >= 2 ) then
132       select case ( cloud_cv_options )
133          case ( 2 )
134 !hcl-check array index of len_scaling
135             mz = be % v6 % mz
136             if ( mz > 0 .and. len_scaling6(1) > 0.0 ) then
137                call da_transform_through_rf_adj(grid, mz, be % v6 % rf_alpha, be % v6 % val, vv % v6)
138             end if
139             mz = be % v7 % mz
140             if ( mz > 0 .and. len_scaling7(1) > 0.0 ) then
141                call da_transform_through_rf_adj(grid, mz, be % v7 % rf_alpha, be % v7 % val, vv % v7)
142             end if
143             mz = be % v8 % mz
144             if ( mz > 0 .and. len_scaling8(1) > 0.0 ) then
145                call da_transform_through_rf_adj(grid, mz, be % v8 % rf_alpha, be % v8 % val, vv % v8)
146             end if
147             mz = be % v9 % mz
148             if ( mz > 0 .and. len_scaling9(1) > 0.0 ) then
149                call da_transform_through_rf_adj(grid, mz, be % v9 % rf_alpha, be % v9 % val, vv % v9)
150             end if
151             mz = be % v10 % mz
152             if ( mz > 0 .and. len_scaling10(1) > 0.0 ) then
153                call da_transform_through_rf_adj(grid, mz, be % v10 % rf_alpha, be % v10 % val, vv % v10)
154             end if
155          case ( 3 )
156             scaling = .true.
157             mz = be % v6 % mz
158             if ( mz > 0 .and. len_scaling6(1) > 0.0 ) then
159                call da_transform_through_rf_adj(grid, mz, be % v6 % rf_alpha, be % v6 % val, vv % v6, scaling)
160             end if
161             mz = be % v7 % mz
162             if ( mz > 0 .and. len_scaling7(1) > 0.0 ) then
163                call da_transform_through_rf_adj(grid, mz, be % v7 % rf_alpha, be % v7 % val, vv % v7, scaling)
164             end if
165             mz = be % v8 % mz
166             if ( mz > 0 .and. len_scaling8(1) > 0.0 ) then
167                call da_transform_through_rf_adj(grid, mz, be % v8 % rf_alpha, be % v8 % val, vv % v8, scaling)
168             end if
169             mz = be % v9 % mz
170             if ( mz > 0 .and. len_scaling9(1) > 0.0 ) then
171                call da_transform_through_rf_adj(grid, mz, be % v9 % rf_alpha, be % v9 % val, vv % v9, scaling)
172             end if
173             mz = be % v10 % mz
174             if ( mz > 0 .and. len_scaling10(1) > 0.0 ) then
175                call da_transform_through_rf_adj(grid, mz, be % v10 % rf_alpha, be % v10 % val, vv % v10, scaling)
176             end if
177       end select
178    end if
180    ! [2.7] Transform w control variables
182    if ( use_rf ) then
183       if ( .not. use_cv_w ) then
184          vv % v11 = 0.0
185       else
186          mz = be % v11 % mz
187          if ( mz > 0 .and. len_scaling11(1) > 0.0 ) then
188             if ( cloud_cv_options == 2 ) then
189                call da_transform_through_rf_adj(grid, mz, be % v11 % rf_alpha, be % v11 % val, vv % v11)
190             else if ( cloud_cv_options == 3 ) then
191                scaling = .true.
192                call da_transform_through_rf_adj(grid, mz, be % v11 % rf_alpha, be % v11 % val, vv % v11, scaling)
193             end if
194          end if
195       end if
196    end if
198    ! [2.8] Transform alpha control variable
200    ne = be % ne
201    if (ne > 0) then
202       mz = be % alpha % mz
203       if( do_normalize )then
204          do n = 1, ne
205             call da_transform_rescale(mz,be%alpha%sd,vv%alpha(:,:,:,n))
206          end do
207       endif
208       if( use_rf )then
209          do n = 1, ne
210             if ( anal_type_hybrid_dual_res ) then
211                call da_transform_through_rf_adj_dual_res(grid % intermediate_grid, mz, be % alpha % rf_alpha, &
212                                                          be % alpha % val, vv % alpha(:,:,:,n))
213             else
214                call da_transform_through_rf_adj(grid, mz, be % alpha % rf_alpha, be % alpha % val, vv % alpha(:,:,:,n))
215             endif
216          end do
217       else
218          do n = 1, ne
219             s(4)=s(2)+nij(0,0,2)*nij(0,1,2)*mz-1
220             call da_transform_through_wavelet_adj(grid,mz,be%alpha%wsd,cv(s(2):s(4)),vv%alpha(:,:,:,n))
221             s(2)=s(4)+1
222          end do
223       endif
224    endif
226 #if (WRF_CHEM == 1)
227 !  if (present(vchem) .and. iv%info(chemic_surf)%nlocal > 0) then
228   if (present(vchem)) then
229 !!!   do ic = 1, num_chem-1
230    do ic = PARAM_FIRST_SCALAR, num_chem
231       mz = be % v12 (ic-1) % mz
232       if( use_rf .and. mz > 0 ) then
233          call da_transform_through_rf_adj(grid, mz, be % v12 (ic-1) % rf_alpha, &
234                                                 be % v12 (ic-1) % val, &
235                                                 vchem % chem_ic (:,:,:,ic) )
236       elseif( .not. use_rf ) then
237          call da_error(__FILE__,__LINE__,(/"no da_transform_through_wavelet for chem_ic"/))
238       endif
239    end do
240   end if
241 #endif
243    if( use_rf )then
244       !-------------------------------------------------------------------------
245       ! [1.0] Fill 1D cv array from 3-dimensional vv arrays.
246       !-------------------------------------------------------------------------
247 #if (WRF_CHEM == 1)
248       if (present(vchem)) then
249          call da_vv_to_cv( vv, grid%xp, be%cv_mz, be%ncv_mz, cv_size, cv, be%cv_mz_chemic, vchem)
250       else
251 #endif
252          call da_vv_to_cv( vv, grid%xp, be%cv_mz, be%ncv_mz, cv_size, cv)
253 #if (WRF_CHEM == 1)
254       end if
255 #endif
256    endif
258    if (trace_use) call da_trace_exit("da_transform_vtovv_adj")
260 endsubroutine da_transform_vtovv_adj