1 |
C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_set_unpack_yz.F,v 1.18 2012/08/10 19:38:58 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "CTRL_OPTIONS.h" |
5 |
|
6 |
subroutine ctrl_set_unpack_yz( |
7 |
& cunit, ivartype, fname, masktype, weighttype, |
8 |
& weightfld, nwetglobal, mythid) |
9 |
|
10 |
c ================================================================== |
11 |
c SUBROUTINE ctrl_set_unpack_yz |
12 |
c ================================================================== |
13 |
c |
14 |
c o Unpack the control vector such that land points are filled in. |
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 changes from Armin 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 "CTRL_OBCS.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 |
character*( 80) weighttype |
46 |
_RL weightfld( nr,nobcs ) |
47 |
integer nwetglobal(nr,nobcs) |
48 |
integer mythid |
49 |
|
50 |
#ifndef EXCLUDE_CTRL_PACK |
51 |
c == local variables == |
52 |
|
53 |
logical lxxadxx |
54 |
|
55 |
integer bi,bj |
56 |
integer ip,jp |
57 |
integer i,j,k |
58 |
integer ii,jj,kk |
59 |
integer irec,iobcs,nrec_nl |
60 |
integer itlo,ithi |
61 |
integer jtlo,jthi |
62 |
integer jmin,jmax |
63 |
integer imin,imax |
64 |
|
65 |
integer cbuffindex |
66 |
|
67 |
real*4 cbuff ( nsx*npx*sny*nsy*npy ) |
68 |
real*4 globfldtmp2( nsx,npx,sny,nsy,npy ) |
69 |
real*4 globfldtmp3( nsx,npx,sny,nsy,npy ) |
70 |
_RL globfldyz( nsx,npx,sny,nsy,npy,nr ) |
71 |
_RL globfld3d( snx,nsx,npx,sny,nsy,npy,nr ) |
72 |
_RL globmskyz( nsx,npx,sny,nsy,npy,nr,nobcs ) |
73 |
|
74 |
integer reclen, irectrue |
75 |
integer cunit2, cunit3 |
76 |
character*(80) cfile2, cfile3 |
77 |
|
78 |
#ifdef CTRL_UNPACK_PRECISE |
79 |
integer il |
80 |
character*(80) weightname |
81 |
_RL weightfldyz( nsx,npx,sny,nsy,npy,nr,nobcs ) |
82 |
|
83 |
c == external == |
84 |
|
85 |
integer ilnblnk |
86 |
external ilnblnk |
87 |
#endif |
88 |
|
89 |
cc == end of interface == |
90 |
|
91 |
jtlo = 1 |
92 |
jthi = nsy |
93 |
itlo = 1 |
94 |
ithi = nsx |
95 |
jmin = 1 |
96 |
jmax = sny |
97 |
imin = 1 |
98 |
imax = snx |
99 |
|
100 |
lxxadxx = .TRUE. |
101 |
|
102 |
c Initialise temporary file |
103 |
do k = 1,nr |
104 |
do jp = 1,nPy |
105 |
do bj = jtlo,jthi |
106 |
do j = jmin,jmax |
107 |
do ip = 1,nPx |
108 |
do bi = itlo,ithi |
109 |
globfldyz (bi,ip,j,bj,jp,k) = 0. _d 0 |
110 |
globfldtmp2(bi,ip,j,bj,jp) = 0. |
111 |
globfldtmp3(bi,ip,j,bj,jp) = 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 |
if ( doPackDiag ) then |
142 |
write(cfile2(1:80),'(80a)') ' ' |
143 |
write(cfile3(1:80),'(80a)') ' ' |
144 |
if ( lxxadxx ) then |
145 |
write(cfile2(1:80),'(a,I2.2,a,I4.4,a)') |
146 |
& 'diag_unpack_nondim_ctrl_', |
147 |
& ivartype, '_', optimcycle, '.bin' |
148 |
write(cfile3(1:80),'(a,I2.2,a,I4.4,a)') |
149 |
& 'diag_unpack_dimens_ctrl_', |
150 |
& ivartype, '_', optimcycle, '.bin' |
151 |
else |
152 |
write(cfile2(1:80),'(a,I2.2,a,I4.4,a)') |
153 |
& 'diag_unpack_nondim_grad_', |
154 |
& ivartype, '_', optimcycle, '.bin' |
155 |
write(cfile3(1:80),'(a,I2.2,a,I4.4,a)') |
156 |
& 'diag_unpack_dimens_grad_', |
157 |
& ivartype, '_', optimcycle, '.bin' |
158 |
endif |
159 |
|
160 |
reclen = nsx*npx*sny*nsy*npy*4 |
161 |
call mdsfindunit( cunit2, mythid ) |
162 |
open( cunit2, file=cfile2, status='unknown', |
163 |
& access='direct', recl=reclen ) |
164 |
call mdsfindunit( cunit3, mythid ) |
165 |
open( cunit3, file=cfile3, status='unknown', |
166 |
& access='direct', recl=reclen ) |
167 |
endif |
168 |
|
169 |
do iobcs=1,nobcs |
170 |
call MDSREADFIELD_YZ_GL( |
171 |
& masktype, ctrlprec, 'RL', |
172 |
& Nr, globmskyz(1,1,1,1,1,1,iobcs), iobcs, mythid) |
173 |
#ifdef CTRL_UNPACK_PRECISE |
174 |
il=ilnblnk( weighttype) |
175 |
write(weightname(1:80),'(80a)') ' ' |
176 |
write(weightname(1:80),'(a)') weighttype(1:il) |
177 |
call MDSREADFIELD_YZ_GL( |
178 |
& weightname, ctrlprec, 'RL', |
179 |
& Nr, weightfldyz(1,1,1,1,1,1,iobcs), iobcs, mythid) |
180 |
#endif /* CTRL_UNPACK_PRECISE */ |
181 |
enddo |
182 |
|
183 |
if ( useSingleCPUio ) then |
184 |
C MDSWRITEFIELD_YZ_GL does not know about useSingleCPUio, so the faster |
185 |
C method that works for .not.useSingleCPUio cannot be used |
186 |
nrec_nl = 0 |
187 |
else |
188 |
nrec_nl = int(ncvarrecs(ivartype)/Nx) |
189 |
endif |
190 |
do irec = 1, nrec_nl |
191 |
c And now back-calculate what iobcs should be. |
192 |
do i=1,snx |
193 |
iobcs= mod((irec-1)*snx+i-1,nobcs)+1 |
194 |
|
195 |
read(cunit) filencvarindex(ivartype) |
196 |
if (filencvarindex(ivartype) .NE. ncvarindex(ivartype)) |
197 |
& then |
198 |
print *, 'ctrl_set_unpack_yz:WARNING: wrong ncvarindex ', |
199 |
& filencvarindex(ivartype), ncvarindex(ivartype) |
200 |
STOP 'in S/R ctrl_set_unpack_yz' |
201 |
endif |
202 |
read(cunit) filej |
203 |
read(cunit) filei |
204 |
do k = 1, Nr |
205 |
irectrue = (irec-1)*nobcs*nr + (iobcs-1)*nr + k |
206 |
cbuffindex = nwetglobal(k,iobcs) |
207 |
if ( cbuffindex .gt. 0 ) then |
208 |
read(cunit) filencbuffindex |
209 |
if (filencbuffindex .NE. cbuffindex) then |
210 |
print *, 'WARNING: wrong cbuffindex ', |
211 |
& filencbuffindex, cbuffindex |
212 |
STOP 'in S/R ctrl_set_unpack_yz' |
213 |
endif |
214 |
read(cunit) filek |
215 |
if (filek .NE. k) then |
216 |
print *, 'WARNING: wrong k ', |
217 |
& filek, k |
218 |
STOP 'in S/R ctrl_set_unpack_yz' |
219 |
endif |
220 |
read(cunit) (cbuff(ii), ii=1,cbuffindex) |
221 |
endif |
222 |
cbuffindex = 0 |
223 |
do jp = 1,nPy |
224 |
do bj = jtlo,jthi |
225 |
do j = jmin,jmax |
226 |
do ip = 1,nPx |
227 |
do bi = itlo,ithi |
228 |
ii=mod((i-1)*nr*sny+(k-1)*sny+j-1,snx)+1 |
229 |
jj=mod(((i-1)*nr*sny+(k-1)*sny+j-1)/snx,sny)+1 |
230 |
kk=int((i-1)*nr*sny+(k-1)*sny+j-1)/(snx*sny)+1 |
231 |
if ( globmskyz(bi,ip,j,bj,jp,k,iobcs) .ne. 0. ) then |
232 |
cbuffindex = cbuffindex + 1 |
233 |
globfld3d(ii,bi,ip,jj,bj,jp,kk) = |
234 |
& cbuff(cbuffindex) |
235 |
cph( |
236 |
globfldtmp2(bi,ip,jj,bj,jp) = |
237 |
& cbuff(cbuffindex) |
238 |
cph) |
239 |
#ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO |
240 |
globfld3d(ii,bi,ip,jj,bj,jp,kk) = |
241 |
& globfld3d(ii,bi,ip,jj,bj,jp,kk)/ |
242 |
# ifdef CTRL_UNPACK_PRECISE |
243 |
& sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs)) |
244 |
# else |
245 |
& sqrt(weightfld(k,iobcs)) |
246 |
# endif |
247 |
#endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */ |
248 |
else |
249 |
globfld3d(ii,bi,ip,jj,bj,jp,kk) = 0. _d 0 |
250 |
endif |
251 |
cph( |
252 |
globfldtmp3(bi,ip,jj,bj,jp) = |
253 |
& globfld3d(ii,bi,ip,jj,bj,jp,kk) |
254 |
cph) |
255 |
enddo |
256 |
enddo |
257 |
enddo |
258 |
enddo |
259 |
enddo |
260 |
c |
261 |
if ( doPackDiag ) then |
262 |
write(cunit2,rec=irectrue) globfldtmp2 |
263 |
write(cunit3,rec=irectrue) globfldtmp3 |
264 |
endif |
265 |
c |
266 |
c -- end of k loop -- |
267 |
enddo |
268 |
c -- end of i loop -- |
269 |
enddo |
270 |
|
271 |
call MDSWRITEFIELD_3d_GL( fname, ctrlprec, 'RL', |
272 |
& Nr, globfld3d, irec, |
273 |
& optimcycle, mythid) |
274 |
|
275 |
c -- end of irec loop -- |
276 |
enddo |
277 |
|
278 |
do irec = nrec_nl*nx+1,ncvarrecs(ivartype) |
279 |
c And now back-calculate what iobcs should be. |
280 |
iobcs= mod(irec-1,nobcs)+1 |
281 |
|
282 |
read(cunit) filencvarindex(ivartype) |
283 |
if (filencvarindex(ivartype) .NE. ncvarindex(ivartype)) |
284 |
& then |
285 |
print *, 'ctrl_set_unpack_yz:WARNING: wrong ncvarindex ', |
286 |
& filencvarindex(ivartype), ncvarindex(ivartype) |
287 |
STOP 'in S/R ctrl_set_unpack_yz' |
288 |
endif |
289 |
read(cunit) filej |
290 |
read(cunit) filei |
291 |
do k = 1, Nr |
292 |
irectrue = (irec-1)*nobcs*nr + (iobcs-1)*nr + k |
293 |
cbuffindex = nwetglobal(k,iobcs) |
294 |
if ( cbuffindex .gt. 0 ) then |
295 |
read(cunit) filencbuffindex |
296 |
if (filencbuffindex .NE. cbuffindex) then |
297 |
print *, 'WARNING: wrong cbuffindex ', |
298 |
& filencbuffindex, cbuffindex |
299 |
STOP 'in S/R ctrl_set_unpack_yz' |
300 |
endif |
301 |
read(cunit) filek |
302 |
if (filek .NE. k) then |
303 |
print *, 'WARNING: wrong k ', |
304 |
& filek, k |
305 |
STOP 'in S/R ctrl_set_unpack_yz' |
306 |
endif |
307 |
read(cunit) (cbuff(ii), ii=1,cbuffindex) |
308 |
endif |
309 |
cbuffindex = 0 |
310 |
do jp = 1,nPy |
311 |
do bj = jtlo,jthi |
312 |
do j = jmin,jmax |
313 |
do ip = 1,nPx |
314 |
do bi = itlo,ithi |
315 |
if ( globmskyz(bi,ip,j,bj,jp,k,iobcs) .ne. 0. ) then |
316 |
cbuffindex = cbuffindex + 1 |
317 |
globfldyz(bi,ip,j,bj,jp,k) = cbuff(cbuffindex) |
318 |
cph( |
319 |
globfldtmp2(bi,ip,j,bj,jp) = cbuff(cbuffindex) |
320 |
cph) |
321 |
#ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO |
322 |
globfldyz(bi,ip,j,bj,jp,k) = |
323 |
& globfldyz(bi,ip,j,bj,jp,k)/ |
324 |
# ifdef CTRL_UNPACK_PRECISE |
325 |
& sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs)) |
326 |
# else |
327 |
& sqrt(weightfld(k,iobcs)) |
328 |
# endif |
329 |
#endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */ |
330 |
else |
331 |
globfldyz(bi,ip,j,bj,jp,k) = 0. _d 0 |
332 |
endif |
333 |
cph( |
334 |
globfldtmp3(bi,ip,j,bj,jp) = |
335 |
& globfldyz(bi,ip,j,bj,jp,k) |
336 |
cph) |
337 |
enddo |
338 |
enddo |
339 |
enddo |
340 |
enddo |
341 |
enddo |
342 |
c |
343 |
c -- end of k loop |
344 |
enddo |
345 |
|
346 |
call MDSWRITEFIELD_YZ_GL( fname, ctrlprec, 'RL', |
347 |
& Nr, globfldyz, irec, |
348 |
& optimcycle, mythid) |
349 |
|
350 |
c -- end of irec loop -- |
351 |
enddo |
352 |
|
353 |
_END_MASTER( mythid ) |
354 |
|
355 |
#endif |
356 |
|
357 |
return |
358 |
end |
359 |
|