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

  ViewVC Help
Powered by ViewVC 1.1.22