/[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.3 - (hide annotations) (download)
Tue Jun 24 16:07:07 2003 UTC (20 years, 11 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint51, checkpoint51d_post, checkpoint51b_pre, checkpoint51b_post, checkpoint51c_post, checkpoint51a_post
Changes since 1.2: +189 -15 lines
Merging for c51 vs. e34

1 heimbach 1.2
2     #include "CTRL_CPPOPTIONS.h"
3    
4    
5     subroutine ctrl_set_pack_yz(
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_yz
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 added :
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     integer ii
63     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     _RL cbuff ( nsx*npx*sny*nsy*npy )
77     _RL globfldyz ( nsx,npx,sny,nsy,npy,nr )
78 heimbach 1.3 _RL globfld3d ( snx,nsx,npx,sny,nsy,npy,nr )
79     _RL globmskyz ( nsx,npx,sny,nsy,npy,nr,nobcs )
80     #ifdef CTRL_PACK_PRECISE
81     _RL weightfldyz( nsx,npx,sny,nsy,npy,nr,nobcs )
82     #endif
83 heimbach 1.2
84     c == external ==
85    
86     integer ilnblnk
87     external ilnblnk
88    
89     c == end of interface ==
90    
91     #ifndef ALLOW_ECCO_OPTIMIZATION
92     optimcycle = 0
93     #endif
94    
95     jtlo = 1
96     jthi = nsy
97     itlo = 1
98     ithi = nsx
99     jmin = 1
100     jmax = sny
101     imin = 1
102     imax = snx
103    
104     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     globfldyz(bi,ip,j,bj,jp,k) = 0. _d 0
112 heimbach 1.3 do iobcs=1,nobcs
113     globmskyz(bi,ip,j,bj,jp,k,iobcs) = 0. _d 0
114     enddo
115     enddo
116     enddo
117     enddo
118     enddo
119     enddo
120     enddo
121     c Initialise temporary file
122     do k = 1,nr
123     do jp = 1,nPy
124     do bj = jtlo,jthi
125     do j = jmin,jmax
126     do ip = 1,nPx
127     do bi = itlo,ithi
128     do i = imin,imax
129     globfld3d(i,bi,ip,j,bj,jp,k) = 0. _d 0
130     enddo
131 heimbach 1.2 enddo
132     enddo
133     enddo
134     enddo
135     enddo
136     enddo
137    
138     c-- Only the master thread will do I/O.
139     _BEGIN_MASTER( mythid )
140    
141 heimbach 1.3 do iobcs=1,nobcs
142     call MDSREADFIELD_YZ_GL(
143     & masktype, ctrlprec, 'RL',
144     & Nr, globmskyz(1,1,1,1,1,1,iobcs), iobcs,mythid)
145     #ifdef CTRL_PACK_PRECISE
146     il=ilnblnk( weighttype)
147     write(weightname(1:80),'(80a)') ' '
148     write(weightname(1:80),'(a)') weighttype(1:il)
149     call MDSREADFIELD_YZ_GL(
150     & weightname, ctrlprec, 'RL',
151     & Nr, weightfldyz(1,1,1,1,1,1,iobcs), iobcs, mythid)
152     CGG One special exception: barotropic velocity should be nondimensionalized
153     cgg differently. Probably introduce new variable.
154     if (iobcs .eq. 3 .or. iobcs .eq. 4) then
155     k = 1
156     do jp = 1,nPy
157     do bj = jtlo,jthi
158     do j = jmin,jmax
159     do ip = 1,nPx
160     do bi = itlo,ithi
161     weightfldyz(bi,ip,j,bj,jp,k,iobcs) = wbaro
162     enddo
163     enddo
164     enddo
165     enddo
166     enddo
167     endif
168     #endif
169     enddo
170 heimbach 1.2
171 heimbach 1.3 nrec_nl=int(ncvarrecs(ivartype)/snx)
172     do irec = 1, nrec_nl
173 heimbach 1.2 cgg do iobcs = 1, nobcs
174     cgg Need to solve for what iobcs would have been.
175    
176 heimbach 1.3 call MDSREADFIELD_3D_GL( fname, ctrlprec, 'RL',
177     & nr, globfld3D, irec, mythid)
178    
179     do i=1,snx
180     iobcs= mod((irec-1)*snx+i-1,nobcs)+1
181    
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     #ifdef NO_CONTROL_BAROTROPIC_VELOCITY
192     if (.not. lxxadxx) then
193     cgg Get rid of any sensitivity to barotropic velocity.
194     globfld3d(i,bi,ip,j,bj,jp,k) = 0.
195     endif
196     #endif
197     enddo
198     enddo
199     enddo
200     enddo
201     enddo
202     endif
203    
204     write(cunit) ncvarindex(ivartype)
205     write(cunit) 1
206     write(cunit) 1
207     do k = 1,nr
208     cbuffindex = 0
209     do jp = 1,nPy
210     do bj = jtlo,jthi
211     do ip = 1,nPx
212     do bi = itlo,ithi
213     do j = jmin,jmax
214     if (globmskyz(bi,ip,j,bj,jp,k,iobcs) .ne. 0. ) then
215     cbuffindex = cbuffindex + 1
216     #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
217     if (lxxadxx) then
218     cbuff(cbuffindex) =
219     & globfld3d(i,bi,ip,j,bj,jp,k) *
220     #ifdef CTRL_PACK_PRECISE
221     & sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs))
222     #else
223     & sqrt(weightfld(k,iobcs))
224     #endif
225     else
226     cbuff(cbuffindex) =
227     & globfld3d(i,bi,ip,j,bj,jp,k) /
228     #ifdef CTRL_PACK_PRECISE
229     & sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs))
230     #else
231     & sqrt(weightfld(k,iobcs))
232     #endif
233     endif
234     #else /* ALLOW_NONDIMENSIONAL_CONTROL_IO undef */
235     cbuff(cbuffindex) = globfld3d(i,bi,ip,j,bj,jp,k)
236     #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     c -- end of k loop --
251     enddo
252     c -- end of iobcs loop --
253     cgg enddo
254     c -- end of i loop --
255     enddo
256     c -- end of irec loop --
257     enddo
258    
259     do irec = nrec_nl*snx+1, ncvarrecs(ivartype)
260     cgg do iobcs = 1, nobcs
261     cgg Need to solve for what iobcs would have been.
262     iobcs= mod(irec-1,nobcs)+1
263 heimbach 1.2
264     call MDSREADFIELD_YZ_GL( fname, ctrlprec, 'RL',
265     & nr, globfldyz, irec, mythid)
266    
267 heimbach 1.3 CGG One special exception: barotropic velocity should be nondimensionalized
268     cgg differently. Probably introduce new variable.
269     if (iobcs .eq. 3 .or. iobcs .eq. 4) then
270     k = 1
271     do jp = 1,nPy
272     do bj = jtlo,jthi
273     do j = jmin,jmax
274     do ip = 1,nPx
275     do bi = itlo,ithi
276     #ifdef NO_CONTROL_BAROTROPIC_VELOCITY
277     if (.not. lxxadxx) then
278     cgg Get rid of any sensitivity to barotropic velocity.
279     globfldyz(bi,ip,j,bj,jp,k) = 0.
280     endif
281     #endif
282     enddo
283     enddo
284     enddo
285     enddo
286     enddo
287     endif
288    
289 heimbach 1.2 write(cunit) ncvarindex(ivartype)
290     write(cunit) 1
291     write(cunit) 1
292     do k = 1,nr
293     cbuffindex = 0
294     do jp = 1,nPy
295     do bj = jtlo,jthi
296     do ip = 1,nPx
297     do bi = itlo,ithi
298     do j = jmin,jmax
299 heimbach 1.3 if (globmskyz(bi,ip,j,bj,jp,k,iobcs) .ne. 0. ) then
300 heimbach 1.2 cbuffindex = cbuffindex + 1
301     #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
302     if (lxxadxx) then
303     cbuff(cbuffindex) =
304     & globfldyz(bi,ip,j,bj,jp,k) *
305 heimbach 1.3 #ifdef CTRL_PACK_PRECISE
306     & sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs))
307     #else
308 heimbach 1.2 & sqrt(weightfld(k,iobcs))
309 heimbach 1.3 #endif
310 heimbach 1.2 else
311     cbuff(cbuffindex) =
312     & globfldyz(bi,ip,j,bj,jp,k) /
313 heimbach 1.3 #ifdef CTRL_PACK_PRECISE
314     & sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs))
315     #else
316 heimbach 1.2 & sqrt(weightfld(k,iobcs))
317 heimbach 1.3 #endif
318 heimbach 1.2 endif
319 heimbach 1.3 #else /* ALLOW_NONDIMENSIONAL_CONTROL_IO undef */
320 heimbach 1.2 cbuff(cbuffindex) = globfldyz(bi,ip,j,bj,jp,k)
321 heimbach 1.3 #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
322 heimbach 1.2 endif
323     enddo
324     enddo
325     enddo
326     enddo
327     enddo
328     c --> check cbuffindex.
329     if ( cbuffindex .gt. 0) then
330     write(cunit) cbuffindex
331     write(cunit) k
332     write(cunit) (cbuff(ii), ii=1,cbuffindex)
333     endif
334 heimbach 1.3 c
335     c -- end of k loop --
336 heimbach 1.2 enddo
337     c -- end of iobcs loop --
338     cgg enddo
339     c -- end of irec loop --
340     enddo
341    
342     _END_MASTER( mythid )
343    
344     return
345     end
346    

  ViewVC Help
Powered by ViewVC 1.1.22