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