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

Annotation of /MITgcm/pkg/ctrl/ctrl_set_pack_xz.F

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


Revision 1.14 - (hide annotations) (download)
Mon Mar 22 02:16:43 2010 UTC (14 years, 2 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62w, checkpoint62z, checkpoint62y, checkpoint62x, checkpoint63
Changes since 1.13: +2 -2 lines
finish removing unbalanced quote (single or double) in commented line

1 jmc 1.14 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_set_pack_xz.F,v 1.13 2009/09/30 16:03:20 mlosch Exp $
2 jmc 1.11 C $Name: $
3 heimbach 1.2
4     #include "CTRL_CPPOPTIONS.h"
5    
6     subroutine ctrl_set_pack_xz(
7 heimbach 1.3 & cunit, ivartype, fname, masktype,weighttype,
8 heimbach 1.2 & 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 heimbach 1.3 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 jmc 1.14 c merged changes from Armin to replace write of
22 heimbach 1.3 c nr * globfld2d by 1 * globfld3d
23     c (ad hoc fix to speed up global I/O)
24     c
25 heimbach 1.2 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 "optim.h"
38    
39     c == routine arguments ==
40    
41     integer cunit
42     integer ivartype
43     character*( 80) fname
44     character*( 9) masktype
45 heimbach 1.3 character*( 80) weighttype
46 heimbach 1.2 _RL weightfld( nr,nobcs )
47     logical lxxadxx
48     integer mythid
49    
50 heimbach 1.10 #ifndef EXCLUDE_CTRL_PACK
51 heimbach 1.2 c == local variables ==
52    
53     integer bi,bj
54     integer ip,jp
55     integer i,j,k
56 heimbach 1.4 integer ii,jj,kk
57 heimbach 1.2 integer il
58 heimbach 1.3 integer irec,iobcs,nrec_nl
59 heimbach 1.2 integer itlo,ithi
60     integer jtlo,jthi
61     integer jmin,jmax
62     integer imin,imax
63    
64     integer cbuffindex
65     cgg(
66     integer igg
67 heimbach 1.12 integer reclen, irectrue
68     integer cunit2, cunit3
69 heimbach 1.2 _RL gg
70 heimbach 1.3 character*(80) weightname
71 heimbach 1.12 character*(80) cfile2, cfile3
72    
73 heimbach 1.2 cgg)
74    
75 heimbach 1.12 real*4 cbuff ( snx*nsx*npx*nsy*npy )
76     real*4 globfldtmp2( snx,nsx,npx,nsy,npy )
77     real*4 globfldtmp3( snx,nsx,npx,nsy,npy )
78 heimbach 1.2 _RL globfldxz ( snx,nsx,npx,nsy,npy,nr )
79 heimbach 1.3 _RL globfld3d ( snx,nsx,npx,sny,nsy,npy,nr )
80     _RL globmskxz ( snx,nsx,npx,nsy,npy,nr,nobcs )
81     #ifdef CTRL_PACK_PRECISE
82     _RL weightfldxz( snx,nsx,npx,nsy,npy,nr,nobcs )
83     #endif
84 heimbach 1.2
85     c == external ==
86    
87     integer ilnblnk
88     external ilnblnk
89    
90     c == end of interface ==
91    
92     jtlo = 1
93     jthi = nsy
94     itlo = 1
95     ithi = nsx
96     jmin = 1
97     jmax = sny
98     imin = 1
99     imax = snx
100    
101     c Initialise temporary file
102     do k = 1,nr
103     do jp = 1,nPy
104     do bj = jtlo,jthi
105     do ip = 1,nPx
106     do bi = itlo,ithi
107     do i = imin,imax
108 heimbach 1.12 globfldxz (i,bi,ip,bj,jp,k) = 0. _d 0
109     globfldtmp2(i,bi,ip,bj,jp) = 0.
110     globfldtmp3(i,bi,ip,bj,jp) = 0.
111 heimbach 1.3 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 heimbach 1.2 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 heimbach 1.12 if ( doPackDiag ) then
141     write(cfile2(1:80),'(80a)') ' '
142     write(cfile3(1:80),'(80a)') ' '
143     if ( lxxadxx ) then
144     write(cfile2(1:80),'(a,I2.2,a,I4.4,a)')
145     & 'diag_pack_nonout_ctrl_',
146     & ivartype, '_', optimcycle, '.bin'
147     write(cfile3(1:80),'(a,I2.2,a,I4.4,a)')
148     & 'diag_pack_dimout_ctrl_',
149     & ivartype, '_', optimcycle, '.bin'
150     else
151     write(cfile2(1:80),'(a,I2.2,a,I4.4,a)')
152     & 'diag_pack_nonout_grad_',
153     & ivartype, '_', optimcycle, '.bin'
154     write(cfile3(1:80),'(a,I2.2,a,I4.4,a)')
155     & 'diag_pack_dimout_grad_',
156     & ivartype, '_', optimcycle, '.bin'
157     endif
158    
159     reclen = snx*nsx*npx*nsy*npy*4
160     call mdsfindunit( cunit2, mythid )
161     open( cunit2, file=cfile2, status='unknown',
162     & access='direct', recl=reclen )
163     call mdsfindunit( cunit3, mythid )
164     open( cunit3, file=cfile3, status='unknown',
165     & access='direct', recl=reclen )
166     endif
167    
168 heimbach 1.3 do iobcs = 1, nobcs
169     call MDSREADFIELD_XZ_GL(
170     & masktype, ctrlprec, 'RL',
171     & Nr, globmskxz(1,1,1,1,1,1,iobcs), iobcs, mythid)
172     #ifdef CTRL_PACK_PRECISE
173     il=ilnblnk( weighttype)
174     write(weightname(1:80),'(80a)') ' '
175     write(weightname(1:80),'(a)') weighttype(1:il)
176     call MDSREADFIELD_XZ_GL(
177     & weightname, ctrlprec, 'RL',
178     & Nr, weightfldxz(1,1,1,1,1,1,iobcs), iobcs, mythid)
179     CGG One special exception: barotropic velocity should be nondimensionalized
180     cgg differently. Probably introduce new variable.
181     if (iobcs .eq. 3 .or. iobcs .eq. 4) then
182     k = 1
183     do jp = 1,nPy
184     do bj = jtlo,jthi
185     do ip = 1,nPx
186     do bi = itlo,ithi
187     do i = imin,imax
188 heimbach 1.9 cph weightfldxz(i,bi,ip,bj,jp,k,iobcs) = wbaro
189 heimbach 1.3 enddo
190     enddo
191     enddo
192     enddo
193     enddo
194     endif
195     #endif
196     enddo
197    
198 mlosch 1.13 if ( useSingleCPUio ) then
199     C MDSREADFIELD_XZ_GL does not know about useSingleCPUio, so the faster
200     C method that works for .not.useSingleCPUio cannot be used
201     nrec_nl = 0
202     else
203     nrec_nl = int(ncvarrecs(ivartype)/Ny)
204     endif
205 heimbach 1.3 do irec = 1, nrec_nl
206     call MDSREADFIELD_3D_GL( fname, ctrlprec, 'RL',
207     & nr, globfld3d, irec, mythid)
208     do j=1,sny
209     iobcs= mod((irec-1)*sny+j-1,nobcs)+1
210    
211     CGG One special exception: barotropic velocity should be nondimensionalized
212     cgg differently. Probably introduce new variable.
213     if (iobcs .eq. 3 .or. iobcs .eq. 4) then
214     k = 1
215     do jp = 1,nPy
216     do bj = jtlo,jthi
217     do ip = 1,nPx
218     do bi = itlo,ithi
219     do i = imin,imax
220     #ifdef NO_CONTROL_BAROTROPIC_VELOCITY
221     if (.not. lxxadxx) then
222     cgg Get rid of any sensitivity to barotropic velocity.
223     globfld3d(i,bi,ip,j,bj,jp,k) = 0.
224     endif
225     #endif
226     enddo
227     enddo
228     enddo
229     enddo
230     enddo
231     endif
232    
233     write(cunit) ncvarindex(ivartype)
234     write(cunit) 1
235     write(cunit) 1
236     do k = 1,nr
237 heimbach 1.12 irectrue = (irec-1)*nobcs*nr + (iobcs-1)*nr + k
238 heimbach 1.3 cbuffindex = 0
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 heimbach 1.7 jj=mod((j-1)*nr+k-1,sny)+1
245     kk=int((j-1)*nr+K-1)/sny+1
246 heimbach 1.3 if (globmskxz(i,bi,ip,bj,jp,k,iobcs) .ne. 0. ) then
247     cbuffindex = cbuffindex + 1
248 heimbach 1.12 cph(
249     globfldtmp3(i,bi,ip,bj,jp) =
250     & globfld3d(i,bi,ip,jj,bj,jp,kk)
251     cph)
252 heimbach 1.3 #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
253     if (lxxadxx) then
254 jmc 1.11 cbuff(cbuffindex) =
255 heimbach 1.4 & globfld3d(i,bi,ip,jj,bj,jp,kk) *
256 heimbach 1.3 # ifdef CTRL_PACK_PRECISE
257     & sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
258     # else
259     & sqrt(weightfld(k,iobcs))
260     # endif
261     else
262 jmc 1.11 cbuff(cbuffindex) =
263 heimbach 1.4 & globfld3d(i,bi,ip,jj,bj,jp,kk) /
264 heimbach 1.3 # ifdef CTRL_PACK_PRECISE
265     & sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
266     # else
267     & sqrt(weightfld(k,iobcs))
268     # endif
269     endif
270 heimbach 1.12 cph(
271     globfldtmp2(i,bi,ip,bj,jp) = cbuff(cbuffindex)
272     cph)
273 heimbach 1.3 #else /* ALLOW_NONDIMENSIONAL_CONTROL_IO undef */
274 heimbach 1.4 cbuff(cbuffindex) = globfld3d(i,bi,ip,jj,bj,jp,kk)
275 heimbach 1.3 #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
276     endif
277     enddo
278     enddo
279     enddo
280     enddo
281     enddo
282     c --> check cbuffindex.
283     if ( cbuffindex .gt. 0) then
284     write(cunit) cbuffindex
285     write(cunit) k
286     write(cunit) (cbuff(ii), ii=1,cbuffindex)
287     endif
288 heimbach 1.12 c
289     if ( doPackDiag ) then
290     write(cunit2,rec=irectrue) globfldtmp2
291     write(cunit3,rec=irectrue) globfldtmp3
292     endif
293     c
294 heimbach 1.3 c -- end of k loop --
295     enddo
296     c -- end of j loop --
297     enddo
298     c -- end of irec loop --
299     enddo
300    
301 mlosch 1.13 do irec = nrec_nl*ny+1, ncvarrecs(ivartype)
302 heimbach 1.2 cgg do iobcs = 1, nobcs
303     cgg Need to solve for what iobcs would have been.
304 heimbach 1.3 iobcs= mod(irec-1,nobcs)+1
305 heimbach 1.2
306     call MDSREADFIELD_XZ_GL( fname, ctrlprec, 'RL',
307     & nr, globfldxz, irec, mythid)
308    
309 heimbach 1.3 CGG One special exception: barotropic velocity should be nondimensionalized
310     cgg differently. Probably introduce new variable.
311     if (iobcs .eq. 3 .or. iobcs .eq. 4) then
312     k = 1
313     do jp = 1,nPy
314     do bj = jtlo,jthi
315     do ip = 1,nPx
316     do bi = itlo,ithi
317     do i = imin,imax
318     #ifdef NO_CONTROL_BAROTROPIC_VELOCITY
319     if (.not. lxxadxx) then
320     cgg Get rid of any sensitivity to barotropic velocity.
321     globfldxz(i,bi,ip,bj,jp,k) = 0.
322     endif
323     #endif
324     enddo
325     enddo
326     enddo
327     enddo
328     enddo
329     endif
330    
331 heimbach 1.2 write(cunit) ncvarindex(ivartype)
332     write(cunit) 1
333     write(cunit) 1
334     do k = 1,nr
335 heimbach 1.12 irectrue = (irec-1)*nobcs*nr + (iobcs-1)*nr + k
336 heimbach 1.2 cbuffindex = 0
337     do jp = 1,nPy
338     do bj = jtlo,jthi
339     do ip = 1,nPx
340     do bi = itlo,ithi
341     do i = imin,imax
342 heimbach 1.3 if (globmskxz(i,bi,ip,bj,jp,k,iobcs) .ne. 0. ) then
343 heimbach 1.2 cbuffindex = cbuffindex + 1
344 heimbach 1.12 cph(
345     globfldtmp3(i,bi,ip,bj,jp) =
346     & globfldxz(i,bi,ip,bj,jp,k)
347     cph)
348 heimbach 1.2 #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
349     if (lxxadxx) then
350 jmc 1.11 cbuff(cbuffindex) =
351 heimbach 1.2 & globfldxz(i,bi,ip,bj,jp,k) *
352 heimbach 1.3 # ifdef CTRL_PACK_PRECISE
353     & sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
354     # else
355 heimbach 1.2 & sqrt(weightfld(k,iobcs))
356 heimbach 1.3 # endif
357 heimbach 1.2 else
358 jmc 1.11 cbuff(cbuffindex) =
359 heimbach 1.2 & globfldxz(i,bi,ip,bj,jp,k) /
360 heimbach 1.3 # ifdef CTRL_PACK_PRECISE
361     & sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
362     # else
363 heimbach 1.2 & sqrt(weightfld(k,iobcs))
364 heimbach 1.3 # endif
365 heimbach 1.2 endif
366 heimbach 1.12 cph(
367     globfldtmp2(i,bi,ip,bj,jp) = cbuff(cbuffindex)
368     cph)
369     #else /* ALLOW_NONDIMENSIONAL_CONTROL_IO undef */
370 heimbach 1.2 cbuff(cbuffindex) = globfldxz(i,bi,ip,bj,jp,k)
371 heimbach 1.12 #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
372 heimbach 1.2 endif
373     enddo
374     enddo
375     enddo
376     enddo
377     enddo
378     c --> check cbuffindex.
379     if ( cbuffindex .gt. 0) then
380     write(cunit) cbuffindex
381     write(cunit) k
382     write(cunit) (cbuff(ii), ii=1,cbuffindex)
383     endif
384 heimbach 1.3 c
385 heimbach 1.12 if ( doPackDiag ) then
386     write(cunit2,rec=irectrue) globfldtmp2
387     write(cunit3,rec=irectrue) globfldtmp3
388     endif
389     c
390 heimbach 1.3 c -- end of k loop --
391 heimbach 1.2 enddo
392     c -- end of iobcs loop --
393     cgg enddo
394     c -- end of irec loop --
395     enddo
396    
397     _END_MASTER( mythid )
398    
399 heimbach 1.10 #endif
400    
401 heimbach 1.2 return
402     end
403    
404    
405    
406    
407    

  ViewVC Help
Powered by ViewVC 1.1.22