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

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

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


Revision 1.15 - (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
Changes since 1.14: +2 -2 lines
finish removing unbalanced quote (single or double) in commented line

1 jmc 1.15 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_set_unpack_yz.F,v 1.14 2009/09/30 16:03:20 mlosch Exp $
2 jmc 1.12 C $Name: $
3 heimbach 1.2
4     #include "CTRL_CPPOPTIONS.h"
5    
6     subroutine ctrl_set_unpack_yz(
7 heimbach 1.3 & cunit, ivartype, fname, masktype, weighttype,
8 heimbach 1.2 & weightfld, nwetglobal, mythid)
9    
10     c ==================================================================
11     c SUBROUTINE ctrl_set_unpack_yz
12     c ==================================================================
13     c
14 heimbach 1.3 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 jmc 1.15 c merged changes from Armin to replace write of
21 heimbach 1.3 c nr * globfld2d by 1 * globfld3d
22     c (ad hoc fix to speed up global I/O)
23 heimbach 1.2 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 heimbach 1.3 character*( 80) weighttype
45 heimbach 1.2 _RL weightfld( nr,nobcs )
46     integer nwetglobal(nr,nobcs)
47     integer mythid
48    
49 heimbach 1.11 #ifndef EXCLUDE_CTRL_PACK
50 heimbach 1.2 c == local variables ==
51    
52 heimbach 1.13 logical lxxadxx
53    
54 heimbach 1.2 integer bi,bj
55     integer ip,jp
56     integer i,j,k
57 heimbach 1.7 integer ii,jj,kk
58 heimbach 1.2 integer il
59 heimbach 1.3 integer irec,iobcs,nrec_nl
60 heimbach 1.2 integer itlo,ithi
61     integer jtlo,jthi
62     integer jmin,jmax
63     integer imin,imax
64    
65     integer cbuffindex
66    
67 heimbach 1.13 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 heimbach 1.2 _RL globfldyz( nsx,npx,sny,nsy,npy,nr )
71 heimbach 1.3 _RL globfld3d( snx,nsx,npx,sny,nsy,npy,nr )
72     _RL globmskyz( nsx,npx,sny,nsy,npy,nr,nobcs )
73     #ifdef CTRL_UNPACK_PRECISE
74     _RL weightfldyz( nsx,npx,sny,nsy,npy,nr,nobcs )
75     #endif
76 heimbach 1.2
77     cgg(
78 heimbach 1.13 integer reclen, irectrue
79     integer cunit2, cunit3
80 heimbach 1.2 integer igg
81     _RL gg
82 heimbach 1.3 character*(80) weightname
83 heimbach 1.13 character*(80) cfile2, cfile3
84 heimbach 1.2 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 heimbach 1.13 lxxadxx = .TRUE.
103    
104 heimbach 1.2 c Initialise temporary file
105     do k = 1,nr
106     do jp = 1,nPy
107     do bj = jtlo,jthi
108     do j = jmin,jmax
109     do ip = 1,nPx
110     do bi = itlo,ithi
111 heimbach 1.13 globfldyz (bi,ip,j,bj,jp,k) = 0. _d 0
112     globfldtmp2(bi,ip,j,bj,jp) = 0.
113     globfldtmp3(bi,ip,j,bj,jp) = 0.
114 heimbach 1.3 do iobcs=1,nobcs
115     globmskyz(bi,ip,j,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 heimbach 1.2 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 heimbach 1.13 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 = nsx*npx*sny*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 heimbach 1.3 do iobcs=1,nobcs
172     call MDSREADFIELD_YZ_GL(
173     & masktype, ctrlprec, 'RL',
174     & Nr, globmskyz(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_YZ_GL(
180     & weightname, ctrlprec, 'RL',
181     & Nr, weightfldyz(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 j = jmin,jmax
189     do ip = 1,nPx
190     do bi = itlo,ithi
191 heimbach 1.10 cph weightfldyz(bi,ip,j,bj,jp,k,iobcs) = wbaro
192 heimbach 1.3 enddo
193     enddo
194     enddo
195     enddo
196     enddo
197     endif
198     #endif
199     enddo
200    
201 mlosch 1.14 if ( useSingleCPUio ) then
202     C MDSWRITEFIELD_YZ_GL does not know about useSingleCPUio, so the faster
203     C method that works for .not.useSingleCPUio cannot be used
204     nrec_nl = 0
205     else
206     nrec_nl = int(ncvarrecs(ivartype)/Nx)
207     endif
208 heimbach 1.3 do irec = 1, nrec_nl
209 heimbach 1.2 cgg do iobcs = 1, nobcs
210     cgg And now back-calculate what iobcs should be.
211 heimbach 1.3 do i=1,snx
212     iobcs= mod((irec-1)*snx+i-1,nobcs)+1
213    
214     read(cunit) filencvarindex(ivartype)
215     if (filencvarindex(ivartype) .NE. ncvarindex(ivartype))
216     & then
217     print *, 'ctrl_set_unpack_yz:WARNING: wrong ncvarindex ',
218     & filencvarindex(ivartype), ncvarindex(ivartype)
219     STOP 'in S/R ctrl_unpack'
220     endif
221     read(cunit) filej
222     read(cunit) filei
223     do k = 1, Nr
224 heimbach 1.13 irectrue = (irec-1)*nobcs*nr + (iobcs-1)*nr + k
225 jmc 1.12 cbuffindex = nwetglobal(k,iobcs)
226 heimbach 1.3 if ( cbuffindex .gt. 0 ) then
227     read(cunit) filencbuffindex
228     if (filencbuffindex .NE. cbuffindex) then
229     print *, 'WARNING: wrong cbuffindex ',
230     & filencbuffindex, cbuffindex
231     STOP 'in S/R ctrl_unpack'
232     endif
233     read(cunit) filek
234     if (filek .NE. k) then
235     print *, 'WARNING: wrong k ',
236     & filek, k
237     STOP 'in S/R ctrl_unpack'
238     endif
239     read(cunit) (cbuff(ii), ii=1,cbuffindex)
240     endif
241     cbuffindex = 0
242     do jp = 1,nPy
243     do bj = jtlo,jthi
244     do j = jmin,jmax
245     do ip = 1,nPx
246     do bi = itlo,ithi
247 heimbach 1.7 ii=mod((i-1)*nr*sny+(k-1)*sny+j-1,snx)+1
248     jj=mod(((i-1)*nr*sny+(k-1)*sny+j-1)/snx,sny)+1
249     kk=int((i-1)*nr*sny+(k-1)*sny+j-1)/(snx*sny)+1
250 heimbach 1.3 if ( globmskyz(bi,ip,j,bj,jp,k,iobcs) .ne. 0. ) then
251     cbuffindex = cbuffindex + 1
252 jmc 1.12 globfld3d(ii,bi,ip,jj,bj,jp,kk) =
253 heimbach 1.4 & cbuff(cbuffindex)
254 heimbach 1.13 cph(
255     globfldtmp2(bi,ip,jj,bj,jp) =
256     & cbuff(cbuffindex)
257     cph)
258 heimbach 1.3 #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
259 jmc 1.12 globfld3d(ii,bi,ip,jj,bj,jp,kk) =
260 heimbach 1.7 & globfld3d(ii,bi,ip,jj,bj,jp,kk)/
261 heimbach 1.3 # ifdef CTRL_UNPACK_PRECISE
262     & sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs))
263     # else
264     & sqrt(weightfld(k,iobcs))
265     # endif
266     #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
267     else
268 heimbach 1.7 globfld3d(ii,bi,ip,jj,bj,jp,kk) = 0. _d 0
269 heimbach 1.3 endif
270 heimbach 1.13 cph(
271     globfldtmp3(bi,ip,jj,bj,jp) =
272     & globfld3d(ii,bi,ip,jj,bj,jp,kk)
273     cph)
274 heimbach 1.3 enddo
275     enddo
276     enddo
277     enddo
278     enddo
279     c
280 heimbach 1.13 if ( doPackDiag ) then
281     write(cunit2,rec=irectrue) globfldtmp2
282     write(cunit3,rec=irectrue) globfldtmp3
283     endif
284     c
285 heimbach 1.3 c -- end of k loop --
286     enddo
287     c -- end of i loop --
288     enddo
289 heimbach 1.2
290 heimbach 1.3 call MDSWRITEFIELD_3d_GL( fname, ctrlprec, 'RL',
291     & Nr, globfld3d, irec,
292     & optimcycle, mythid)
293    
294     c -- end of iobcs loop -- This loop has been removed.
295     cgg enddo
296     c -- end of irec loop --
297     enddo
298    
299 mlosch 1.14 do irec = nrec_nl*nx+1,ncvarrecs(ivartype)
300 heimbach 1.3 iobcs= mod(irec-1,nobcs)+1
301 heimbach 1.2
302     read(cunit) filencvarindex(ivartype)
303     if (filencvarindex(ivartype) .NE. ncvarindex(ivartype))
304     & then
305     print *, 'ctrl_set_unpack_yz:WARNING: wrong ncvarindex ',
306     & filencvarindex(ivartype), ncvarindex(ivartype)
307     STOP 'in S/R ctrl_unpack'
308     endif
309     read(cunit) filej
310     read(cunit) filei
311     do k = 1, Nr
312 heimbach 1.13 irectrue = (irec-1)*nobcs*nr + (iobcs-1)*nr + k
313 jmc 1.12 cbuffindex = nwetglobal(k,iobcs)
314 heimbach 1.2 if ( cbuffindex .gt. 0 ) then
315     read(cunit) filencbuffindex
316     if (filencbuffindex .NE. cbuffindex) then
317     print *, 'WARNING: wrong cbuffindex ',
318     & filencbuffindex, cbuffindex
319     STOP 'in S/R ctrl_unpack'
320     endif
321     read(cunit) filek
322     if (filek .NE. k) then
323     print *, 'WARNING: wrong k ',
324     & filek, k
325     STOP 'in S/R ctrl_unpack'
326     endif
327     read(cunit) (cbuff(ii), ii=1,cbuffindex)
328     endif
329     cbuffindex = 0
330     do jp = 1,nPy
331     do bj = jtlo,jthi
332     do j = jmin,jmax
333     do ip = 1,nPx
334     do bi = itlo,ithi
335 heimbach 1.3 if ( globmskyz(bi,ip,j,bj,jp,k,iobcs) .ne. 0. ) then
336 heimbach 1.2 cbuffindex = cbuffindex + 1
337     globfldyz(bi,ip,j,bj,jp,k) = cbuff(cbuffindex)
338 heimbach 1.13 cph(
339     globfldtmp2(bi,ip,j,bj,jp) = cbuff(cbuffindex)
340     cph)
341 heimbach 1.2 #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
342 jmc 1.12 globfldyz(bi,ip,j,bj,jp,k) =
343 heimbach 1.2 & globfldyz(bi,ip,j,bj,jp,k)/
344 heimbach 1.3 # ifdef CTRL_UNPACK_PRECISE
345     & sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs))
346     # else
347 heimbach 1.2 & sqrt(weightfld(k,iobcs))
348 heimbach 1.3 # endif
349     #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
350 heimbach 1.2 else
351     globfldyz(bi,ip,j,bj,jp,k) = 0. _d 0
352     endif
353 heimbach 1.13 cph(
354     globfldtmp3(bi,ip,j,bj,jp) =
355     & globfldyz(bi,ip,j,bj,jp,k)
356     cph)
357 heimbach 1.2 enddo
358     enddo
359     enddo
360     enddo
361     enddo
362     c
363 heimbach 1.3 c -- end of k loop
364 heimbach 1.2 enddo
365 jmc 1.12
366 heimbach 1.2 call MDSWRITEFIELD_YZ_GL( fname, ctrlprec, 'RL',
367     & Nr, globfldyz, irec,
368     & optimcycle, mythid)
369    
370     c -- end of iobcs loop -- This loop has been removed.
371     cgg enddo
372     c -- end of irec loop --
373     enddo
374    
375     _END_MASTER( mythid )
376    
377 heimbach 1.11 #endif
378    
379 heimbach 1.2 return
380     end
381    

  ViewVC Help
Powered by ViewVC 1.1.22