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