da_cv_to_vv.inc

References to this file elsewhere.
1 subroutine da_cv_to_vv (cv_size, rcv, xp, mzs, vv)
2 
3    !---------------------------------------------------------------------------
4    ! Purpose: Fill (local) vv arrays from 1D (local) cv array.
5    !---------------------------------------------------------------------------
6 
7    implicit none
8 
9    integer,        intent(in)    :: cv_size     ! Length of CV array.
10    real,           intent(in)    :: rcv(1:cv_size) ! Control variables v.
11    type (xpose_type), intent(in) :: xp          ! Dimensions and xpose buffers.
12    integer,        intent(in)    :: mzs(:) ! Background error structure levels.
13    type (vp_type), intent(inout) :: vv          ! Grdipt/EOF cv_array.
14 
15    integer   :: is,ie       ! Local grid range in y coordinate.
16    integer   :: js,je       ! Local grid range in x coordinate.
17    integer   :: ix,jy       ! Local grid horizontal dimensions.
18    integer   :: mz          ! Max vertical coordinate for v1 through v5 arrays.
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 
22    if (trace_use) call da_trace_entry("da_cv_to_vv")
23 
24    is = xp % its
25    ie = xp % ite
26    js = xp % jts
27    je = xp % jte
28    ix = ie-is+1
29    jy = je-js+1
30    cv_e = 0
31 
32    !--------------------------------------------------------------------------
33    ! [1] Transfer components of Jb control variable:
34    !--------------------------------------------------------------------------
35 
36    ! Fill v1
37    mz = mzs(1)
38    if (mz > 0) then
39       ijm = ix * jy * mz
40       cv_s = cv_e + 1
41       cv_e = cv_s + ijm - 1
42       vv % v1(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
43    end if
44 
45    ! Fill v2
46    mz = mzs(2)
47    if (mz > 0) then
48       ijm = ix * jy * mz
49       cv_s = cv_e + 1
50       cv_e = cv_s + ijm - 1
51       vv % v2(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
52    end if
53 
54    ! Fill v3
55    mz = mzs(3)
56    if (mz > 0) then
57       ijm = ix * jy * mz
58       cv_s = cv_e + 1
59       cv_e = cv_s + ijm - 1
60       vv % v3(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
61    end if
62 
63    ! Fill v4
64    mz = mzs(4)
65    if (mz > 0) then
66      ijm = ix * jy * mz
67      cv_s = cv_e + 1
68      cv_e = cv_s + ijm - 1
69      vv % v4(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
70    end if
71 
72    ! Fill v5
73    mz = mzs(5)
74    if (mz == 1) then ! Can only be 0 or 1 (2D ps_u field)
75      ijm = ix * jy * mz
76      cv_s = cv_e + 1
77      cv_e = cv_s + ijm - 1
78      vv % v5(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy,mz/))
79    end if
80 
81    !--------------------------------------------------------------------------
82    ! [2] Transfer components of Je control variable:
83    !--------------------------------------------------------------------------
84 
85    mz = mzs(6)
86    if (mz > 0) then
87       ijm = ix * jy * mz
88       cv_s = cv_e + 1
89       cv_e = cv_s + ijm - 1
90       vv % alpha(is:ie,js:je,1:mz) = RESHAPE(rcv(cv_s:cv_e),(/ix, jy, mz/))
91    end if
92 
93    if (trace_use) call da_trace_exit("da_cv_to_vv")
94 
95 end subroutine da_cv_to_vv
96 
97