/[MITgcm]/MITgcm/pkg/ctrl/ctrl_set_unpack_xz.F
ViewVC logotype

Contents of /MITgcm/pkg/ctrl/ctrl_set_unpack_xz.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.13 - (show annotations) (download)
Wed Jan 23 22:38:43 2008 UTC (16 years, 4 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint60, checkpoint61, checkpoint59q, checkpoint59p, checkpoint59r, checkpoint59o, checkpoint59n, checkpoint61f, checkpoint61n, checkpoint61q, checkpoint61e, checkpoint61g, checkpoint61d, checkpoint61b, checkpoint61c, checkpoint61a, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61j, checkpoint61k, checkpoint61h, checkpoint61i, checkpoint61v, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p
Changes since 1.12: +71 -3 lines
Mehr fuer die Luetten
(this time unpack).

1 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_set_unpack_xz.F,v 1.12 2007/10/09 00:00:01 jmc Exp $
2 C $Name: $
3
4 #include "CTRL_CPPOPTIONS.h"
5
6 subroutine ctrl_set_unpack_xz(
7 & cunit, ivartype, fname, masktype, weighttype,
8 & weightfld, nwetglobal, mythid)
9
10 c ==================================================================
11 c SUBROUTINE ctrl_set_unpack_xz
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 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 "optim.h"
37
38 c == routine arguments ==
39
40 integer cunit
41 integer ivartype
42 character*( 80) fname
43 character* (9) masktype
44 character*( 80) weighttype
45 _RL weightfld( nr,nobcs )
46 integer nwetglobal(nr,nobcs)
47 integer mythid
48
49 #ifndef EXCLUDE_CTRL_PACK
50 c == local variables ==
51
52 logical lxxadxx
53
54 integer bi,bj
55 integer ip,jp
56 integer i,j,k
57 integer ii,jj,kk
58 integer il
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 ( snx*nsx*npx*nsy*npy )
68 real*4 globfldtmp2( snx,nsx,npx,nsy,npy )
69 real*4 globfldtmp3( snx,nsx,npx,nsy,npy )
70 _RL globfldxz( snx,nsx,npx,nsy,npy,nr )
71 _RL globfld3d( snx,nsx,npx,sny,nsy,npy,nr )
72 _RL globmskxz( snx,nsx,npx,nsy,npy,nr,nobcs )
73 #ifdef CTRL_UNPACK_PRECISE
74 _RL weightfldxz( snx,nsx,npx,nsy,npy,nr,nobcs )
75 #endif
76
77 cgg(
78 integer reclen, irectrue
79 integer cunit2, cunit3
80 integer igg
81 _RL gg
82 character*(80) weightname
83 character*(80) cfile2, cfile3
84 cgg)
85
86 c == external ==
87
88 integer ilnblnk
89 external ilnblnk
90
91 cc == end of interface ==
92
93 jtlo = 1
94 jthi = nsy
95 itlo = 1
96 ithi = nsx
97 jmin = 1
98 jmax = sny
99 imin = 1
100 imax = snx
101
102 lxxadxx = .TRUE.
103
104 c Initialise temporary file
105 do k = 1,nr
106 do jp = 1,nPy
107 do bj = jtlo,jthi
108 do ip = 1,nPx
109 do bi = itlo,ithi
110 do i = imin,imax
111 globfldxz (i,bi,ip,bj,jp,k) = 0. _d 0
112 globfldtmp2(i,bi,ip,bj,jp) = 0.
113 globfldtmp3(i,bi,ip,bj,jp) = 0.
114 do iobcs=1,nobcs
115 globmskxz(i,bi,ip,bj,jp,k,iobcs) = 0. _d 0
116 enddo
117 enddo
118 enddo
119 enddo
120 enddo
121 enddo
122 enddo
123 c Initialise temporary file
124 do k = 1,nr
125 do jp = 1,nPy
126 do bj = jtlo,jthi
127 do j = jmin,jmax
128 do ip = 1,nPx
129 do bi = itlo,ithi
130 do i = imin,imax
131 globfld3d(i,bi,ip,j,bj,jp,k) = 0. _d 0
132 enddo
133 enddo
134 enddo
135 enddo
136 enddo
137 enddo
138 enddo
139
140 c-- Only the master thread will do I/O.
141 _BEGIN_MASTER( mythid )
142
143 if ( doPackDiag ) then
144 write(cfile2(1:80),'(80a)') ' '
145 write(cfile3(1:80),'(80a)') ' '
146 if ( lxxadxx ) then
147 write(cfile2(1:80),'(a,I2.2,a,I4.4,a)')
148 & 'diag_unpack_nondim_ctrl_',
149 & ivartype, '_', optimcycle, '.bin'
150 write(cfile3(1:80),'(a,I2.2,a,I4.4,a)')
151 & 'diag_unpack_dimens_ctrl_',
152 & ivartype, '_', optimcycle, '.bin'
153 else
154 write(cfile2(1:80),'(a,I2.2,a,I4.4,a)')
155 & 'diag_unpack_nondim_grad_',
156 & ivartype, '_', optimcycle, '.bin'
157 write(cfile3(1:80),'(a,I2.2,a,I4.4,a)')
158 & 'diag_unpack_dimens_grad_',
159 & ivartype, '_', optimcycle, '.bin'
160 endif
161
162 reclen = snx*nsx*npx*nsy*npy*4
163 call mdsfindunit( cunit2, mythid )
164 open( cunit2, file=cfile2, status='unknown',
165 & access='direct', recl=reclen )
166 call mdsfindunit( cunit3, mythid )
167 open( cunit3, file=cfile3, status='unknown',
168 & access='direct', recl=reclen )
169 endif
170
171 do iobcs=1,nobcs
172 call MDSREADFIELD_XZ_GL(
173 & masktype, ctrlprec, 'RL',
174 & Nr, globmskxz(1,1,1,1,1,1,iobcs), iobcs,mythid)
175 #ifdef CTRL_UNPACK_PRECISE
176 il=ilnblnk( weighttype)
177 write(weightname(1:80),'(80a)') ' '
178 write(weightname(1:80),'(a)') weighttype(1:il)
179 call MDSREADFIELD_XZ_GL(
180 & weightname, ctrlprec, 'RL',
181 & Nr, weightfldxz(1,1,1,1,1,1,iobcs), iobcs, mythid)
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 ip = 1,nPx
189 do bi = itlo,ithi
190 do i = imin,imax
191 cph weightfldxz(i,bi,ip,bj,jp,k,iobcs) = wbaro
192 enddo
193 enddo
194 enddo
195 enddo
196 enddo
197 endif
198 #endif /* CTRL_UNPACK_PRECISE */
199 enddo
200
201 nrec_nl=int(ncvarrecs(ivartype)/sny)
202 do irec = 1, nrec_nl
203 cgg do iobcs = 1, nobcs
204 cgg And now back-calculate what iobcs should be.
205 do j=1,sny
206 iobcs= mod((irec-1)*sny+j-1,nobcs)+1
207
208 read(cunit) filencvarindex(ivartype)
209 if (filencvarindex(ivartype) .NE. ncvarindex(ivartype))
210 & then
211 print *, 'ctrl-set_unpack:xz:WARNING: wrong ncvarindex ',
212 & filencvarindex(ivartype), ncvarindex(ivartype)
213 STOP 'in S/R ctrl_unpack'
214 endif
215 read(cunit) filej
216 read(cunit) filei
217 do k = 1, Nr
218 irectrue = (irec-1)*nobcs*nr + (iobcs-1)*nr + k
219 cbuffindex = nwetglobal(k,iobcs)
220 if ( cbuffindex .gt. 0 ) then
221 read(cunit) filencbuffindex
222 if (filencbuffindex .NE. cbuffindex) then
223 print *, 'WARNING: wrong cbuffindex ',
224 & filencbuffindex, cbuffindex
225 STOP 'in S/R ctrl_unpack'
226 endif
227 read(cunit) filek
228 if (filek .NE. k) then
229 print *, 'WARNING: wrong k ',
230 & filek, k
231 STOP 'in S/R ctrl_unpack'
232 endif
233 read(cunit) (cbuff(ii), ii=1,cbuffindex)
234 endif
235
236 cbuffindex = 0
237 jj=mod((j-1)*nr+k-1,sny)+1
238 kk=int((j-1)*nr+k-1)/sny+1
239 do jp = 1,nPy
240 do bj = jtlo,jthi
241 do ip = 1,nPx
242 do bi = itlo,ithi
243 do i = imin,imax
244 if ( globmskxz(i,bi,ip,bj,jp,k,iobcs) .ne. 0. ) then
245 cbuffindex = cbuffindex + 1
246 globfld3d(i,bi,ip,jj,bj,jp,kk) =
247 & cbuff(cbuffindex)
248 cph(
249 globfldtmp2(i,bi,ip,bj,jp) =
250 & cbuff(cbuffindex)
251 cph)
252 #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
253 globfld3d(i,bi,ip,jj,bj,jp,kk) =
254 & globfld3d(i,bi,ip,jj,bj,jp,kk)/
255 # ifdef CTRL_UNPACK_PRECISE
256 & sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
257 # else
258 & sqrt(weightfld(k,iobcs))
259 # endif
260 #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
261 else
262 globfld3d(i,bi,ip,jj,bj,jp,kk) = 0. _d 0
263 endif
264 cph(
265 globfldtmp3(i,bi,ip,bj,jp) =
266 & globfld3d(i,bi,ip,jj,bj,jp,kk)
267 cph)
268 enddo
269 enddo
270 enddo
271 enddo
272 enddo
273 c
274 if ( doPackDiag ) then
275 write(cunit2,rec=irectrue) globfldtmp2
276 write(cunit3,rec=irectrue) globfldtmp3
277 endif
278 c
279 c -- end of k loop --
280 enddo
281 c -- end of j loop --
282 enddo
283
284 call MDSWRITEFIELD_3D_GL( fname, ctrlprec, 'RL',
285 & Nr, globfld3d, irec,
286 & optimcycle, mythid)
287
288 c -- end of iobcs loop -- This loop removed. 3-28-02.
289 cgg enddo
290 c -- end of irec loop --
291 enddo
292
293 do irec = nrec_nl*sny+1, ncvarrecs(ivartype)
294 cgg do iobcs = 1, nobcs
295 cgg And now back-calculate what iobcs should be.
296 iobcs= mod(irec-1,nobcs)+1
297
298 read(cunit) filencvarindex(ivartype)
299 if (filencvarindex(ivartype) .NE. ncvarindex(ivartype))
300 & then
301 print *, 'ctrl-set_unpack:xz:WARNING: wrong ncvarindex ',
302 & filencvarindex(ivartype), ncvarindex(ivartype)
303 STOP 'in S/R ctrl_unpack'
304 endif
305 read(cunit) filej
306 read(cunit) filei
307 do k = 1, Nr
308 irectrue = (irec-1)*nobcs*nr + (iobcs-1)*nr + k
309 cbuffindex = nwetglobal(k,iobcs)
310 if ( cbuffindex .gt. 0 ) then
311 read(cunit) filencbuffindex
312 if (filencbuffindex .NE. cbuffindex) then
313 print *, 'WARNING: wrong cbuffindex ',
314 & filencbuffindex, cbuffindex
315 STOP 'in S/R ctrl_unpack'
316 endif
317 read(cunit) filek
318 if (filek .NE. k) then
319 print *, 'WARNING: wrong k ',
320 & filek, k
321 STOP 'in S/R ctrl_unpack'
322 endif
323 read(cunit) (cbuff(ii), ii=1,cbuffindex)
324 endif
325
326 cbuffindex = 0
327 do jp = 1,nPy
328 do bj = jtlo,jthi
329 do ip = 1,nPx
330 do bi = itlo,ithi
331 do i = imin,imax
332 if ( globmskxz(i,bi,ip,bj,jp,k,iobcs) .ne. 0. ) then
333 cbuffindex = cbuffindex + 1
334 globfldxz(i,bi,ip,bj,jp,k) = cbuff(cbuffindex)
335 cph(
336 globfldtmp2(i,bi,ip,bj,jp) = cbuff(cbuffindex)
337 cph)
338 #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
339 globfldxz(i,bi,ip,bj,jp,k) =
340 & globfldxz(i,bi,ip,bj,jp,k)/
341 # ifdef CTRL_UNPACK_PRECISE
342 & sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
343 # else
344 & sqrt(weightfld(k,iobcs))
345 # endif
346 #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
347 else
348 globfldxz(i,bi,ip,bj,jp,k) = 0. _d 0
349 endif
350 cph(
351 globfldtmp3(i,bi,ip,bj,jp) =
352 & globfldxz(i,bi,ip,bj,jp,k)
353 cph)
354 enddo
355 enddo
356 enddo
357 enddo
358 enddo
359 c
360 if ( doPackDiag ) then
361 write(cunit2,rec=irectrue) globfldtmp2
362 write(cunit3,rec=irectrue) globfldtmp3
363 endif
364 c
365 c -- end of k loop --
366 enddo
367
368 call MDSWRITEFIELD_XZ_GL( fname, ctrlprec, 'RL',
369 & Nr, globfldxz, irec,
370 & optimcycle, mythid)
371
372 c -- end of iobcs loop -- This loop removed. 3-28-02.
373 cgg enddo
374 c -- end of irec loop --
375 enddo
376
377 _END_MASTER( mythid )
378
379 #endif
380
381 return
382 end
383
384
385
386
387

  ViewVC Help
Powered by ViewVC 1.1.22