/[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.11 - (hide annotations) (download)
Thu Jun 14 18:55:36 2007 UTC (16 years, 11 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59h
Changes since 1.10: +3 -0 lines
Exclude global arrays if we dont need/want them
(thought we had checked this in a while ago, but apparently not)

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     #include "optim.h"
35    
36     c == routine arguments ==
37    
38     integer cunit
39     integer ivartype
40     character*( 80) fname
41     character* (9) masktype
42 heimbach 1.3 character*( 80) weighttype
43 heimbach 1.2 _RL weightfld( nr,nobcs )
44     integer nwetglobal(nr,nobcs)
45     integer mythid
46    
47 heimbach 1.11 #ifndef EXCLUDE_CTRL_PACK
48 heimbach 1.2 c == local variables ==
49    
50     integer bi,bj
51     integer ip,jp
52     integer i,j,k
53 heimbach 1.4 integer ii,jj,kk
54 heimbach 1.2 integer il
55 heimbach 1.3 integer irec,iobcs,nrec_nl
56 heimbach 1.2 integer itlo,ithi
57     integer jtlo,jthi
58     integer jmin,jmax
59     integer imin,imax
60    
61     integer cbuffindex
62    
63 heimbach 1.7 real*4 cbuff ( snx*nsx*npx*nsy*npy )
64 heimbach 1.2 _RL globfldxz( snx,nsx,npx,nsy,npy,nr )
65 heimbach 1.3 _RL globfld3d( snx,nsx,npx,sny,nsy,npy,nr )
66     _RL globmskxz( snx,nsx,npx,nsy,npy,nr,nobcs )
67     #ifdef CTRL_UNPACK_PRECISE
68     _RL weightfldxz( snx,nsx,npx,nsy,npy,nr,nobcs )
69     #endif
70 heimbach 1.2
71     cgg(
72     integer igg
73     _RL gg
74 heimbach 1.3 character*(80) weightname
75 heimbach 1.2 cgg)
76    
77     c == external ==
78    
79     integer ilnblnk
80     external ilnblnk
81    
82     cc == end of interface ==
83    
84     jtlo = 1
85     jthi = nsy
86     itlo = 1
87     ithi = nsx
88     jmin = 1
89     jmax = sny
90     imin = 1
91     imax = snx
92    
93     c Initialise temporary file
94     do k = 1,nr
95     do jp = 1,nPy
96     do bj = jtlo,jthi
97     do ip = 1,nPx
98     do bi = itlo,ithi
99     do i = imin,imax
100     globfldxz(i,bi,ip,bj,jp,k) = 0. _d 0
101 heimbach 1.3 do iobcs=1,nobcs
102     globmskxz(i,bi,ip,bj,jp,k,iobcs) = 0. _d 0
103     enddo
104     enddo
105     enddo
106     enddo
107     enddo
108     enddo
109     enddo
110     c Initialise temporary file
111     do k = 1,nr
112     do jp = 1,nPy
113     do bj = jtlo,jthi
114     do j = jmin,jmax
115     do ip = 1,nPx
116     do bi = itlo,ithi
117     do i = imin,imax
118     globfld3d(i,bi,ip,j,bj,jp,k) = 0. _d 0
119     enddo
120 heimbach 1.2 enddo
121     enddo
122     enddo
123     enddo
124     enddo
125     enddo
126    
127     c-- Only the master thread will do I/O.
128     _BEGIN_MASTER( mythid )
129    
130 heimbach 1.3 do iobcs=1,nobcs
131     call MDSREADFIELD_XZ_GL(
132     & masktype, ctrlprec, 'RL',
133     & Nr, globmskxz(1,1,1,1,1,1,iobcs), iobcs,mythid)
134     #ifdef CTRL_UNPACK_PRECISE
135     il=ilnblnk( weighttype)
136     write(weightname(1:80),'(80a)') ' '
137     write(weightname(1:80),'(a)') weighttype(1:il)
138     call MDSREADFIELD_XZ_GL(
139     & weightname, ctrlprec, 'RL',
140     & Nr, weightfldxz(1,1,1,1,1,1,iobcs), iobcs, mythid)
141     CGG One special exception: barotropic velocity should be nondimensionalized
142     cgg differently. Probably introduce new variable.
143     if (iobcs .eq. 3 .or. iobcs .eq. 4) then
144     k = 1
145     do jp = 1,nPy
146     do bj = jtlo,jthi
147     do ip = 1,nPx
148     do bi = itlo,ithi
149     do i = imin,imax
150 heimbach 1.10 cph weightfldxz(i,bi,ip,bj,jp,k,iobcs) = wbaro
151 heimbach 1.3 enddo
152     enddo
153     enddo
154     enddo
155     enddo
156     endif
157     #endif /* CTRL_UNPACK_PRECISE */
158     enddo
159    
160     nrec_nl=int(ncvarrecs(ivartype)/sny)
161     do irec = 1, nrec_nl
162 heimbach 1.2 cgg do iobcs = 1, nobcs
163 heimbach 1.3 cgg And now back-calculate what iobcs should be.
164     do j=1,sny
165     iobcs= mod((irec-1)*sny+j-1,nobcs)+1
166    
167     read(cunit) filencvarindex(ivartype)
168     if (filencvarindex(ivartype) .NE. ncvarindex(ivartype))
169     & then
170     print *, 'ctrl-set_unpack:xz:WARNING: wrong ncvarindex ',
171     & filencvarindex(ivartype), ncvarindex(ivartype)
172     STOP 'in S/R ctrl_unpack'
173     endif
174     read(cunit) filej
175     read(cunit) filei
176     do k = 1, Nr
177     cbuffindex = nwetglobal(k,iobcs)
178     if ( cbuffindex .gt. 0 ) then
179     read(cunit) filencbuffindex
180     if (filencbuffindex .NE. cbuffindex) then
181     print *, 'WARNING: wrong cbuffindex ',
182     & filencbuffindex, cbuffindex
183     STOP 'in S/R ctrl_unpack'
184     endif
185     read(cunit) filek
186     if (filek .NE. k) then
187     print *, 'WARNING: wrong k ',
188     & filek, k
189     STOP 'in S/R ctrl_unpack'
190     endif
191     read(cunit) (cbuff(ii), ii=1,cbuffindex)
192     endif
193     cbuffindex = 0
194 heimbach 1.7 jj=mod((j-1)*nr+k-1,sny)+1
195     kk=int((j-1)*nr+k-1)/sny+1
196 heimbach 1.3 do jp = 1,nPy
197     do bj = jtlo,jthi
198     do ip = 1,nPx
199     do bi = itlo,ithi
200     do i = imin,imax
201     if ( globmskxz(i,bi,ip,bj,jp,k,iobcs) .ne. 0. ) then
202     cbuffindex = cbuffindex + 1
203 heimbach 1.4 globfld3d(i,bi,ip,jj,bj,jp,kk) =
204     & cbuff(cbuffindex)
205 heimbach 1.3 #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
206 heimbach 1.4 globfld3d(i,bi,ip,jj,bj,jp,kk) =
207     & globfld3d(i,bi,ip,jj,bj,jp,kk)/
208 heimbach 1.3 # ifdef CTRL_UNPACK_PRECISE
209     & sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
210     # else
211     & sqrt(weightfld(k,iobcs))
212     # endif
213     #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
214     else
215 heimbach 1.4 globfld3d(i,bi,ip,jj,bj,jp,kk) = 0. _d 0
216 heimbach 1.3 endif
217     enddo
218     enddo
219     enddo
220     enddo
221     enddo
222     c
223     c -- end of k loop --
224     enddo
225     c -- end of j loop --
226     enddo
227    
228     call MDSWRITEFIELD_3D_GL( fname, ctrlprec, 'RL',
229     & Nr, globfld3d, irec,
230     & optimcycle, mythid)
231    
232     c -- end of iobcs loop -- This loop removed. 3-28-02.
233     cgg enddo
234     c -- end of irec loop --
235     enddo
236 heimbach 1.2
237 heimbach 1.3 do irec = nrec_nl*sny+1, ncvarrecs(ivartype)
238     cgg do iobcs = 1, nobcs
239     cgg And now back-calculate what iobcs should be.
240     iobcs= mod(irec-1,nobcs)+1
241 heimbach 1.2
242     read(cunit) filencvarindex(ivartype)
243     if (filencvarindex(ivartype) .NE. ncvarindex(ivartype))
244     & then
245     print *, 'ctrl-set_unpack:xz:WARNING: wrong ncvarindex ',
246     & filencvarindex(ivartype), ncvarindex(ivartype)
247     STOP 'in S/R ctrl_unpack'
248     endif
249     read(cunit) filej
250     read(cunit) filei
251     do k = 1, Nr
252     cbuffindex = nwetglobal(k,iobcs)
253     if ( cbuffindex .gt. 0 ) then
254     read(cunit) filencbuffindex
255     if (filencbuffindex .NE. cbuffindex) then
256     print *, 'WARNING: wrong cbuffindex ',
257     & filencbuffindex, cbuffindex
258     STOP 'in S/R ctrl_unpack'
259     endif
260     read(cunit) filek
261     if (filek .NE. k) then
262     print *, 'WARNING: wrong k ',
263     & filek, k
264     STOP 'in S/R ctrl_unpack'
265     endif
266     read(cunit) (cbuff(ii), ii=1,cbuffindex)
267     endif
268     cbuffindex = 0
269     do jp = 1,nPy
270     do bj = jtlo,jthi
271     do ip = 1,nPx
272     do bi = itlo,ithi
273     do i = imin,imax
274 heimbach 1.3 if ( globmskxz(i,bi,ip,bj,jp,k,iobcs) .ne. 0. ) then
275 heimbach 1.2 cbuffindex = cbuffindex + 1
276     globfldxz(i,bi,ip,bj,jp,k) = cbuff(cbuffindex)
277     #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
278     globfldxz(i,bi,ip,bj,jp,k) =
279     & globfldxz(i,bi,ip,bj,jp,k)/
280 heimbach 1.3 # ifdef CTRL_UNPACK_PRECISE
281     & sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
282     # else
283 heimbach 1.2 & sqrt(weightfld(k,iobcs))
284 heimbach 1.3 # endif
285     #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
286 heimbach 1.2 else
287     globfldxz(i,bi,ip,bj,jp,k) = 0. _d 0
288     endif
289     enddo
290     enddo
291     enddo
292     enddo
293     enddo
294     c
295 heimbach 1.3 c -- end of k loop --
296 heimbach 1.2 enddo
297    
298     call MDSWRITEFIELD_XZ_GL( fname, ctrlprec, 'RL',
299     & Nr, globfldxz, irec,
300     & optimcycle, mythid)
301    
302     c -- end of iobcs loop -- This loop removed. 3-28-02.
303     cgg enddo
304     c -- end of irec loop --
305     enddo
306    
307     _END_MASTER( mythid )
308    
309 heimbach 1.11 #endif
310    
311 heimbach 1.2 return
312     end
313    
314    
315    
316    
317    

  ViewVC Help
Powered by ViewVC 1.1.22