Merge remote-tracking branch 'origin/release-v4.5.2'
[WRF.git] / var / da / da_par_util / da_cv_to_vv.inc
blob11b076d530a9d55cbc9b75bee4faf5beea58f25a
1 subroutine da_cv_to_vv (cv_size, rcv, mzs, nmzs, vv &
2 #if (WRF_CHEM == 1)
3                         ,mzs_chem, vchem &
4 #endif
5                         )
6    !---------------------------------------------------------------------------
7    ! Purpose: Fill (local) vv arrays from 1D (local) cv array.
8    !---------------------------------------------------------------------------
10    implicit none
12    integer,        intent(in)    :: cv_size     ! Length of CV array.
13    real,           intent(in)    :: rcv(1:cv_size) ! Control variables v.
14    integer,        intent(in)    :: nmzs        ! number of cv for mz
15    integer,        intent(in)    :: mzs(nmzs)   ! Background error structure levels.
16    type (vp_type), intent(inout) :: vv          ! Grdipt/EOF cv_array.
17 #if (WRF_CHEM == 1)
18    integer, optional, intent(in)    :: mzs_chem(:)   ! chem Background error structure levels.
19    type (xchem_type), optional, intent(inout) :: vchem   ! Grdipt/EOF cv_array.
20 #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
32 #endif
34    if (trace_use) call da_trace_entry("da_cv_to_vv")
36    if( use_rf )then
37       is = its
38       ie = ite
39       js = jts
40       je = jte
41    else
42       call da_error(__FILE__,__LINE__,(/"This subroutine should not be called for use_rf = .false."/))
43    endif
44    ix = ie-is+1
45    jy = je-js+1
46    cv_e = 0
48    !--------------------------------------------------------------------------
49    ! [1] Transfer components of Jb control variable:
50    !--------------------------------------------------------------------------
52    ! Fill v1
53    mz = mzs(1)
54    if (mz > 0) then
55       ijm = ix * jy * mz
56       cv_s = cv_e + 1
57       cv_e = cv_s + ijm - 1
58       vv % v1(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
59    end if
61    ! Fill v2
62    mz = mzs(2)
63    if (mz > 0) then
64       ijm = ix * jy * mz
65       cv_s = cv_e + 1
66       cv_e = cv_s + ijm - 1
67       vv % v2(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
68    end if
70    ! Fill v3
71    mz = mzs(3)
72    if (mz > 0) then
73       ijm = ix * jy * mz
74       cv_s = cv_e + 1
75       cv_e = cv_s + ijm - 1
76       vv % v3(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
77    end if
79    ! Fill v4
80    mz = mzs(4)
81    if (mz > 0) then
82       ijm = ix * jy * mz
83       cv_s = cv_e + 1
84       cv_e = cv_s + ijm - 1
85       vv % v4(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
86    end if
88    ! Fill v5
89    mz = mzs(5)
90    if (mz == 1) then ! Can only be 0 or 1 (2D ps_u field)
91       ijm = ix * jy * mz
92       cv_s = cv_e + 1
93       cv_e = cv_s + ijm - 1
94       vv % v5(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy,mz/))
95    end if
97    if ( nmzs > 10 ) then
98       ! Fill v6
99       mz = mzs(6)
100       if (mz > 0) then
101          ijm = ix * jy * mz
102          cv_s = cv_e + 1
103          cv_e = cv_s + ijm - 1
104          vv % v6(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
105       end if
106       ! Fill v7
107       mz = mzs(7)
108       if (mz > 0) then
109          ijm = ix * jy * mz
110          cv_s = cv_e + 1
111          cv_e = cv_s + ijm - 1
112          vv % v7(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
113       end if
114       ! Fill v8
115       mz = mzs(8)
116       if (mz > 0) then
117          ijm = ix * jy * mz
118          cv_s = cv_e + 1
119          cv_e = cv_s + ijm - 1
120          vv % v8(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
121       end if
122       ! Fill v9
123       mz = mzs(9)
124       if (mz > 0) then
125          ijm = ix * jy * mz
126          cv_s = cv_e + 1
127          cv_e = cv_s + ijm - 1
128          vv % v9(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
129       end if
130       ! Fill v10
131       mz = mzs(10)
132       if (mz > 0) then
133          ijm = ix * jy * mz
134          cv_s = cv_e + 1
135          cv_e = cv_s + ijm - 1
136          vv % v10(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
137       end if
138    end if ! cloud cv
140    if ( nmzs == 8 .or. nmzs == 13 ) then
141       ! Fill v11
142       mz = mzs(nmzs-2)
143       if (mz > 0) then
144          ijm = ix * jy * mz
145          cv_s = cv_e + 1
146          cv_e = cv_s + ijm - 1
147          vv % v11(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
148       end if
149    end if ! w
151    !--------------------------------------------------------------------------
152    ! [2] Transfer components of Je control variable:
153    !--------------------------------------------------------------------------
154    mz = mzs(nmzs-1)
155    ne = mzs(nmzs)
156    if ( ne > 0 ) then
157       ix = ite_int - its_int + 1
158       jy = jte_int - jts_int + 1
159       ijmn = ix * jy * mz * ne
160       cv_s = cv_e + 1
161       cv_e = cv_s + ijmn - 1
162 !     vv % alpha(is:ie,js:je,1:mz,1:ne) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz, ne/))
163       vv % alpha(its_int:ite_int,jts_int:jte_int,1:mz,1:ne) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz, ne/))
164    end if
166 #if (WRF_CHEM == 1)
167    if (present(mzs_chem) .and. present(vchem)) then
168       ix = ie-is+1
169       jy = je-js+1
171       if ( num_chem >= PARAM_FIRST_SCALAR ) then
172          do ic = PARAM_FIRST_SCALAR, num_chem
173             mz = mzs_chem(ic-1)
174             ijm = ix * jy * mz
175             cv_s = cv_e + 1
176             cv_e = cv_s + ijm - 1
177             vchem % chem_ic(is:ie,js:je,1:mz,ic) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
178          end do
179       end if
180    end if
182 #endif
184    if (trace_use) call da_trace_exit("da_cv_to_vv")
186 end subroutine da_cv_to_vv