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

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

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


Revision 1.15 - (hide annotations) (download)
Tue Jul 19 12:44:36 2011 UTC (12 years, 10 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint63g, checkpoint63p, checkpoint63q, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63a, checkpoint63b, checkpoint63c
Changes since 1.14: +7 -76 lines
remove unused variables and some obsolete code

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

  ViewVC Help
Powered by ViewVC 1.1.22