1 subroutine da_transform_vtovv(grid, cv_size, be, cv, vv &
7 !-----------------------------------------------------------------------
9 !-----------------------------------------------------------------------
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(in) :: cv(cv_size) ! control variables.
17 type(vp_type), intent(inout) :: vv ! Grid point/EOF equivalent.
20 type(xchem_type), optional, intent(inout) :: vchem ! Grid point/EOF equivalent.
24 integer :: s(4) ! Index bounds into arrays.
25 integer :: n ! Loop counter.
26 integer :: mz ! Vertical truncations.
27 integer :: ne ! Ensemble size.
33 if (trace_use) call da_trace_entry("da_transform_vtovv")
36 !-------------------------------------------------------------------------
37 ! [1.0] Fill vv arrays from 1-dimensional cv array.
38 !-------------------------------------------------------------------------
40 if (present(vchem)) then
41 call da_cv_to_vv(cv_size, cv, be%cv_mz, be%ncv_mz, vv, be%cv_mz_chemic, vchem=vchem)
44 call da_cv_to_vv(cv_size, cv, be%cv_mz, be%ncv_mz, vv)
49 if( .not. use_rf .or. do_normalize ) s(1:2)=1
51 !-------------------------------------------------------------------------
52 ! [2.0] Perform VToVV Transform:
53 !-------------------------------------------------------------------------
55 ! [2.1] Transform 1st control variable:
59 if( use_rf .and. mz > 0 .and. len_scaling1(1) /= 0.0) then
60 call da_transform_through_rf(grid, mz, be % v1 % rf_alpha, be % v1 % val, vv % v1)
62 s(4)=s(2)+nij(0,0,2)*nij(0,1,2)*mz-1
63 call da_transform_through_wavelet(grid,mz,be%wsd(:,:,s(1):s(3)),cv(s(2):s(4)),vv%v1)
66 ! print'(a,": be%v1%mz=",I0)',__FILE__,mz
68 if( do_normalize )call da_transform_rescale(mz,be%sd(:,:,s(1):s(3)),vv%v1)
71 ! [2.2] Transform 2nd control variable:
75 if( use_rf .and. mz > 0 .and. len_scaling2(1) /= 0.0) then
76 call da_transform_through_rf(grid, mz, be % v2 % rf_alpha, be % v2 % val, vv % v2)
78 s(4)=s(2)+nij(0,0,2)*nij(0,1,2)*mz-1
79 call da_transform_through_wavelet(grid,mz,be%wsd(:,:,s(1):s(3)),cv(s(2):s(4)),vv%v2)
82 ! print'(a,": be%v2%mz=",I0)',__FILE__,mz
84 if( do_normalize )call da_transform_rescale(mz,be%sd(:,:,s(1):s(3)),vv%v2)
87 ! [2.3] Transform 3rd control variable
91 if( use_rf .and. mz > 0 .and. len_scaling3(1) /= 0.0) then
92 call da_transform_through_rf(grid, mz, be % v3 % rf_alpha,be % v3 % val, vv % v3)
94 s(4)=s(2)+nij(0,0,2)*nij(0,1,2)*mz-1
95 call da_transform_through_wavelet(grid,mz,be%wsd(:,:,s(1):s(3)),cv(s(2):s(4)),vv%v3)
98 ! print'(a,": be%v3%mz=",I0)',__FILE__,mz
100 if( do_normalize )call da_transform_rescale(mz,be%sd(:,:,s(1):s(3)),vv%v3)
103 ! [2.4] Transform 4th control variable
107 if( use_rf .and. mz > 0 .and. len_scaling4(1) /= 0.0) then
108 call da_transform_through_rf(grid, mz, be % v4 % rf_alpha, be % v4 % val, vv % v4)
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(grid,mz,be%wsd(:,:,s(1):s(3)),cv(s(2):s(4)),vv%v4)
114 ! print'(a,": be%v4%mz=",I0)',__FILE__,mz
116 if( do_normalize )call da_transform_rescale(mz,be%sd(:,:,s(1):s(3)),vv%v4)
119 ! [2.5] Transform 5th control variable
123 if( use_rf .and. mz > 0 .and. len_scaling5(1) /= 0.0) then
124 call da_transform_through_rf(grid, mz, be % v5 % rf_alpha, be % v5 % val, vv % v5)
125 elseif( mz > 0 ) then
126 s(4)=s(2)+nij(0,0,2)*nij(0,1,2)*mz-1
127 call da_transform_through_wavelet(grid,mz,be%wsd(:,:,s(1):s(3)),cv(s(2):s(4)),vv%v5)
130 ! print'(a,": be%v5%mz=",I0)',__FILE__,mz
132 if( do_normalize )call da_transform_rescale(mz,be%sd(:,:,s(1):s(3)),vv%v5)
135 if ( .not. use_rf .and. cloud_cv_options > 0 ) then
136 call da_error(__FILE__,__LINE__,(/"no da_transform_through_wavelet for v6-v11"/))
139 if ( use_rf .and. cloud_cv_options <= 1 ) then
148 ! [2.6] Transform 6th-10th cloud control variables
150 if ( use_rf .and. cloud_cv_options >= 2 ) then
151 select case ( cloud_cv_options )
153 !hcl-check array index of len_scaling
155 if ( mz > 0 .and. len_scaling6(1) > 0.0 ) then
156 call da_transform_through_rf(grid, mz, be % v6 % rf_alpha, be % v6 % val, vv % v6)
159 if ( mz > 0 .and. len_scaling7(1) > 0.0 ) then
160 call da_transform_through_rf(grid, mz, be % v7 % rf_alpha, be % v7 % val, vv % v7)
163 if ( mz > 0 .and. len_scaling8(1) > 0.0 ) then
164 call da_transform_through_rf(grid, mz, be % v8 % rf_alpha, be % v8 % val, vv % v8)
167 if ( mz > 0 .and. len_scaling9(1) > 0.0 ) then
168 call da_transform_through_rf(grid, mz, be % v9 % rf_alpha, be % v9 % val, vv % v9)
171 if ( mz > 0 .and. len_scaling10(1) > 0.0 ) then
172 call da_transform_through_rf(grid, mz, be % v10 % rf_alpha, be % v10 % val, vv % v10)
177 if ( mz > 0 .and. len_scaling6(1) > 0.0 ) then
178 call da_transform_through_rf(grid, mz, be % v6 % rf_alpha, be % v6 % val, vv % v6, scaling)
181 if ( mz > 0 .and. len_scaling7(1) > 0.0 ) then
182 call da_transform_through_rf(grid, mz, be % v7 % rf_alpha, be % v7 % val, vv % v7, scaling)
185 if ( mz > 0 .and. len_scaling8(1) > 0.0 ) then
186 call da_transform_through_rf(grid, mz, be % v8 % rf_alpha, be % v8 % val, vv % v8, scaling)
189 if ( mz > 0 .and. len_scaling9(1) > 0.0 ) then
190 call da_transform_through_rf(grid, mz, be % v9 % rf_alpha, be % v9 % val, vv % v9, scaling)
193 if ( mz > 0 .and. len_scaling10(1) > 0.0 ) then
194 call da_transform_through_rf(grid, mz, be % v10 % rf_alpha, be % v10 % val, vv % v10, scaling)
199 ! [2.7] Transform w control variables
202 if ( .not. use_cv_w ) then
206 if ( mz > 0 .and. len_scaling11(1) > 0.0 ) then
207 if ( cloud_cv_options == 2 ) then
208 call da_transform_through_rf(grid, mz, be % v11 % rf_alpha, be % v11 % val, vv % v11)
209 else if ( cloud_cv_options == 3 ) then
211 call da_transform_through_rf(grid, mz, be % v11 % rf_alpha, be % v11 % val, vv % v11, scaling)
217 ! [2.8] Transform alpha control variable
224 if ( anal_type_hybrid_dual_res ) then
225 call da_transform_through_rf_dual_res(grid%intermediate_grid, mz, be % alpha % rf_alpha, &
226 be % alpha % val, vv % alpha(:,:,:,n) )
228 call da_transform_through_rf(grid, mz, be % alpha % rf_alpha, be % alpha % val, vv % alpha(:,:,:,n) )
233 s(4)=s(2)+nij(0,0,2)*nij(0,1,2)*mz-1
234 call da_transform_through_wavelet(grid,mz,be%alpha%wsd,cv(s(2):s(4)),vv%alpha(:,:,:,n))
238 if( do_normalize )then
240 call da_transform_rescale(mz,be%alpha%sd,vv%alpha(:,:,:,n))
246 ! [2.10] Transform chem_ic control variables
248 if (present(vchem)) then
249 do ic = PARAM_FIRST_SCALAR, num_chem
250 mz = be % v12 (ic-1) % mz
251 if( use_rf .and. mz > 0 ) then
252 call da_transform_through_rf(grid, mz, be % v12 (ic-1) % rf_alpha, &
253 be % v12 (ic-1) % val, &
254 vchem % chem_ic (:,:,:,ic))
256 elseif( .not. use_rf ) then
257 call da_error(__FILE__,__LINE__,(/"no da_transform_through_wavelet for chem_ic"/))
264 if (trace_use) call da_trace_exit("da_transform_vtovv")
266 endsubroutine da_transform_vtovv