wrf svn trunk commit r4103
[wrffire.git] / wrfv2_fire / var / da / da_par_util / da_cv_to_vv.inc
blobd9676e29a91a05df14af921829957395645be7b7
1 subroutine da_cv_to_vv (cv_size, rcv, mzs, vv)
3    !---------------------------------------------------------------------------
4    ! Purpose: Fill (local) vv arrays from 1D (local) cv array.
5    !---------------------------------------------------------------------------
7    implicit none
9    integer,        intent(in)    :: cv_size     ! Length of CV array.
10    real,           intent(in)    :: rcv(1:cv_size) ! Control variables v.
11    integer,        intent(in)    :: mzs(:) ! Background error structure levels.
12    type (vp_type), intent(inout) :: vv          ! Grdipt/EOF cv_array.
14    integer   :: is,ie       ! Local grid range in y coordinate.
15    integer   :: js,je       ! Local grid range in x coordinate.
16    integer   :: ix,jy       ! Local grid horizontal dimensions.
17    integer   :: mz          ! Max vertical coordinate for v1 through v5 arrays.
18    integer   :: ne          ! Ensemble size.
19    integer   :: cv_s,cv_e   ! Starting and ending indices into CV array.
20    integer   :: ijm         ! Size of interior of v1 through v5 arrays.
21    integer   :: ijmn        ! Size of interior of alpha cv arrays.
23    if (trace_use) call da_trace_entry("da_cv_to_vv")
25    is = its
26    ie = ite
27    js = jts
28    je = jte
29    ix = ie-is+1
30    jy = je-js+1
31    cv_e = 0
33    !--------------------------------------------------------------------------
34    ! [1] Transfer components of Jb control variable:
35    !--------------------------------------------------------------------------
37    ! Fill v1
38    mz = mzs(1)
39    if (mz > 0) then
40       ijm = ix * jy * mz
41       cv_s = cv_e + 1
42       cv_e = cv_s + ijm - 1
43       vv % v1(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
44    end if
46    ! Fill v2
47    mz = mzs(2)
48    if (mz > 0) then
49       ijm = ix * jy * mz
50       cv_s = cv_e + 1
51       cv_e = cv_s + ijm - 1
52       vv % v2(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
53    end if
55    ! Fill v3
56    mz = mzs(3)
57    if (mz > 0) then
58       ijm = ix * jy * mz
59       cv_s = cv_e + 1
60       cv_e = cv_s + ijm - 1
61       vv % v3(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
62    end if
64    ! Fill v4
65    mz = mzs(4)
66    if (mz > 0) then
67      ijm = ix * jy * mz
68      cv_s = cv_e + 1
69      cv_e = cv_s + ijm - 1
70      vv % v4(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
71    end if
73    ! Fill v5
74    mz = mzs(5)
75    if (mz == 1) then ! Can only be 0 or 1 (2D ps_u field)
76      ijm = ix * jy * mz
77      cv_s = cv_e + 1
78      cv_e = cv_s + ijm - 1
79      vv % v5(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy,mz/))
80    end if
82    !--------------------------------------------------------------------------
83    ! [2] Transfer components of Je control variable:
84    !--------------------------------------------------------------------------
86    mz = mzs(6)
87    ne = mzs(7)
88    if ( ne > 0 ) then
89       ijmn = ix * jy * mz * ne
90       cv_s = cv_e + 1
91       cv_e = cv_s + ijmn - 1
92       vv % alpha(is:ie,js:je,1:mz,1:ne) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz, ne/))
93    end if
95    if (trace_use) call da_trace_exit("da_cv_to_vv")
97 end subroutine da_cv_to_vv