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