1 subroutine da_cv_to_vv (cv_size, rcv, mzs, nmzs, vv &
6 !---------------------------------------------------------------------------
7 ! Purpose: Fill (local) vv arrays from 1D (local) cv array.
8 !---------------------------------------------------------------------------
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.
18 integer, optional, intent(in) :: mzs_chem(:) ! chem Background error structure levels.
19 type (xchem_type), optional, intent(inout) :: vchem ! Grdipt/EOF cv_array.
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.
34 if (trace_use) call da_trace_entry("da_cv_to_vv")
42 call da_error(__FILE__,__LINE__,(/"This subroutine should not be called for use_rf = .false."/))
48 !--------------------------------------------------------------------------
49 ! [1] Transfer components of Jb control variable:
50 !--------------------------------------------------------------------------
58 vv % v1(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
67 vv % v2(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
76 vv % v3(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
85 vv % v4(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
90 if (mz == 1) then ! Can only be 0 or 1 (2D ps_u field)
94 vv % v5(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy,mz/))
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/))
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/))
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/))
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/))
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/))
140 if ( nmzs == 8 .or. nmzs == 13 ) then
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/))
151 !--------------------------------------------------------------------------
152 ! [2] Transfer components of Je control variable:
153 !--------------------------------------------------------------------------
157 ix = ite_int - its_int + 1
158 jy = jte_int - jts_int + 1
159 ijmn = ix * jy * mz * ne
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/))
167 if (present(mzs_chem) .and. present(vchem)) then
171 if ( num_chem >= PARAM_FIRST_SCALAR ) then
172 do ic = PARAM_FIRST_SCALAR, num_chem
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/))
184 if (trace_use) call da_trace_exit("da_cv_to_vv")
186 end subroutine da_cv_to_vv