1 |
jmc |
1.11 |
C $Header: $ |
2 |
|
|
C $Name: $ |
3 |
heimbach |
1.2 |
|
4 |
|
|
#include "CTRL_CPPOPTIONS.h" |
5 |
|
|
|
6 |
|
|
subroutine ctrl_set_pack_yz( |
7 |
heimbach |
1.3 |
& cunit, ivartype, fname, masktype, weighttype, |
8 |
heimbach |
1.2 |
& weightfld, lxxadxx, mythid) |
9 |
|
|
|
10 |
|
|
c ================================================================== |
11 |
|
|
c SUBROUTINE ctrl_set_pack_yz |
12 |
|
|
c ================================================================== |
13 |
|
|
c |
14 |
|
|
c o Compress the control vector such that only ocean points are |
15 |
|
|
c written to file. |
16 |
|
|
c |
17 |
heimbach |
1.3 |
c o Open boundary packing added : |
18 |
|
|
c gebbie@mit.edu, 18-Mar-2003 |
19 |
|
|
c |
20 |
|
|
c changed: heimbach@mit.edu 17-Jun-2003 |
21 |
|
|
c merged Armin's changes to replace write of |
22 |
|
|
c nr * globfld2d by 1 * globfld3d |
23 |
|
|
c (ad hoc fix to speed up global I/O) |
24 |
|
|
c |
25 |
heimbach |
1.2 |
c ================================================================== |
26 |
|
|
|
27 |
|
|
implicit none |
28 |
|
|
|
29 |
|
|
c == global variables == |
30 |
|
|
|
31 |
|
|
#include "EEPARAMS.h" |
32 |
|
|
#include "SIZE.h" |
33 |
|
|
#include "PARAMS.h" |
34 |
|
|
#include "GRID.h" |
35 |
|
|
|
36 |
|
|
#include "ctrl.h" |
37 |
|
|
#include "optim.h" |
38 |
|
|
|
39 |
|
|
c == routine arguments == |
40 |
|
|
|
41 |
|
|
integer cunit |
42 |
|
|
integer ivartype |
43 |
|
|
character*( 80) fname |
44 |
|
|
character*( 9) masktype |
45 |
heimbach |
1.3 |
character*( 80) weighttype |
46 |
heimbach |
1.2 |
_RL weightfld( nr,nobcs ) |
47 |
|
|
logical lxxadxx |
48 |
|
|
integer mythid |
49 |
|
|
|
50 |
heimbach |
1.10 |
#ifndef EXCLUDE_CTRL_PACK |
51 |
heimbach |
1.2 |
c == local variables == |
52 |
|
|
|
53 |
|
|
integer bi,bj |
54 |
|
|
integer ip,jp |
55 |
|
|
integer i,j,k |
56 |
heimbach |
1.7 |
integer ii,jj,kk |
57 |
heimbach |
1.2 |
integer il |
58 |
heimbach |
1.3 |
integer irec,iobcs,nrec_nl |
59 |
heimbach |
1.2 |
integer itlo,ithi |
60 |
|
|
integer jtlo,jthi |
61 |
|
|
integer jmin,jmax |
62 |
|
|
integer imin,imax |
63 |
|
|
|
64 |
|
|
integer cbuffindex |
65 |
|
|
cgg( |
66 |
|
|
integer igg |
67 |
|
|
_RL gg |
68 |
heimbach |
1.3 |
character*(80) weightname |
69 |
heimbach |
1.2 |
cgg) |
70 |
heimbach |
1.7 |
real*4 cbuff ( nsx*npx*sny*nsy*npy ) |
71 |
heimbach |
1.2 |
_RL globfldyz ( nsx,npx,sny,nsy,npy,nr ) |
72 |
heimbach |
1.3 |
_RL globfld3d ( snx,nsx,npx,sny,nsy,npy,nr ) |
73 |
|
|
_RL globmskyz ( nsx,npx,sny,nsy,npy,nr,nobcs ) |
74 |
|
|
#ifdef CTRL_PACK_PRECISE |
75 |
|
|
_RL weightfldyz( nsx,npx,sny,nsy,npy,nr,nobcs ) |
76 |
|
|
#endif |
77 |
heimbach |
1.2 |
|
78 |
|
|
c == external == |
79 |
|
|
|
80 |
|
|
integer ilnblnk |
81 |
|
|
external ilnblnk |
82 |
|
|
|
83 |
|
|
c == end of interface == |
84 |
|
|
|
85 |
|
|
jtlo = 1 |
86 |
|
|
jthi = nsy |
87 |
|
|
itlo = 1 |
88 |
|
|
ithi = nsx |
89 |
|
|
jmin = 1 |
90 |
|
|
jmax = sny |
91 |
|
|
imin = 1 |
92 |
|
|
imax = snx |
93 |
|
|
|
94 |
|
|
c Initialise temporary file |
95 |
|
|
do k = 1,nr |
96 |
|
|
do jp = 1,nPy |
97 |
|
|
do bj = jtlo,jthi |
98 |
|
|
do j = jmin,jmax |
99 |
|
|
do ip = 1,nPx |
100 |
|
|
do bi = itlo,ithi |
101 |
|
|
globfldyz(bi,ip,j,bj,jp,k) = 0. _d 0 |
102 |
heimbach |
1.3 |
do iobcs=1,nobcs |
103 |
|
|
globmskyz(bi,ip,j,bj,jp,k,iobcs) = 0. _d 0 |
104 |
|
|
enddo |
105 |
|
|
enddo |
106 |
|
|
enddo |
107 |
|
|
enddo |
108 |
|
|
enddo |
109 |
|
|
enddo |
110 |
|
|
enddo |
111 |
|
|
c Initialise temporary file |
112 |
|
|
do k = 1,nr |
113 |
|
|
do jp = 1,nPy |
114 |
|
|
do bj = jtlo,jthi |
115 |
|
|
do j = jmin,jmax |
116 |
|
|
do ip = 1,nPx |
117 |
|
|
do bi = itlo,ithi |
118 |
|
|
do i = imin,imax |
119 |
|
|
globfld3d(i,bi,ip,j,bj,jp,k) = 0. _d 0 |
120 |
|
|
enddo |
121 |
heimbach |
1.2 |
enddo |
122 |
|
|
enddo |
123 |
|
|
enddo |
124 |
|
|
enddo |
125 |
|
|
enddo |
126 |
|
|
enddo |
127 |
|
|
|
128 |
|
|
c-- Only the master thread will do I/O. |
129 |
|
|
_BEGIN_MASTER( mythid ) |
130 |
|
|
|
131 |
heimbach |
1.3 |
do iobcs=1,nobcs |
132 |
jmc |
1.11 |
call MDSREADFIELD_YZ_GL( |
133 |
heimbach |
1.3 |
& masktype, ctrlprec, 'RL', |
134 |
|
|
& Nr, globmskyz(1,1,1,1,1,1,iobcs), iobcs,mythid) |
135 |
|
|
#ifdef CTRL_PACK_PRECISE |
136 |
|
|
il=ilnblnk( weighttype) |
137 |
|
|
write(weightname(1:80),'(80a)') ' ' |
138 |
|
|
write(weightname(1:80),'(a)') weighttype(1:il) |
139 |
|
|
call MDSREADFIELD_YZ_GL( |
140 |
|
|
& weightname, ctrlprec, 'RL', |
141 |
|
|
& Nr, weightfldyz(1,1,1,1,1,1,iobcs), iobcs, mythid) |
142 |
|
|
CGG One special exception: barotropic velocity should be nondimensionalized |
143 |
|
|
cgg differently. Probably introduce new variable. |
144 |
|
|
if (iobcs .eq. 3 .or. iobcs .eq. 4) then |
145 |
|
|
k = 1 |
146 |
|
|
do jp = 1,nPy |
147 |
|
|
do bj = jtlo,jthi |
148 |
|
|
do j = jmin,jmax |
149 |
|
|
do ip = 1,nPx |
150 |
|
|
do bi = itlo,ithi |
151 |
heimbach |
1.9 |
cph weightfldyz(bi,ip,j,bj,jp,k,iobcs) = wbaro |
152 |
heimbach |
1.3 |
enddo |
153 |
|
|
enddo |
154 |
|
|
enddo |
155 |
|
|
enddo |
156 |
|
|
enddo |
157 |
|
|
endif |
158 |
|
|
#endif |
159 |
|
|
enddo |
160 |
heimbach |
1.2 |
|
161 |
heimbach |
1.3 |
nrec_nl=int(ncvarrecs(ivartype)/snx) |
162 |
|
|
do irec = 1, nrec_nl |
163 |
heimbach |
1.2 |
cgg do iobcs = 1, nobcs |
164 |
|
|
cgg Need to solve for what iobcs would have been. |
165 |
|
|
|
166 |
heimbach |
1.3 |
call MDSREADFIELD_3D_GL( fname, ctrlprec, 'RL', |
167 |
|
|
& nr, globfld3D, irec, mythid) |
168 |
|
|
|
169 |
|
|
do i=1,snx |
170 |
|
|
iobcs= mod((irec-1)*snx+i-1,nobcs)+1 |
171 |
|
|
|
172 |
|
|
CGG One special exception: barotropic velocity should be nondimensionalized |
173 |
|
|
cgg differently. Probably introduce new variable. |
174 |
|
|
if (iobcs .eq. 3 .or. iobcs .eq. 4) then |
175 |
|
|
k = 1 |
176 |
|
|
do jp = 1,nPy |
177 |
|
|
do bj = jtlo,jthi |
178 |
|
|
do j = jmin,jmax |
179 |
|
|
do ip = 1,nPx |
180 |
|
|
do bi = itlo,ithi |
181 |
|
|
#ifdef NO_CONTROL_BAROTROPIC_VELOCITY |
182 |
|
|
if (.not. lxxadxx) then |
183 |
|
|
cgg Get rid of any sensitivity to barotropic velocity. |
184 |
|
|
globfld3d(i,bi,ip,j,bj,jp,k) = 0. |
185 |
|
|
endif |
186 |
|
|
#endif |
187 |
|
|
enddo |
188 |
|
|
enddo |
189 |
|
|
enddo |
190 |
|
|
enddo |
191 |
|
|
enddo |
192 |
|
|
endif |
193 |
|
|
|
194 |
|
|
write(cunit) ncvarindex(ivartype) |
195 |
|
|
write(cunit) 1 |
196 |
|
|
write(cunit) 1 |
197 |
|
|
do k = 1,nr |
198 |
|
|
cbuffindex = 0 |
199 |
|
|
do jp = 1,nPy |
200 |
|
|
do bj = jtlo,jthi |
201 |
|
|
do ip = 1,nPx |
202 |
|
|
do bi = itlo,ithi |
203 |
|
|
do j = jmin,jmax |
204 |
heimbach |
1.7 |
ii=mod((i-1)*nr*sny+(k-1)*sny+j-1,snx)+1 |
205 |
|
|
jj=mod(((i-1)*nr*sny+(k-1)*sny+j-1)/snx,sny)+1 |
206 |
|
|
kk=int((i-1)*nr*sny+(k-1)*sny+j-1)/(snx*sny)+1 |
207 |
heimbach |
1.3 |
if (globmskyz(bi,ip,j,bj,jp,k,iobcs) .ne. 0. ) then |
208 |
|
|
cbuffindex = cbuffindex + 1 |
209 |
|
|
#ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO |
210 |
|
|
if (lxxadxx) then |
211 |
jmc |
1.11 |
cbuff(cbuffindex) = |
212 |
heimbach |
1.7 |
& globfld3d(ii,bi,ip,jj,bj,jp,kk) * |
213 |
heimbach |
1.3 |
#ifdef CTRL_PACK_PRECISE |
214 |
|
|
& sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs)) |
215 |
|
|
#else |
216 |
|
|
& sqrt(weightfld(k,iobcs)) |
217 |
|
|
#endif |
218 |
|
|
else |
219 |
jmc |
1.11 |
cbuff(cbuffindex) = |
220 |
heimbach |
1.7 |
& globfld3d(ii,bi,ip,jj,bj,jp,kk) / |
221 |
heimbach |
1.3 |
#ifdef CTRL_PACK_PRECISE |
222 |
|
|
& sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs)) |
223 |
|
|
#else |
224 |
|
|
& sqrt(weightfld(k,iobcs)) |
225 |
|
|
#endif |
226 |
|
|
endif |
227 |
|
|
#else /* ALLOW_NONDIMENSIONAL_CONTROL_IO undef */ |
228 |
heimbach |
1.7 |
cbuff(cbuffindex) = globfld3d(ii,bi,ip,jj,bj,jp,kk) |
229 |
heimbach |
1.3 |
#endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */ |
230 |
|
|
endif |
231 |
|
|
enddo |
232 |
|
|
enddo |
233 |
|
|
enddo |
234 |
|
|
enddo |
235 |
|
|
enddo |
236 |
|
|
c --> check cbuffindex. |
237 |
|
|
if ( cbuffindex .gt. 0) then |
238 |
|
|
write(cunit) cbuffindex |
239 |
|
|
write(cunit) k |
240 |
|
|
write(cunit) (cbuff(ii), ii=1,cbuffindex) |
241 |
|
|
endif |
242 |
|
|
c |
243 |
|
|
c -- end of k loop -- |
244 |
|
|
enddo |
245 |
|
|
c -- end of iobcs loop -- |
246 |
|
|
cgg enddo |
247 |
|
|
c -- end of i loop -- |
248 |
|
|
enddo |
249 |
|
|
c -- end of irec loop -- |
250 |
|
|
enddo |
251 |
|
|
|
252 |
|
|
do irec = nrec_nl*snx+1, ncvarrecs(ivartype) |
253 |
|
|
cgg do iobcs = 1, nobcs |
254 |
|
|
cgg Need to solve for what iobcs would have been. |
255 |
|
|
iobcs= mod(irec-1,nobcs)+1 |
256 |
heimbach |
1.2 |
|
257 |
|
|
call MDSREADFIELD_YZ_GL( fname, ctrlprec, 'RL', |
258 |
|
|
& nr, globfldyz, irec, mythid) |
259 |
|
|
|
260 |
heimbach |
1.3 |
CGG One special exception: barotropic velocity should be nondimensionalized |
261 |
|
|
cgg differently. Probably introduce new variable. |
262 |
|
|
if (iobcs .eq. 3 .or. iobcs .eq. 4) then |
263 |
|
|
k = 1 |
264 |
|
|
do jp = 1,nPy |
265 |
|
|
do bj = jtlo,jthi |
266 |
|
|
do j = jmin,jmax |
267 |
|
|
do ip = 1,nPx |
268 |
|
|
do bi = itlo,ithi |
269 |
|
|
#ifdef NO_CONTROL_BAROTROPIC_VELOCITY |
270 |
|
|
if (.not. lxxadxx) then |
271 |
|
|
cgg Get rid of any sensitivity to barotropic velocity. |
272 |
|
|
globfldyz(bi,ip,j,bj,jp,k) = 0. |
273 |
|
|
endif |
274 |
|
|
#endif |
275 |
|
|
enddo |
276 |
|
|
enddo |
277 |
|
|
enddo |
278 |
|
|
enddo |
279 |
|
|
enddo |
280 |
|
|
endif |
281 |
|
|
|
282 |
heimbach |
1.2 |
write(cunit) ncvarindex(ivartype) |
283 |
|
|
write(cunit) 1 |
284 |
|
|
write(cunit) 1 |
285 |
|
|
do k = 1,nr |
286 |
|
|
cbuffindex = 0 |
287 |
|
|
do jp = 1,nPy |
288 |
|
|
do bj = jtlo,jthi |
289 |
|
|
do ip = 1,nPx |
290 |
|
|
do bi = itlo,ithi |
291 |
|
|
do j = jmin,jmax |
292 |
heimbach |
1.3 |
if (globmskyz(bi,ip,j,bj,jp,k,iobcs) .ne. 0. ) then |
293 |
heimbach |
1.2 |
cbuffindex = cbuffindex + 1 |
294 |
|
|
#ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO |
295 |
|
|
if (lxxadxx) then |
296 |
jmc |
1.11 |
cbuff(cbuffindex) = |
297 |
heimbach |
1.2 |
& globfldyz(bi,ip,j,bj,jp,k) * |
298 |
heimbach |
1.3 |
#ifdef CTRL_PACK_PRECISE |
299 |
|
|
& sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs)) |
300 |
|
|
#else |
301 |
heimbach |
1.2 |
& sqrt(weightfld(k,iobcs)) |
302 |
heimbach |
1.3 |
#endif |
303 |
heimbach |
1.2 |
else |
304 |
jmc |
1.11 |
cbuff(cbuffindex) = |
305 |
heimbach |
1.2 |
& globfldyz(bi,ip,j,bj,jp,k) / |
306 |
heimbach |
1.3 |
#ifdef CTRL_PACK_PRECISE |
307 |
|
|
& sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs)) |
308 |
|
|
#else |
309 |
heimbach |
1.2 |
& sqrt(weightfld(k,iobcs)) |
310 |
heimbach |
1.3 |
#endif |
311 |
heimbach |
1.2 |
endif |
312 |
heimbach |
1.3 |
#else /* ALLOW_NONDIMENSIONAL_CONTROL_IO undef */ |
313 |
heimbach |
1.2 |
cbuff(cbuffindex) = globfldyz(bi,ip,j,bj,jp,k) |
314 |
heimbach |
1.3 |
#endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */ |
315 |
heimbach |
1.2 |
endif |
316 |
|
|
enddo |
317 |
|
|
enddo |
318 |
|
|
enddo |
319 |
|
|
enddo |
320 |
|
|
enddo |
321 |
|
|
c --> check cbuffindex. |
322 |
|
|
if ( cbuffindex .gt. 0) then |
323 |
|
|
write(cunit) cbuffindex |
324 |
|
|
write(cunit) k |
325 |
|
|
write(cunit) (cbuff(ii), ii=1,cbuffindex) |
326 |
|
|
endif |
327 |
heimbach |
1.3 |
c |
328 |
|
|
c -- end of k loop -- |
329 |
heimbach |
1.2 |
enddo |
330 |
|
|
c -- end of iobcs loop -- |
331 |
|
|
cgg enddo |
332 |
|
|
c -- end of irec loop -- |
333 |
|
|
enddo |
334 |
|
|
|
335 |
|
|
_END_MASTER( mythid ) |
336 |
|
|
|
337 |
heimbach |
1.10 |
#endif |
338 |
|
|
|
339 |
heimbach |
1.2 |
return |
340 |
|
|
end |
341 |
|
|
|