Merge remote-tracking branch 'origin/release-v4.6.1'
[WRF.git] / var / da / da_par_util / da_vv_to_cv.inc
blob4d37c11b20bc80525977a4002ecbc855c2d508a9
1 subroutine da_vv_to_cv(vv, xp, mzs, nmzs, cv_size, rcv &
2 #if (WRF_CHEM == 1)
3                         , mzs_chem, vchem &
4 #endif
5                         )
6    !---------------------------------------------------------------------------
7    ! Purpose: Fill (local) 1D cv array from 3D (local) vv arrays.
8    !---------------------------------------------------------------------------
10    implicit none
12    type (vp_type), intent(in)    :: vv          ! Grdipt/EOF rcv.
13    type (xpose_type), intent(in) :: xp          ! Dimensions and xpose buffers.
14    integer,        intent(in)    :: nmzs        ! number of cv for mz
15    integer,        intent(in)    :: mzs(nmzs)   ! Background error structure levels.
16    integer,        intent(in)    :: cv_size     ! Length of CV array.
17    real,           intent(inout) :: rcv(1:cv_size) ! Control variables v.
18 #if (WRF_CHEM == 1)
19    integer, optional, intent(in)    :: mzs_chem(:)   ! chem Background error structure levels.
20    type (xchem_type), optional, intent(inout) :: vchem   ! Grdipt/EOF rcv for chem
21 #endif
22    integer   :: is,ie       ! Local grid range in y coordinate.
23    integer   :: js,je       ! Local grid range in x coordinate.
24    integer   :: ix,jy       ! Local grid horizontal dimensions.
25    integer   :: mz          ! Max vertical coordinate for v1 through v5 arrays.
26    integer   :: ne          ! Ensemble size.
27    integer   :: cv_s,cv_e   ! Starting and ending indices into CV array.
28    integer   :: ijm         ! Size of interior of v1 through v5 arrays.
29    integer   :: ijmn        ! Size of interior of alpha cv arrays.
30 #if (WRF_CHEM == 1)
31    integer   :: ic, offset,tt
32 #endif
34 #ifdef NORESHAPE
35    integer   :: i,j,k,ijk,m,n
36 #endif
38    if (trace_use) call da_trace_entry("da_vv_to_cv")
40    if( use_rf )then
41       is = xp % its
42       ie = xp % ite
43       js = xp % jts
44       je = xp % jte
45    else
46       call da_error(__FILE__,__LINE__,(/"This subroutine should not be called for use_rf = .false."/))
47    endif
48    ix = ie-is+1
49    jy = je-js+1
50    cv_e = 0
51 #ifdef NORESHAPE
52    ijk=0
53 #endif
55    ! Store v1
56    mz = mzs(1)
57    if (mz > 0) then
58       ijm = ix * jy * mz
59       cv_s = cv_e + 1
60       cv_e = cv_s + ijm - 1
61 #ifdef NORESHAPE
62       do k=1,mz
63          do j=js,je
64             do i=is,ie
65                ijk = ijk + 1
66                rcv(ijk) = vv%v1(i,j,k)
67             end do
68         end do
69      end do
70 #else
71      rcv(cv_s:cv_e) = RESHAPE(vv % v1(is:ie,js:je,1:mz), (/ijm/))
72 #endif
73    end if
74    ! Store v2
75    mz = mzs(2)
76    if (mz > 0) then
77       ijm = ix * jy * mz
78       cv_s = cv_e + 1
79       cv_e = cv_s + ijm - 1
80 #ifdef NORESHAPE
81       do k=1,mz
82          do j=js,je
83             do i=is,ie
84                ijk = ijk + 1
85                rcv(ijk) = vv%v2(i,j,k)
86             end do
87          end do
88       end do
89 #else
90       rcv(cv_s:cv_e) = RESHAPE(vv % v2(is:ie,js:je,1:mz), (/ijm/))
91 #endif
92    end if
93    ! Store v3
94    mz = mzs(3)
95    if (mz > 0) then
96       ijm = ix * jy * mz
97       cv_s = cv_e + 1
98       cv_e = cv_s + ijm - 1
100 #ifdef NORESHAPE
101       do k=1,mz
102          do j=js,je
103             do i=is,ie
104                ijk = ijk + 1
105                rcv(ijk) = vv%v3(i,j,k)
106             end do
107          end do
108       end do
109 #else
110       rcv(cv_s:cv_e) = RESHAPE(vv % v3(is:ie,js:je,1:mz), (/ijm/))
111 #endif
112    end if
113    ! Store v4
114    mz = mzs(4)
115    if (mz > 0) then
116       ijm = ix * jy * mz
117       cv_s = cv_e + 1
118       cv_e = cv_s + ijm - 1
119 #ifdef NORESHAPE
120       do k=1,mz
121          do j=js,je
122             do i=is,ie
123                ijk = ijk + 1
124                rcv(ijk) = vv%v4(i,j,k)
125             end do
126          end do
127       end do
128 #else
129       rcv(cv_s:cv_e) = RESHAPE(vv % v4(is:ie,js:je,1:mz), (/ijm/))
130 #endif
131    end if
132    ! Store v5
133    mz = mzs(5)
134    if (mz > 0) then
135       ijm = ix * jy * mz
136       cv_s = cv_e + 1
137       cv_e = cv_s + ijm - 1
138 #ifdef NORESHAPE
139       do k=1,mz
140          do j=js,je
141             do i=is,ie
142                ijk = ijk + 1
143                rcv(ijk) = vv%v5(i,j,k)
144             end do
145          end do
146       end do
147 #else
148       rcv(cv_s:cv_e) = RESHAPE(vv % v5(is:ie,js:je,1:mz), (/ijm/))
149 #endif
150    end if
151    if ( nmzs > 10 ) then
152       ! Store v6
153       mz = mzs(6)
154       if (mz > 0) then
155          ijm = ix * jy * mz
156          cv_s = cv_e + 1
157          cv_e = cv_s + ijm - 1
158 #ifdef NORESHAPE
159          do k=1,mz
160             do j=js,je
161                do i=is,ie
162                   ijk = ijk + 1
163                   rcv(ijk) = vv%v6(i,j,k)
164                end do
165             end do
166          end do
167 #else
168          rcv(cv_s:cv_e) = RESHAPE(vv % v6(is:ie,js:je,1:mz), (/ijm/))
169 #endif
170       end if
171       ! Store v7
172       mz = mzs(7)
173       if (mz > 0) then
174          ijm = ix * jy * mz
175          cv_s = cv_e + 1
176          cv_e = cv_s + ijm - 1
177 #ifdef NORESHAPE
178          do k=1,mz
179             do j=js,je
180                do i=is,ie
181                   ijk = ijk + 1
182                   rcv(ijk) = vv%v7(i,j,k)
183                end do
184             end do
185          end do
186 #else
187          rcv(cv_s:cv_e) = RESHAPE(vv % v7(is:ie,js:je,1:mz), (/ijm/))
188 #endif
189       end if
190       ! Store v8
191       mz = mzs(8)
192       if (mz > 0) then
193          ijm = ix * jy * mz
194          cv_s = cv_e + 1
195          cv_e = cv_s + ijm - 1
196 #ifdef NORESHAPE
197          do k=1,mz
198             do j=js,je
199                do i=is,ie
200                   ijk = ijk + 1
201                   rcv(ijk) = vv%v8(i,j,k)
202                end do
203             end do
204          end do
205 #else
206          rcv(cv_s:cv_e) = RESHAPE(vv % v8(is:ie,js:je,1:mz), (/ijm/))
207 #endif
208       end if
209       ! Store v9
210       mz = mzs(9)
211       if (mz > 0) then
212          ijm = ix * jy * mz
213          cv_s = cv_e + 1
214          cv_e = cv_s + ijm - 1
215 #ifdef NORESHAPE
216          do k=1,mz
217             do j=js,je
218                do i=is,ie
219                   ijk = ijk + 1
220                   rcv(ijk) = vv%v9(i,j,k)
221                end do
222             end do
223          end do
224 #else
225          rcv(cv_s:cv_e) = RESHAPE(vv % v9(is:ie,js:je,1:mz), (/ijm/))
226 #endif
227       end if
228       ! Store v10
229       mz = mzs(10)
230       if (mz > 0) then
231          ijm = ix * jy * mz
232          cv_s = cv_e + 1
233          cv_e = cv_s + ijm - 1
234 #ifdef NORESHAPE
235          do k=1,mz
236             do j=js,je
237                do i=is,ie
238                   ijk = ijk + 1
239                   rcv(ijk) = vv%v10(i,j,k)
240                end do
241             end do
242          end do
243 #else
244          rcv(cv_s:cv_e) = RESHAPE(vv % v10(is:ie,js:je,1:mz), (/ijm/))
245 #endif
246       end if
247    end if ! cloud_cv
248    if ( nmzs == 8 .or. nmzs == 13 ) then
249       ! Store v11
250       mz = mzs(nmzs-2)
251       if (mz > 0) then
252          ijm = ix * jy * mz
253          cv_s = cv_e + 1
254          cv_e = cv_s + ijm - 1
255 #ifdef NORESHAPE
256          do k=1,mz
257             do j=js,je
258                do i=is,ie
259                   ijk = ijk + 1
260                   rcv(ijk) = vv%v11(i,j,k)
261                end do
262             end do
263          end do
264 #else
265          rcv(cv_s:cv_e) = RESHAPE(vv % v11(is:ie,js:je,1:mz), (/ijm/))
266 #endif
267       end if
268    end if ! w
270    ! Store alpha:
271    mz = mzs(nmzs-1)
272    ne = mzs(nmzs)
273    if ( ne > 0 ) then
274       ix = ite_int - its_int + 1
275       jy = jte_int - jts_int + 1
276       ijmn = ix * jy * mz * ne
277       cv_s = cv_e + 1
278       cv_e = cv_s + ijmn - 1
279 #ifdef NORESHAPE
280       do n = 1, ne
281          do m = 1, mz
282             do j = jts_int, jte_int !js, je
283                do i = its_int, ite_int ! is, ie
284                   ijk = ijk + 1
285                   rcv(ijk) = vv%alpha(i,j,m,n)
286                end do
287             end do
288          end do
289       end do
290 #else
291 !     rcv(cv_s:cv_e) = RESHAPE(vv % alpha(is:ie,js:je,1:mz,1:ne), (/ijmn/))
292       rcv(cv_s:cv_e) = RESHAPE(vv % alpha(its_int:ite_int,jts_int:jte_int,1:mz,1:ne), (/ijmn/))
293 #endif
294    end if
296 #if (WRF_CHEM == 1)
297    if (present(mzs_chem) .and. present(vchem)) then
298       ix = ie-is+1
299       jy = je-js+1
300       offset = 1
302       if ( num_chem >= PARAM_FIRST_SCALAR ) then
303          do ic = PARAM_FIRST_SCALAR, num_chem
304             mz = mzs_chem(ic-1)
305             ijm = ix * jy * mz
306             cv_s = cv_e + 1
307             cv_e = cv_s + ijm - 1
308 #ifdef NORESHAPE
309             do k=1,mz
310                do j=js,je
311                   do i=is,ie
312                      ijk = ijk + 1
313                      rcv(ijk) = vchem % chem_ic(i,j,k,ic)
314                   end do
315                end do
316             end do
317 #else
318             rcv(cv_s:cv_e) = RESHAPE(vchem % chem_ic(is:ie,js:je,1:mz,ic))
319 #endif
320          end do
321       end if
323    end if
325 #endif
327    if (trace_use) call da_trace_exit("da_vv_to_cv")
329 end subroutine da_vv_to_cv