da_vv_to_cv.inc

References to this file elsewhere.
1 subroutine da_vv_to_cv(vv, xp, mzs, cv_size, rcv)
2 
3    !---------------------------------------------------------------------------
4    ! Purpose: Fill (local) 1D cv array from 3D (local) vv arrays.
5    !---------------------------------------------------------------------------
6 
7    implicit none
8 
9    type (vp_type), intent(in)    :: vv          ! Grdipt/EOF rcv.
10    type (xpose_type), intent(in) :: xp          ! Dimensions and xpose buffers.
11    integer,        intent(in)    :: mzs(:) ! Background error structure levels.
12    integer,        intent(in)    :: cv_size     ! Length of CV array.
13    real,           intent(out)   :: rcv(1:cv_size) ! Control variables v.
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 #ifdef NORESHAPE
23    integer   :: i,j,k,ijk
24 #endif
25 
26    is = xp % its
27    ie = xp % ite
28    js = xp % jts
29    je = xp % jte
30    ix = ie-is+1
31    jy = je-js+1
32    cv_e = 0
33 
34    ! Store v1
35    mz = mzs(1)
36    if (mz > 0) then
37       ijm = ix * jy * mz
38       cv_s = cv_e + 1
39       cv_e = cv_s + ijm - 1
40 #ifdef NORESHAPE
41       ijk = 0
42       do k=1,mz
43          do j=js,je
44             do i=is,ie
45                ijk = ijk + 1
46                rcv(ijk) = vv%v1(i,j,k)
47             end do
48         end do
49      end do
50 #else
51      rcv(cv_s:cv_e) = RESHAPE(vv % v1(is:ie,js:je,1:mz), (/ijm/))
52 #endif
53    end if
54 
55    ! Store v2
56    mz = mzs(2)
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%v2(i,j,k)
67             end do
68         end do
69      end do
70 #else
71      rcv(cv_s:cv_e) = RESHAPE(vv % v2(is:ie,js:je,1:mz), (/ijm/))
72 #endif
73    end if
74 
75    ! Store v3
76    mz = mzs(3)
77    if (mz > 0) then
78       ijm = ix * jy * mz
79       cv_s = cv_e + 1
80       cv_e = cv_s + ijm - 1
81 
82 #ifdef NORESHAPE
83       do k=1,mz
84          do j=js,je
85             do i=is,ie
86                ijk = ijk + 1
87                rcv(ijk) = vv%v3(i,j,k)
88             end do
89          end do
90       end do
91 #else
92       rcv(cv_s:cv_e) = RESHAPE(vv % v3(is:ie,js:je,1:mz), (/ijm/))
93 #endif
94    end if
95 
96    ! Store v4
97    mz = mzs(4)
98    if (mz > 0) then
99       ijm = ix * jy * mz
100       cv_s = cv_e + 1
101       cv_e = cv_s + ijm - 1
102 #ifdef NORESHAPE
103       do k=1,mz
104          do j=js,je
105             do i=is,ie
106                ijk = ijk + 1
107                rcv(ijk) = vv%v4(i,j,k)
108             end do
109          end do
110       end do
111 #else
112       rcv(cv_s:cv_e) = RESHAPE(vv % v4(is:ie,js:je,1:mz), (/ijm/))
113 #endif
114    end if
115 
116    ! Store v5
117    mz = mzs(5)
118    if (mz > 0) then
119       ijm = ix * jy * mz
120       cv_s = cv_e + 1
121       cv_e = cv_s + ijm - 1
122 #ifdef NORESHAPE
123       do k=1,mz
124          do j=js,je
125             do i=is,ie
126                ijk = ijk + 1
127                rcv(ijk) = vv%v5(i,j,k)
128             end do
129          end do
130       end do
131 #else
132       rcv(cv_s:cv_e) = RESHAPE(vv % v5(is:ie,js:je,1:mz), (/ijm/))
133 #endif
134    end if
135 
136    ! Store alpha:
137    mz = mzs(6)
138    if (mz > 0) then
139       ijm = ix * jy * mz
140       cv_s = cv_e + 1
141       cv_e = cv_s + ijm - 1
142 #ifdef NORESHAPE
143       do k=1,mz
144          do j=js,je
145             do i=is,ie
146                ijk = ijk + 1
147                rcv(ijk) = vv%alpha(i,j,k)
148             end do
149          end do
150       end do
151 #else
152       rcv(cv_s:cv_e) = RESHAPE(vv % alpha(is:ie,js:je,1:mz), (/ijm/))
153 #endif
154    end if
155 
156 end subroutine da_vv_to_cv
157 
158