/[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.4 - (hide annotations) (download)
Thu Jul 24 22:00:18 2003 UTC (20 years, 11 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint51k_post, checkpoint51l_post, checkpoint51f_post, checkpoint51j_post, checkpoint51l_pre, checkpoint51h_pre, branchpoint-genmake2, checkpoint51i_post, checkpoint51i_pre, checkpoint51e_post, checkpoint51f_pre, checkpoint51g_post, checkpoint51m_post
Branch point for: branch-genmake2, tg2-branch
Changes since 1.3: +6 -4 lines
bug fixes for 3d packing and I/O of sliced (xz/yz) fields
to increase I/O performance.

1 heimbach 1.2
2     #include "CTRL_CPPOPTIONS.h"
3    
4    
5     subroutine ctrl_set_pack_xz(
6 heimbach 1.3 & cunit, ivartype, fname, masktype,weighttype,
7 heimbach 1.2 & weightfld, lxxadxx, mythid)
8    
9     c ==================================================================
10     c SUBROUTINE ctrl_set_pack_xz
11     c ==================================================================
12     c
13     c o Compress the control vector such that only ocean points are
14     c written to file.
15     c
16 heimbach 1.3 c o Open boundary packing finalized :
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 heimbach 1.2 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 "cost.h"
37    
38     #ifdef ALLOW_ECCO_OPTIMIZATION
39     #include "optim.h"
40     #endif
41    
42     c == routine arguments ==
43    
44     integer cunit
45     integer ivartype
46     character*( 80) fname
47     character*( 9) masktype
48 heimbach 1.3 character*( 80) weighttype
49 heimbach 1.2 _RL weightfld( nr,nobcs )
50     logical lxxadxx
51     integer mythid
52    
53     c == local variables ==
54    
55     #ifndef ALLOW_ECCO_OPTIMIZATION
56     integer optimcycle
57     #endif
58    
59     integer bi,bj
60     integer ip,jp
61     integer i,j,k
62 heimbach 1.4 integer ii,jj,kk
63 heimbach 1.2 integer il
64 heimbach 1.3 integer irec,iobcs,nrec_nl
65 heimbach 1.2 integer itlo,ithi
66     integer jtlo,jthi
67     integer jmin,jmax
68     integer imin,imax
69    
70     integer cbuffindex
71     cgg(
72     integer igg
73     _RL gg
74 heimbach 1.3 character*(80) weightname
75 heimbach 1.2 cgg)
76    
77     _RL cbuff ( snx*nsx*npx*nsy*npy )
78     _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     #ifndef ALLOW_ECCO_OPTIMIZATION
93     optimcycle = 0
94     #endif
95    
96     jtlo = 1
97     jthi = nsy
98     itlo = 1
99     ithi = nsx
100     jmin = 1
101     jmax = sny
102     imin = 1
103     imax = snx
104    
105     c Initialise temporary file
106     do k = 1,nr
107     do jp = 1,nPy
108     do bj = jtlo,jthi
109     do ip = 1,nPx
110     do bi = itlo,ithi
111     do i = imin,imax
112     globfldxz(i,bi,ip,bj,jp,k) = 0. _d 0
113 heimbach 1.3 do iobcs=1,nobcs
114     globmskxz(i,bi,ip,bj,jp,k,iobcs) = 0. _d 0
115     enddo
116     enddo
117     enddo
118     enddo
119     enddo
120     enddo
121     enddo
122     c Initialise temporary file
123     do k = 1,nr
124     do jp = 1,nPy
125     do bj = jtlo,jthi
126     do j = jmin,jmax
127     do ip = 1,nPx
128     do bi = itlo,ithi
129     do i = imin,imax
130     globfld3d(i,bi,ip,j,bj,jp,k) = 0. _d 0
131     enddo
132 heimbach 1.2 enddo
133     enddo
134     enddo
135     enddo
136     enddo
137     enddo
138    
139     c-- Only the master thread will do I/O.
140     _BEGIN_MASTER( mythid )
141    
142 heimbach 1.3 do iobcs = 1, nobcs
143     call MDSREADFIELD_XZ_GL(
144     & masktype, ctrlprec, 'RL',
145     & Nr, globmskxz(1,1,1,1,1,1,iobcs), iobcs, mythid)
146     #ifdef CTRL_PACK_PRECISE
147     il=ilnblnk( weighttype)
148     write(weightname(1:80),'(80a)') ' '
149     write(weightname(1:80),'(a)') weighttype(1:il)
150     call MDSREADFIELD_XZ_GL(
151     & weightname, ctrlprec, 'RL',
152     & Nr, weightfldxz(1,1,1,1,1,1,iobcs), iobcs, mythid)
153     CGG One special exception: barotropic velocity should be nondimensionalized
154     cgg differently. Probably introduce new variable.
155     if (iobcs .eq. 3 .or. iobcs .eq. 4) then
156     k = 1
157     do jp = 1,nPy
158     do bj = jtlo,jthi
159     do ip = 1,nPx
160     do bi = itlo,ithi
161     do i = imin,imax
162     weightfldxz(i,bi,ip,bj,jp,k,iobcs) = wbaro
163     enddo
164     enddo
165     enddo
166     enddo
167     enddo
168     endif
169     #endif
170     enddo
171    
172     nrec_nl=int(ncvarrecs(ivartype)/sny)
173     do irec = 1, nrec_nl
174     call MDSREADFIELD_3D_GL( fname, ctrlprec, 'RL',
175     & nr, globfld3d, irec, mythid)
176     do j=1,sny
177     iobcs= mod((irec-1)*sny+j-1,nobcs)+1
178    
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     #ifdef NO_CONTROL_BAROTROPIC_VELOCITY
189     if (.not. lxxadxx) then
190     cgg Get rid of any sensitivity to barotropic velocity.
191     globfld3d(i,bi,ip,j,bj,jp,k) = 0.
192     endif
193     #endif
194     enddo
195     enddo
196     enddo
197     enddo
198     enddo
199     endif
200    
201     write(cunit) ncvarindex(ivartype)
202     write(cunit) 1
203     write(cunit) 1
204     do k = 1,nr
205     cbuffindex = 0
206     do jp = 1,nPy
207     do bj = jtlo,jthi
208     do ip = 1,nPx
209     do bi = itlo,ithi
210     do i = imin,imax
211     if (globmskxz(i,bi,ip,bj,jp,k,iobcs) .ne. 0. ) then
212     cbuffindex = cbuffindex + 1
213 heimbach 1.4 jj=mod((j-1)*nr+k-1,sny)+1
214     kk=int((j-1)*nr+K-1)/sny+1
215 heimbach 1.3 #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
216     if (lxxadxx) then
217     cbuff(cbuffindex) =
218 heimbach 1.4 & globfld3d(i,bi,ip,jj,bj,jp,kk) *
219 heimbach 1.3 # ifdef CTRL_PACK_PRECISE
220     & sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
221     # else
222     & sqrt(weightfld(k,iobcs))
223     # endif
224     else
225     cbuff(cbuffindex) =
226 heimbach 1.4 & globfld3d(i,bi,ip,jj,bj,jp,kk) /
227 heimbach 1.3 # ifdef CTRL_PACK_PRECISE
228     & sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
229     # else
230     & sqrt(weightfld(k,iobcs))
231     # endif
232     endif
233     #else /* ALLOW_NONDIMENSIONAL_CONTROL_IO undef */
234 heimbach 1.4 cbuff(cbuffindex) = globfld3d(i,bi,ip,jj,bj,jp,kk)
235 heimbach 1.3 #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
236     endif
237     enddo
238     enddo
239     enddo
240     enddo
241     enddo
242     c --> check cbuffindex.
243     if ( cbuffindex .gt. 0) then
244     write(cunit) cbuffindex
245     write(cunit) k
246     write(cunit) (cbuff(ii), ii=1,cbuffindex)
247     endif
248     c -- end of k loop --
249     enddo
250     c -- end of j loop --
251     enddo
252     c -- end of irec loop --
253     enddo
254    
255     do irec = nrec_nl*sny+1, ncvarrecs(ivartype)
256 heimbach 1.2 cgg do iobcs = 1, nobcs
257     cgg Need to solve for what iobcs would have been.
258 heimbach 1.3 iobcs= mod(irec-1,nobcs)+1
259 heimbach 1.2
260     call MDSREADFIELD_XZ_GL( fname, ctrlprec, 'RL',
261     & nr, globfldxz, irec, mythid)
262    
263 heimbach 1.3 CGG One special exception: barotropic velocity should be nondimensionalized
264     cgg differently. Probably introduce new variable.
265     if (iobcs .eq. 3 .or. iobcs .eq. 4) then
266     k = 1
267     do jp = 1,nPy
268     do bj = jtlo,jthi
269     do ip = 1,nPx
270     do bi = itlo,ithi
271     do i = imin,imax
272     #ifdef NO_CONTROL_BAROTROPIC_VELOCITY
273     if (.not. lxxadxx) then
274     cgg Get rid of any sensitivity to barotropic velocity.
275     globfldxz(i,bi,ip,bj,jp,k) = 0.
276     endif
277     #endif
278     enddo
279     enddo
280     enddo
281     enddo
282     enddo
283     endif
284    
285 heimbach 1.2 write(cunit) ncvarindex(ivartype)
286     write(cunit) 1
287     write(cunit) 1
288     do k = 1,nr
289     cbuffindex = 0
290     do jp = 1,nPy
291     do bj = jtlo,jthi
292     do ip = 1,nPx
293     do bi = itlo,ithi
294     do i = imin,imax
295 heimbach 1.3 if (globmskxz(i,bi,ip,bj,jp,k,iobcs) .ne. 0. ) then
296 heimbach 1.2 cbuffindex = cbuffindex + 1
297     #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
298     if (lxxadxx) then
299     cbuff(cbuffindex) =
300     & globfldxz(i,bi,ip,bj,jp,k) *
301 heimbach 1.3 # ifdef CTRL_PACK_PRECISE
302     & sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
303     # else
304 heimbach 1.2 & sqrt(weightfld(k,iobcs))
305 heimbach 1.3 # endif
306 heimbach 1.2 else
307     cbuff(cbuffindex) =
308     & globfldxz(i,bi,ip,bj,jp,k) /
309 heimbach 1.3 # ifdef CTRL_PACK_PRECISE
310     & sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
311     # else
312 heimbach 1.2 & sqrt(weightfld(k,iobcs))
313 heimbach 1.3 # endif
314 heimbach 1.2 endif
315     #else
316     cbuff(cbuffindex) = globfldxz(i,bi,ip,bj,jp,k)
317     #endif
318     endif
319     enddo
320     enddo
321     enddo
322     enddo
323     enddo
324     c --> check cbuffindex.
325     if ( cbuffindex .gt. 0) then
326     write(cunit) cbuffindex
327     write(cunit) k
328     write(cunit) (cbuff(ii), ii=1,cbuffindex)
329     endif
330 heimbach 1.3 c
331     c -- end of k loop --
332 heimbach 1.2 enddo
333     c -- end of iobcs loop --
334     cgg enddo
335     c -- end of irec loop --
336     enddo
337    
338     _END_MASTER( mythid )
339    
340     return
341     end
342    
343    
344    
345    
346    

  ViewVC Help
Powered by ViewVC 1.1.22