1 |
C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_set_pack_xz.F,v 1.16 2012/08/10 19:38:58 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "CTRL_OPTIONS.h" |
5 |
|
6 |
subroutine ctrl_set_pack_xz( |
7 |
& cunit, ivartype, fname, masktype,weighttype, |
8 |
& weightfld, lxxadxx, mythid) |
9 |
|
10 |
c ================================================================== |
11 |
c SUBROUTINE ctrl_set_pack_xz |
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 |
c o Open boundary packing finalized : |
18 |
c gebbie@mit.edu, 18-Mar-2003 |
19 |
c |
20 |
c changed: heimbach@mit.edu 17-Jun-2003 |
21 |
c merged changes from Armin 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 |
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 "CTRL_OBCS.h" |
38 |
#include "optim.h" |
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 |
#ifndef EXCLUDE_CTRL_PACK |
52 |
c == local variables == |
53 |
|
54 |
integer bi,bj |
55 |
integer ip,jp |
56 |
integer i,j,k |
57 |
integer ii,jj,kk |
58 |
integer irec,iobcs,nrec_nl |
59 |
integer itlo,ithi |
60 |
integer jtlo,jthi |
61 |
integer jmin,jmax |
62 |
integer imin,imax |
63 |
|
64 |
integer cbuffindex |
65 |
integer reclen, irectrue |
66 |
integer cunit2, cunit3 |
67 |
character*(80) cfile2, cfile3 |
68 |
|
69 |
real*4 cbuff ( snx*nsx*npx*nsy*npy ) |
70 |
real*4 globfldtmp2( snx,nsx,npx,nsy,npy ) |
71 |
real*4 globfldtmp3( snx,nsx,npx,nsy,npy ) |
72 |
_RL globfldxz ( snx,nsx,npx,nsy,npy,nr ) |
73 |
_RL globfld3d ( snx,nsx,npx,sny,nsy,npy,nr ) |
74 |
_RL globmskxz ( snx,nsx,npx,nsy,npy,nr,nobcs ) |
75 |
#ifdef CTRL_PACK_PRECISE |
76 |
integer il |
77 |
character*(80) weightname |
78 |
_RL weightfldxz( snx,nsx,npx,nsy,npy,nr,nobcs ) |
79 |
|
80 |
c == external == |
81 |
|
82 |
integer ilnblnk |
83 |
external ilnblnk |
84 |
#endif |
85 |
|
86 |
c == end of interface == |
87 |
|
88 |
jtlo = 1 |
89 |
jthi = nsy |
90 |
itlo = 1 |
91 |
ithi = nsx |
92 |
jmin = 1 |
93 |
jmax = sny |
94 |
imin = 1 |
95 |
imax = snx |
96 |
|
97 |
c Initialise temporary file |
98 |
do k = 1,nr |
99 |
do jp = 1,nPy |
100 |
do bj = jtlo,jthi |
101 |
do ip = 1,nPx |
102 |
do bi = itlo,ithi |
103 |
do i = imin,imax |
104 |
globfldxz (i,bi,ip,bj,jp,k) = 0. _d 0 |
105 |
globfldtmp2(i,bi,ip,bj,jp) = 0. |
106 |
globfldtmp3(i,bi,ip,bj,jp) = 0. |
107 |
do iobcs=1,nobcs |
108 |
globmskxz(i,bi,ip,bj,jp,k,iobcs) = 0. _d 0 |
109 |
enddo |
110 |
enddo |
111 |
enddo |
112 |
enddo |
113 |
enddo |
114 |
enddo |
115 |
enddo |
116 |
c Initialise temporary file |
117 |
do k = 1,nr |
118 |
do jp = 1,nPy |
119 |
do bj = jtlo,jthi |
120 |
do j = jmin,jmax |
121 |
do ip = 1,nPx |
122 |
do bi = itlo,ithi |
123 |
do i = imin,imax |
124 |
globfld3d(i,bi,ip,j,bj,jp,k) = 0. _d 0 |
125 |
enddo |
126 |
enddo |
127 |
enddo |
128 |
enddo |
129 |
enddo |
130 |
enddo |
131 |
enddo |
132 |
|
133 |
c-- Only the master thread will do I/O. |
134 |
_BEGIN_MASTER( mythid ) |
135 |
|
136 |
if ( doPackDiag ) then |
137 |
write(cfile2(1:80),'(80a)') ' ' |
138 |
write(cfile3(1:80),'(80a)') ' ' |
139 |
if ( lxxadxx ) then |
140 |
write(cfile2(1:80),'(a,I2.2,a,I4.4,a)') |
141 |
& 'diag_pack_nonout_ctrl_', |
142 |
& ivartype, '_', optimcycle, '.bin' |
143 |
write(cfile3(1:80),'(a,I2.2,a,I4.4,a)') |
144 |
& 'diag_pack_dimout_ctrl_', |
145 |
& ivartype, '_', optimcycle, '.bin' |
146 |
else |
147 |
write(cfile2(1:80),'(a,I2.2,a,I4.4,a)') |
148 |
& 'diag_pack_nonout_grad_', |
149 |
& ivartype, '_', optimcycle, '.bin' |
150 |
write(cfile3(1:80),'(a,I2.2,a,I4.4,a)') |
151 |
& 'diag_pack_dimout_grad_', |
152 |
& ivartype, '_', optimcycle, '.bin' |
153 |
endif |
154 |
|
155 |
reclen = snx*nsx*npx*nsy*npy*4 |
156 |
call mdsfindunit( cunit2, mythid ) |
157 |
open( cunit2, file=cfile2, status='unknown', |
158 |
& access='direct', recl=reclen ) |
159 |
call mdsfindunit( cunit3, mythid ) |
160 |
open( cunit3, file=cfile3, status='unknown', |
161 |
& access='direct', recl=reclen ) |
162 |
endif |
163 |
|
164 |
do iobcs = 1, nobcs |
165 |
call MDSREADFIELD_XZ_GL( |
166 |
& masktype, ctrlprec, 'RL', |
167 |
& Nr, globmskxz(1,1,1,1,1,1,iobcs), iobcs, mythid) |
168 |
#ifdef CTRL_PACK_PRECISE |
169 |
il=ilnblnk( weighttype) |
170 |
write(weightname(1:80),'(80a)') ' ' |
171 |
write(weightname(1:80),'(a)') weighttype(1:il) |
172 |
call MDSREADFIELD_XZ_GL( |
173 |
& weightname, ctrlprec, 'RL', |
174 |
& Nr, weightfldxz(1,1,1,1,1,1,iobcs), iobcs, mythid) |
175 |
#endif |
176 |
enddo |
177 |
|
178 |
if ( useSingleCPUio ) then |
179 |
C MDSREADFIELD_XZ_GL does not know about useSingleCPUio, so the faster |
180 |
C method that works for .not.useSingleCPUio cannot be used |
181 |
nrec_nl = 0 |
182 |
else |
183 |
nrec_nl = int(ncvarrecs(ivartype)/Ny) |
184 |
endif |
185 |
do irec = 1, nrec_nl |
186 |
call MDSREADFIELD_3D_GL( fname, ctrlprec, 'RL', |
187 |
& nr, globfld3d, irec, mythid) |
188 |
do j=1,sny |
189 |
iobcs= mod((irec-1)*sny+j-1,nobcs)+1 |
190 |
|
191 |
write(cunit) ncvarindex(ivartype) |
192 |
write(cunit) 1 |
193 |
write(cunit) 1 |
194 |
do k = 1,nr |
195 |
irectrue = (irec-1)*nobcs*nr + (iobcs-1)*nr + k |
196 |
cbuffindex = 0 |
197 |
do jp = 1,nPy |
198 |
do bj = jtlo,jthi |
199 |
do ip = 1,nPx |
200 |
do bi = itlo,ithi |
201 |
do i = imin,imax |
202 |
jj=mod((j-1)*nr+k-1,sny)+1 |
203 |
kk=int((j-1)*nr+K-1)/sny+1 |
204 |
if (globmskxz(i,bi,ip,bj,jp,k,iobcs) .ne. 0. ) then |
205 |
cbuffindex = cbuffindex + 1 |
206 |
cph( |
207 |
globfldtmp3(i,bi,ip,bj,jp) = |
208 |
& globfld3d(i,bi,ip,jj,bj,jp,kk) |
209 |
cph) |
210 |
#ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO |
211 |
if (lxxadxx) then |
212 |
cbuff(cbuffindex) = |
213 |
& globfld3d(i,bi,ip,jj,bj,jp,kk) * |
214 |
# ifdef CTRL_PACK_PRECISE |
215 |
& sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs)) |
216 |
# else |
217 |
& sqrt(weightfld(k,iobcs)) |
218 |
# endif |
219 |
else |
220 |
cbuff(cbuffindex) = |
221 |
& globfld3d(i,bi,ip,jj,bj,jp,kk) / |
222 |
# ifdef CTRL_PACK_PRECISE |
223 |
& sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs)) |
224 |
# else |
225 |
& sqrt(weightfld(k,iobcs)) |
226 |
# endif |
227 |
endif |
228 |
cph( |
229 |
globfldtmp2(i,bi,ip,bj,jp) = cbuff(cbuffindex) |
230 |
cph) |
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 |
247 |
if ( doPackDiag ) then |
248 |
write(cunit2,rec=irectrue) globfldtmp2 |
249 |
write(cunit3,rec=irectrue) globfldtmp3 |
250 |
endif |
251 |
c |
252 |
c -- end of k loop -- |
253 |
enddo |
254 |
c -- end of j loop -- |
255 |
enddo |
256 |
c -- end of irec loop -- |
257 |
enddo |
258 |
|
259 |
do irec = nrec_nl*ny+1, ncvarrecs(ivartype) |
260 |
c Need to solve for what iobcs would have been. |
261 |
iobcs= mod(irec-1,nobcs)+1 |
262 |
|
263 |
call MDSREADFIELD_XZ_GL( fname, ctrlprec, 'RL', |
264 |
& nr, globfldxz, irec, mythid) |
265 |
|
266 |
write(cunit) ncvarindex(ivartype) |
267 |
write(cunit) 1 |
268 |
write(cunit) 1 |
269 |
do k = 1,nr |
270 |
irectrue = (irec-1)*nobcs*nr + (iobcs-1)*nr + k |
271 |
cbuffindex = 0 |
272 |
do jp = 1,nPy |
273 |
do bj = jtlo,jthi |
274 |
do ip = 1,nPx |
275 |
do bi = itlo,ithi |
276 |
do i = imin,imax |
277 |
if (globmskxz(i,bi,ip,bj,jp,k,iobcs) .ne. 0. ) then |
278 |
cbuffindex = cbuffindex + 1 |
279 |
cph( |
280 |
globfldtmp3(i,bi,ip,bj,jp) = |
281 |
& globfldxz(i,bi,ip,bj,jp,k) |
282 |
cph) |
283 |
#ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO |
284 |
if (lxxadxx) then |
285 |
cbuff(cbuffindex) = |
286 |
& globfldxz(i,bi,ip,bj,jp,k) * |
287 |
# ifdef CTRL_PACK_PRECISE |
288 |
& sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs)) |
289 |
# else |
290 |
& sqrt(weightfld(k,iobcs)) |
291 |
# endif |
292 |
else |
293 |
cbuff(cbuffindex) = |
294 |
& globfldxz(i,bi,ip,bj,jp,k) / |
295 |
# ifdef CTRL_PACK_PRECISE |
296 |
& sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs)) |
297 |
# else |
298 |
& sqrt(weightfld(k,iobcs)) |
299 |
# endif |
300 |
endif |
301 |
cph( |
302 |
globfldtmp2(i,bi,ip,bj,jp) = cbuff(cbuffindex) |
303 |
cph) |
304 |
#else /* ALLOW_NONDIMENSIONAL_CONTROL_IO undef */ |
305 |
cbuff(cbuffindex) = globfldxz(i,bi,ip,bj,jp,k) |
306 |
#endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */ |
307 |
endif |
308 |
enddo |
309 |
enddo |
310 |
enddo |
311 |
enddo |
312 |
enddo |
313 |
c --> check cbuffindex. |
314 |
if ( cbuffindex .gt. 0) then |
315 |
write(cunit) cbuffindex |
316 |
write(cunit) k |
317 |
write(cunit) (cbuff(ii), ii=1,cbuffindex) |
318 |
endif |
319 |
c |
320 |
if ( doPackDiag ) then |
321 |
write(cunit2,rec=irectrue) globfldtmp2 |
322 |
write(cunit3,rec=irectrue) globfldtmp3 |
323 |
endif |
324 |
c |
325 |
c -- end of k loop -- |
326 |
enddo |
327 |
c -- end of irec loop -- |
328 |
enddo |
329 |
|
330 |
_END_MASTER( mythid ) |
331 |
|
332 |
#endif |
333 |
|
334 |
return |
335 |
end |
336 |
|
337 |
|
338 |
|
339 |
|
340 |
|