/[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.17 - (hide annotations) (download)
Thu Oct 9 00:49:27 2014 UTC (9 years, 8 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65l, checkpoint65m, checkpoint65f, checkpoint65g, checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65o, HEAD
Changes since 1.16: +2 -1 lines
- pkg/ctrl/CTRL_OBCS.h (new) : regroup all obcs ctrl variables.
- pkg/ctrl/ctrl.h, ctrl_dummy.h, ctrl_weights.h : rm obcs
  ctrl variables (now all in CTRL_OBCS.h).

- pkg/ctrl/ctrl_getobcse.F, ctrl_getobcsn.F, ctrl_getobcss.F,
  ctrl_getobcsw.F, ctrl_getrec.F, ctrl_init.F, ctrl_init_obcs_variables.F,
  ctrl_init_wet.F, ctrl_mask_set_xz.F, ctrl_mask_set_yz.F,
  ctrl_pack.F, ctrl_unpack.F, ctrl_readparms.F,
  ctrl_set_pack_xz.F, ctrl_set_pack_yz.F, ctrl_set_unpack_xz.F,
  ctrl_set_unpack_yz.F : add CPP brackets and CTRL_OBCS.h

- pkg/ctrl/ctrl_pack.F, ctrl_unpack.F : add CPP brackets

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

  ViewVC Help
Powered by ViewVC 1.1.22