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

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

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


Revision 1.5 - (hide annotations) (download)
Thu Oct 23 04:41:40 2003 UTC (20 years, 7 months ago) by edhill
Branch: MAIN
CVS Tags: checkpoint51o_pre, checkpoint51n_post, checkpoint51n_pre, checkpoint51o_post, checkpoint51p_post
Branch point for: checkpoint51n_branch
Changes since 1.4: +4 -0 lines
 o added the [#include "AD_CONFIG.h"] statement to all files that need
   it for adjoint/tl #defines
 o re-worked the build logic in genmake2 to support AD_CONFIG.h
 o removed tools/genmake since it no longer works

1 edhill 1.5 C
2     C $Header: $
3     C $Name: $
4 heimbach 1.2
5 edhill 1.5 #include "AD_CONFIG.h"
6 heimbach 1.2 #include "CTRL_CPPOPTIONS.h"
7    
8    
9     subroutine ctrl_set_unpack_xz(
10 heimbach 1.3 & cunit, ivartype, fname, masktype, weighttype,
11 heimbach 1.2 & weightfld, nwetglobal, mythid)
12    
13     c ==================================================================
14     c SUBROUTINE ctrl_set_unpack_xz
15     c ==================================================================
16     c
17 heimbach 1.3 c o Unpack the control vector such that land points are filled in.
18     c
19     c o Open boundary packing added :
20     c gebbie@mit.edu, 18-Mar-2003
21     c
22     c changed: heimbach@mit.edu 17-Jun-2003
23     c merged Armin's changes to replace write of
24     c nr * globfld2d by 1 * globfld3d
25     c (ad hoc fix to speed up global I/O)
26 heimbach 1.2 c
27     c ==================================================================
28    
29     implicit none
30    
31     c == global variables ==
32    
33     #include "EEPARAMS.h"
34     #include "SIZE.h"
35     #include "PARAMS.h"
36     #include "GRID.h"
37    
38     #include "ctrl.h"
39     #include "cost.h"
40    
41     #ifdef ALLOW_ECCO_OPTIMIZATION
42     #include "optim.h"
43     #endif
44    
45     c == routine arguments ==
46    
47     integer cunit
48     integer ivartype
49     character*( 80) fname
50     character* (9) masktype
51 heimbach 1.3 character*( 80) weighttype
52 heimbach 1.2 _RL weightfld( nr,nobcs )
53     integer nwetglobal(nr,nobcs)
54     integer mythid
55    
56     c == local variables ==
57    
58     #ifndef ALLOW_ECCO_OPTIMIZATION
59     integer optimcycle
60     #endif
61    
62     integer bi,bj
63     integer ip,jp
64     integer i,j,k
65 heimbach 1.4 integer ii,jj,kk
66 heimbach 1.2 integer il
67 heimbach 1.3 integer irec,iobcs,nrec_nl
68 heimbach 1.2 integer itlo,ithi
69     integer jtlo,jthi
70     integer jmin,jmax
71     integer imin,imax
72    
73     integer cbuffindex
74    
75     _RL cbuff ( snx*nsx*npx*nsy*npy )
76     _RL globfldxz( snx,nsx,npx,nsy,npy,nr )
77 heimbach 1.3 _RL globfld3d( snx,nsx,npx,sny,nsy,npy,nr )
78     _RL globmskxz( snx,nsx,npx,nsy,npy,nr,nobcs )
79     #ifdef CTRL_UNPACK_PRECISE
80     _RL weightfldxz( snx,nsx,npx,nsy,npy,nr,nobcs )
81     #endif
82 heimbach 1.2
83     integer filenvartype
84     integer filenvarlength
85     character*(10) fileExpId
86     integer fileOptimCycle
87     integer filencbuffindex
88     _RL fileDummy
89     integer fileIg
90     integer fileJg
91     integer fileI
92     integer fileJ
93     integer filensx
94     integer filensy
95     integer filek
96     integer filencvarindex(maxcvars)
97     integer filencvarrecs(maxcvars)
98     integer filencvarxmax(maxcvars)
99     integer filencvarymax(maxcvars)
100     integer filencvarnrmax(maxcvars)
101     character*( 1) filencvargrd(maxcvars)
102     cgg(
103     integer igg
104     _RL gg
105 heimbach 1.3 character*(80) weightname
106 heimbach 1.2 cgg)
107    
108     c == external ==
109    
110     integer ilnblnk
111     external ilnblnk
112    
113     cc == end of interface ==
114    
115     jtlo = 1
116     jthi = nsy
117     itlo = 1
118     ithi = nsx
119     jmin = 1
120     jmax = sny
121     imin = 1
122     imax = snx
123    
124     c Initialise temporary file
125     do k = 1,nr
126     do jp = 1,nPy
127     do bj = jtlo,jthi
128     do ip = 1,nPx
129     do bi = itlo,ithi
130     do i = imin,imax
131     globfldxz(i,bi,ip,bj,jp,k) = 0. _d 0
132 heimbach 1.3 do iobcs=1,nobcs
133     globmskxz(i,bi,ip,bj,jp,k,iobcs) = 0. _d 0
134     enddo
135     enddo
136     enddo
137     enddo
138     enddo
139     enddo
140     enddo
141     c Initialise temporary file
142     do k = 1,nr
143     do jp = 1,nPy
144     do bj = jtlo,jthi
145     do j = jmin,jmax
146     do ip = 1,nPx
147     do bi = itlo,ithi
148     do i = imin,imax
149     globfld3d(i,bi,ip,j,bj,jp,k) = 0. _d 0
150     enddo
151 heimbach 1.2 enddo
152     enddo
153     enddo
154     enddo
155     enddo
156     enddo
157    
158     #ifndef ALLOW_ECCO_OPTIMIZATION
159     optimcycle = 0
160     #endif
161    
162     c-- Only the master thread will do I/O.
163     _BEGIN_MASTER( mythid )
164    
165 heimbach 1.3 do iobcs=1,nobcs
166     call MDSREADFIELD_XZ_GL(
167     & masktype, ctrlprec, 'RL',
168     & Nr, globmskxz(1,1,1,1,1,1,iobcs), iobcs,mythid)
169     #ifdef CTRL_UNPACK_PRECISE
170     il=ilnblnk( weighttype)
171     write(weightname(1:80),'(80a)') ' '
172     write(weightname(1:80),'(a)') weighttype(1:il)
173     call MDSREADFIELD_XZ_GL(
174     & weightname, ctrlprec, 'RL',
175     & Nr, weightfldxz(1,1,1,1,1,1,iobcs), iobcs, mythid)
176     CGG One special exception: barotropic velocity should be nondimensionalized
177     cgg differently. Probably introduce new variable.
178     if (iobcs .eq. 3 .or. iobcs .eq. 4) then
179     k = 1
180     do jp = 1,nPy
181     do bj = jtlo,jthi
182     do ip = 1,nPx
183     do bi = itlo,ithi
184     do i = imin,imax
185     weightfldxz(i,bi,ip,bj,jp,k,iobcs) = wbaro
186     enddo
187     enddo
188     enddo
189     enddo
190     enddo
191     endif
192     #endif /* CTRL_UNPACK_PRECISE */
193     enddo
194    
195     nrec_nl=int(ncvarrecs(ivartype)/sny)
196     do irec = 1, nrec_nl
197 heimbach 1.2 cgg do iobcs = 1, nobcs
198 heimbach 1.3 cgg And now back-calculate what iobcs should be.
199     do j=1,sny
200     iobcs= mod((irec-1)*sny+j-1,nobcs)+1
201    
202     read(cunit) filencvarindex(ivartype)
203     if (filencvarindex(ivartype) .NE. ncvarindex(ivartype))
204     & then
205     print *, 'ctrl-set_unpack:xz:WARNING: wrong ncvarindex ',
206     & filencvarindex(ivartype), ncvarindex(ivartype)
207     STOP 'in S/R ctrl_unpack'
208     endif
209     read(cunit) filej
210     read(cunit) filei
211     do k = 1, Nr
212     cbuffindex = nwetglobal(k,iobcs)
213     if ( cbuffindex .gt. 0 ) then
214     read(cunit) filencbuffindex
215     if (filencbuffindex .NE. cbuffindex) then
216     print *, 'WARNING: wrong cbuffindex ',
217     & filencbuffindex, cbuffindex
218     STOP 'in S/R ctrl_unpack'
219     endif
220     read(cunit) filek
221     if (filek .NE. k) then
222     print *, 'WARNING: wrong k ',
223     & filek, k
224     STOP 'in S/R ctrl_unpack'
225     endif
226     read(cunit) (cbuff(ii), ii=1,cbuffindex)
227     endif
228     cbuffindex = 0
229     do jp = 1,nPy
230     do bj = jtlo,jthi
231     do ip = 1,nPx
232     do bi = itlo,ithi
233     do i = imin,imax
234     if ( globmskxz(i,bi,ip,bj,jp,k,iobcs) .ne. 0. ) then
235     cbuffindex = cbuffindex + 1
236 heimbach 1.4 jj=mod((j-1)*nr+k-1,sny)+1
237     kk=int((j-1)*nr+k-1)/sny+1
238     globfld3d(i,bi,ip,jj,bj,jp,kk) =
239     & cbuff(cbuffindex)
240 heimbach 1.3 #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
241 heimbach 1.4 globfld3d(i,bi,ip,jj,bj,jp,kk) =
242     & globfld3d(i,bi,ip,jj,bj,jp,kk)/
243 heimbach 1.3 # ifdef CTRL_UNPACK_PRECISE
244     & sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
245     # else
246     & sqrt(weightfld(k,iobcs))
247     # endif
248     #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
249     else
250 heimbach 1.4 globfld3d(i,bi,ip,jj,bj,jp,kk) = 0. _d 0
251 heimbach 1.3 endif
252     enddo
253     enddo
254     enddo
255     enddo
256     enddo
257     c
258     c -- end of k loop --
259     enddo
260     c -- end of j loop --
261     enddo
262    
263     call MDSWRITEFIELD_3D_GL( fname, ctrlprec, 'RL',
264     & Nr, globfld3d, irec,
265     & optimcycle, mythid)
266    
267     c -- end of iobcs loop -- This loop removed. 3-28-02.
268     cgg enddo
269     c -- end of irec loop --
270     enddo
271 heimbach 1.2
272 heimbach 1.3 do irec = nrec_nl*sny+1, ncvarrecs(ivartype)
273     cgg do iobcs = 1, nobcs
274     cgg And now back-calculate what iobcs should be.
275     iobcs= mod(irec-1,nobcs)+1
276 heimbach 1.2
277     read(cunit) filencvarindex(ivartype)
278     if (filencvarindex(ivartype) .NE. ncvarindex(ivartype))
279     & then
280     print *, 'ctrl-set_unpack:xz:WARNING: wrong ncvarindex ',
281     & filencvarindex(ivartype), ncvarindex(ivartype)
282     STOP 'in S/R ctrl_unpack'
283     endif
284     read(cunit) filej
285     read(cunit) filei
286     do k = 1, Nr
287     cbuffindex = nwetglobal(k,iobcs)
288     if ( cbuffindex .gt. 0 ) then
289     read(cunit) filencbuffindex
290     if (filencbuffindex .NE. cbuffindex) then
291     print *, 'WARNING: wrong cbuffindex ',
292     & filencbuffindex, cbuffindex
293     STOP 'in S/R ctrl_unpack'
294     endif
295     read(cunit) filek
296     if (filek .NE. k) then
297     print *, 'WARNING: wrong k ',
298     & filek, k
299     STOP 'in S/R ctrl_unpack'
300     endif
301     read(cunit) (cbuff(ii), ii=1,cbuffindex)
302     endif
303     cbuffindex = 0
304     do jp = 1,nPy
305     do bj = jtlo,jthi
306     do ip = 1,nPx
307     do bi = itlo,ithi
308     do i = imin,imax
309 heimbach 1.3 if ( globmskxz(i,bi,ip,bj,jp,k,iobcs) .ne. 0. ) then
310 heimbach 1.2 cbuffindex = cbuffindex + 1
311     globfldxz(i,bi,ip,bj,jp,k) = cbuff(cbuffindex)
312     #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
313     globfldxz(i,bi,ip,bj,jp,k) =
314     & globfldxz(i,bi,ip,bj,jp,k)/
315 heimbach 1.3 # ifdef CTRL_UNPACK_PRECISE
316     & sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
317     # else
318 heimbach 1.2 & sqrt(weightfld(k,iobcs))
319 heimbach 1.3 # endif
320     #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
321 heimbach 1.2 else
322     globfldxz(i,bi,ip,bj,jp,k) = 0. _d 0
323     endif
324     enddo
325     enddo
326     enddo
327     enddo
328     enddo
329     c
330 heimbach 1.3 c -- end of k loop --
331 heimbach 1.2 enddo
332    
333     call MDSWRITEFIELD_XZ_GL( fname, ctrlprec, 'RL',
334     & Nr, globfldxz, irec,
335     & optimcycle, mythid)
336    
337     c -- end of iobcs loop -- This loop removed. 3-28-02.
338     cgg enddo
339     c -- end of irec loop --
340     enddo
341    
342     _END_MASTER( mythid )
343    
344     return
345     end
346    
347    
348    
349    
350    

  ViewVC Help
Powered by ViewVC 1.1.22