/[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.7 - (hide annotations) (download)
Thu Nov 6 22:05:08 2003 UTC (20 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint52l_pre, hrcube4, hrcube5, checkpoint52d_pre, checkpoint52j_pre, checkpoint54d_post, checkpoint54e_post, checkpoint52l_post, checkpoint52k_post, checkpoint55, checkpoint54, checkpoint53, checkpoint52, checkpoint52f_post, checkpoint54f_post, checkpoint55i_post, checkpoint52i_pre, hrcube_1, hrcube_2, hrcube_3, checkpoint55c_post, checkpoint52e_pre, checkpoint52e_post, checkpoint53d_post, checkpoint52b_pre, checkpoint54b_post, checkpoint52m_post, checkpoint55g_post, checkpoint52b_post, checkpoint52c_post, checkpoint52f_pre, checkpoint55d_post, checkpoint54a_pre, checkpoint53c_post, checkpoint55d_pre, checkpoint55j_post, checkpoint54a_post, checkpoint55h_post, checkpoint55b_post, checkpoint53a_post, checkpoint55f_post, checkpoint52d_post, checkpoint53g_post, checkpoint52a_pre, checkpoint52i_post, checkpoint52h_pre, checkpoint53f_post, checkpoint52j_post, branch-netcdf, checkpoint52n_post, checkpoint53b_pre, checkpoint55a_post, checkpoint53b_post, checkpoint52a_post, ecco_c52_e35, checkpoint53d_pre, checkpoint55e_post, checkpoint54c_post, checkpoint51u_post
Branch point for: netcdf-sm0
Changes since 1.6: +8 -12 lines
o merging from ecco-branch
o cleaned some cross-dependencies and updated CPP options

1 heimbach 1.2
2     #include "CTRL_CPPOPTIONS.h"
3    
4     subroutine ctrl_set_pack_yz(
5 heimbach 1.3 & cunit, ivartype, fname, masktype, weighttype,
6 heimbach 1.2 & weightfld, lxxadxx, mythid)
7    
8     c ==================================================================
9     c SUBROUTINE ctrl_set_pack_yz
10     c ==================================================================
11     c
12     c o Compress the control vector such that only ocean points are
13     c written to file.
14     c
15 heimbach 1.3 c o Open boundary packing added :
16     c gebbie@mit.edu, 18-Mar-2003
17     c
18     c changed: heimbach@mit.edu 17-Jun-2003
19     c merged Armin's changes to replace write of
20     c nr * globfld2d by 1 * globfld3d
21     c (ad hoc fix to speed up global I/O)
22     c
23 heimbach 1.2 c ==================================================================
24    
25     implicit none
26    
27     c == global variables ==
28    
29     #include "EEPARAMS.h"
30     #include "SIZE.h"
31     #include "PARAMS.h"
32     #include "GRID.h"
33    
34     #include "ctrl.h"
35    
36     #ifdef ALLOW_ECCO_OPTIMIZATION
37     #include "optim.h"
38     #endif
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     c == local variables ==
52    
53     #ifndef ALLOW_ECCO_OPTIMIZATION
54     integer optimcycle
55     #endif
56    
57     integer bi,bj
58     integer ip,jp
59     integer i,j,k
60 heimbach 1.7 integer ii,jj,kk
61 heimbach 1.2 integer il
62 heimbach 1.3 integer irec,iobcs,nrec_nl
63 heimbach 1.2 integer itlo,ithi
64     integer jtlo,jthi
65     integer jmin,jmax
66     integer imin,imax
67    
68     integer cbuffindex
69     cgg(
70     integer igg
71     _RL gg
72 heimbach 1.3 character*(80) weightname
73 heimbach 1.2 cgg)
74 heimbach 1.7 real*4 cbuff ( nsx*npx*sny*nsy*npy )
75 heimbach 1.2 _RL globfldyz ( nsx,npx,sny,nsy,npy,nr )
76 heimbach 1.3 _RL globfld3d ( snx,nsx,npx,sny,nsy,npy,nr )
77     _RL globmskyz ( nsx,npx,sny,nsy,npy,nr,nobcs )
78     #ifdef CTRL_PACK_PRECISE
79     _RL weightfldyz( nsx,npx,sny,nsy,npy,nr,nobcs )
80     #endif
81 heimbach 1.2
82     c == external ==
83    
84     integer ilnblnk
85     external ilnblnk
86    
87     c == end of interface ==
88    
89     #ifndef ALLOW_ECCO_OPTIMIZATION
90     optimcycle = 0
91     #endif
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     c Initialise temporary file
103     do k = 1,nr
104     do jp = 1,nPy
105     do bj = jtlo,jthi
106     do j = jmin,jmax
107     do ip = 1,nPx
108     do bi = itlo,ithi
109     globfldyz(bi,ip,j,bj,jp,k) = 0. _d 0
110 heimbach 1.3 do iobcs=1,nobcs
111     globmskyz(bi,ip,j,bj,jp,k,iobcs) = 0. _d 0
112     enddo
113     enddo
114     enddo
115     enddo
116     enddo
117     enddo
118     enddo
119     c Initialise temporary file
120     do k = 1,nr
121     do jp = 1,nPy
122     do bj = jtlo,jthi
123     do j = jmin,jmax
124     do ip = 1,nPx
125     do bi = itlo,ithi
126     do i = imin,imax
127     globfld3d(i,bi,ip,j,bj,jp,k) = 0. _d 0
128     enddo
129 heimbach 1.2 enddo
130     enddo
131     enddo
132     enddo
133     enddo
134     enddo
135    
136     c-- Only the master thread will do I/O.
137     _BEGIN_MASTER( mythid )
138    
139 heimbach 1.3 do iobcs=1,nobcs
140     call MDSREADFIELD_YZ_GL(
141     & masktype, ctrlprec, 'RL',
142     & Nr, globmskyz(1,1,1,1,1,1,iobcs), iobcs,mythid)
143     #ifdef CTRL_PACK_PRECISE
144     il=ilnblnk( weighttype)
145     write(weightname(1:80),'(80a)') ' '
146     write(weightname(1:80),'(a)') weighttype(1:il)
147     call MDSREADFIELD_YZ_GL(
148     & weightname, ctrlprec, 'RL',
149     & Nr, weightfldyz(1,1,1,1,1,1,iobcs), iobcs, mythid)
150     CGG One special exception: barotropic velocity should be nondimensionalized
151     cgg differently. Probably introduce new variable.
152     if (iobcs .eq. 3 .or. iobcs .eq. 4) then
153     k = 1
154     do jp = 1,nPy
155     do bj = jtlo,jthi
156     do j = jmin,jmax
157     do ip = 1,nPx
158     do bi = itlo,ithi
159     weightfldyz(bi,ip,j,bj,jp,k,iobcs) = wbaro
160     enddo
161     enddo
162     enddo
163     enddo
164     enddo
165     endif
166     #endif
167     enddo
168 heimbach 1.2
169 heimbach 1.3 nrec_nl=int(ncvarrecs(ivartype)/snx)
170     do irec = 1, nrec_nl
171 heimbach 1.2 cgg do iobcs = 1, nobcs
172     cgg Need to solve for what iobcs would have been.
173    
174 heimbach 1.3 call MDSREADFIELD_3D_GL( fname, ctrlprec, 'RL',
175     & nr, globfld3D, irec, mythid)
176    
177     do i=1,snx
178     iobcs= mod((irec-1)*snx+i-1,nobcs)+1
179    
180     CGG One special exception: barotropic velocity should be nondimensionalized
181     cgg differently. Probably introduce new variable.
182     if (iobcs .eq. 3 .or. iobcs .eq. 4) then
183     k = 1
184     do jp = 1,nPy
185     do bj = jtlo,jthi
186     do j = jmin,jmax
187     do ip = 1,nPx
188     do bi = itlo,ithi
189     #ifdef NO_CONTROL_BAROTROPIC_VELOCITY
190     if (.not. lxxadxx) then
191     cgg Get rid of any sensitivity to barotropic velocity.
192     globfld3d(i,bi,ip,j,bj,jp,k) = 0.
193     endif
194     #endif
195     enddo
196     enddo
197     enddo
198     enddo
199     enddo
200     endif
201    
202     write(cunit) ncvarindex(ivartype)
203     write(cunit) 1
204     write(cunit) 1
205     do k = 1,nr
206     cbuffindex = 0
207     do jp = 1,nPy
208     do bj = jtlo,jthi
209     do ip = 1,nPx
210     do bi = itlo,ithi
211     do j = jmin,jmax
212 heimbach 1.7 ii=mod((i-1)*nr*sny+(k-1)*sny+j-1,snx)+1
213     jj=mod(((i-1)*nr*sny+(k-1)*sny+j-1)/snx,sny)+1
214     kk=int((i-1)*nr*sny+(k-1)*sny+j-1)/(snx*sny)+1
215 heimbach 1.3 if (globmskyz(bi,ip,j,bj,jp,k,iobcs) .ne. 0. ) then
216     cbuffindex = cbuffindex + 1
217     #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
218     if (lxxadxx) then
219     cbuff(cbuffindex) =
220 heimbach 1.7 & globfld3d(ii,bi,ip,jj,bj,jp,kk) *
221 heimbach 1.3 #ifdef CTRL_PACK_PRECISE
222     & sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs))
223     #else
224     & sqrt(weightfld(k,iobcs))
225     #endif
226     else
227     cbuff(cbuffindex) =
228 heimbach 1.7 & globfld3d(ii,bi,ip,jj,bj,jp,kk) /
229 heimbach 1.3 #ifdef CTRL_PACK_PRECISE
230     & sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs))
231     #else
232     & sqrt(weightfld(k,iobcs))
233     #endif
234     endif
235     #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     c -- end of k loop --
252     enddo
253     c -- end of iobcs loop --
254     cgg enddo
255     c -- end of i loop --
256     enddo
257     c -- end of irec loop --
258     enddo
259    
260     do irec = nrec_nl*snx+1, ncvarrecs(ivartype)
261     cgg do iobcs = 1, nobcs
262     cgg Need to solve for what iobcs would have been.
263     iobcs= mod(irec-1,nobcs)+1
264 heimbach 1.2
265     call MDSREADFIELD_YZ_GL( fname, ctrlprec, 'RL',
266     & nr, globfldyz, irec, mythid)
267    
268 heimbach 1.3 CGG One special exception: barotropic velocity should be nondimensionalized
269     cgg differently. Probably introduce new variable.
270     if (iobcs .eq. 3 .or. iobcs .eq. 4) then
271     k = 1
272     do jp = 1,nPy
273     do bj = jtlo,jthi
274     do j = jmin,jmax
275     do ip = 1,nPx
276     do bi = itlo,ithi
277     #ifdef NO_CONTROL_BAROTROPIC_VELOCITY
278     if (.not. lxxadxx) then
279     cgg Get rid of any sensitivity to barotropic velocity.
280     globfldyz(bi,ip,j,bj,jp,k) = 0.
281     endif
282     #endif
283     enddo
284     enddo
285     enddo
286     enddo
287     enddo
288     endif
289    
290 heimbach 1.2 write(cunit) ncvarindex(ivartype)
291     write(cunit) 1
292     write(cunit) 1
293     do k = 1,nr
294     cbuffindex = 0
295     do jp = 1,nPy
296     do bj = jtlo,jthi
297     do ip = 1,nPx
298     do bi = itlo,ithi
299     do j = jmin,jmax
300 heimbach 1.3 if (globmskyz(bi,ip,j,bj,jp,k,iobcs) .ne. 0. ) then
301 heimbach 1.2 cbuffindex = cbuffindex + 1
302     #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
303     if (lxxadxx) then
304     cbuff(cbuffindex) =
305     & globfldyz(bi,ip,j,bj,jp,k) *
306 heimbach 1.3 #ifdef CTRL_PACK_PRECISE
307     & sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs))
308     #else
309 heimbach 1.2 & sqrt(weightfld(k,iobcs))
310 heimbach 1.3 #endif
311 heimbach 1.2 else
312     cbuff(cbuffindex) =
313     & globfldyz(bi,ip,j,bj,jp,k) /
314 heimbach 1.3 #ifdef CTRL_PACK_PRECISE
315     & sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs))
316     #else
317 heimbach 1.2 & sqrt(weightfld(k,iobcs))
318 heimbach 1.3 #endif
319 heimbach 1.2 endif
320 heimbach 1.3 #else /* ALLOW_NONDIMENSIONAL_CONTROL_IO undef */
321 heimbach 1.2 cbuff(cbuffindex) = globfldyz(bi,ip,j,bj,jp,k)
322 heimbach 1.3 #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
323 heimbach 1.2 endif
324     enddo
325     enddo
326     enddo
327     enddo
328     enddo
329     c --> check cbuffindex.
330     if ( cbuffindex .gt. 0) then
331     write(cunit) cbuffindex
332     write(cunit) k
333     write(cunit) (cbuff(ii), ii=1,cbuffindex)
334     endif
335 heimbach 1.3 c
336     c -- end of k loop --
337 heimbach 1.2 enddo
338     c -- end of iobcs loop --
339     cgg enddo
340     c -- end of irec loop --
341     enddo
342    
343     _END_MASTER( mythid )
344    
345     return
346     end
347    

  ViewVC Help
Powered by ViewVC 1.1.22