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

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

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


Revision 1.17 - (show annotations) (download)
Tue Jul 19 12:44:36 2011 UTC (12 years, 10 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint63g, checkpoint63p, checkpoint63q, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63a, checkpoint63b, checkpoint63c
Changes since 1.16: +9 -34 lines
remove unused variables and some obsolete code

1 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_set_unpack_xz.F,v 1.16 2011/05/05 09:22:57 mlosch Exp $
2 C $Name: $
3
4 #include "CTRL_CPPOPTIONS.h"
5
6 subroutine ctrl_set_unpack_xz(
7 & cunit, ivartype, fname, masktype, weighttype,
8 & weightfld, nwetglobal, mythid)
9
10 c ==================================================================
11 c SUBROUTINE ctrl_set_unpack_xz
12 c ==================================================================
13 c
14 c o Unpack the control vector such that land points are filled in.
15 c
16 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 changes from Armin 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 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 "optim.h"
37
38 c == routine arguments ==
39
40 integer cunit
41 integer ivartype
42 character*( 80) fname
43 character* (9) masktype
44 character*( 80) weighttype
45 _RL weightfld( nr,nobcs )
46 integer nwetglobal(nr,nobcs)
47 integer mythid
48
49 #ifndef EXCLUDE_CTRL_PACK
50 c == local variables ==
51
52 logical lxxadxx
53
54 integer bi,bj
55 integer ip,jp
56 integer i,j,k
57 integer ii,jj,kk
58 integer irec,iobcs,nrec_nl
59 integer itlo,ithi
60 integer jtlo,jthi
61 integer jmin,jmax
62 integer imin,imax
63
64 integer cbuffindex
65
66 real*4 cbuff ( snx*nsx*npx*nsy*npy )
67 real*4 globfldtmp2( snx,nsx,npx,nsy,npy )
68 real*4 globfldtmp3( snx,nsx,npx,nsy,npy )
69 _RL globfldxz( snx,nsx,npx,nsy,npy,nr )
70 _RL globfld3d( snx,nsx,npx,sny,nsy,npy,nr )
71 _RL globmskxz( snx,nsx,npx,nsy,npy,nr,nobcs )
72
73 integer reclen, irectrue
74 integer cunit2, cunit3
75 character*(80) cfile2, cfile3
76
77 #ifdef CTRL_UNPACK_PRECISE
78 integer il
79 character*(80) weightname
80 _RL weightfldxz( snx,nsx,npx,nsy,npy,nr,nobcs )
81
82 c == external ==
83
84 integer ilnblnk
85 external ilnblnk
86 #endif
87
88 cc == end of interface ==
89
90 jtlo = 1
91 jthi = nsy
92 itlo = 1
93 ithi = nsx
94 jmin = 1
95 jmax = sny
96 imin = 1
97 imax = snx
98
99 lxxadxx = .TRUE.
100
101 c Initialise temporary file
102 do k = 1,nr
103 do jp = 1,nPy
104 do bj = jtlo,jthi
105 do ip = 1,nPx
106 do bi = itlo,ithi
107 do i = imin,imax
108 globfldxz (i,bi,ip,bj,jp,k) = 0. _d 0
109 globfldtmp2(i,bi,ip,bj,jp) = 0.
110 globfldtmp3(i,bi,ip,bj,jp) = 0.
111 do iobcs=1,nobcs
112 globmskxz(i,bi,ip,bj,jp,k,iobcs) = 0. _d 0
113 enddo
114 enddo
115 enddo
116 enddo
117 enddo
118 enddo
119 enddo
120 c Initialise temporary file
121 do k = 1,nr
122 do jp = 1,nPy
123 do bj = jtlo,jthi
124 do j = jmin,jmax
125 do ip = 1,nPx
126 do bi = itlo,ithi
127 do i = imin,imax
128 globfld3d(i,bi,ip,j,bj,jp,k) = 0. _d 0
129 enddo
130 enddo
131 enddo
132 enddo
133 enddo
134 enddo
135 enddo
136
137 c-- Only the master thread will do I/O.
138 _BEGIN_MASTER( mythid )
139
140 if ( doPackDiag ) then
141 write(cfile2(1:80),'(80a)') ' '
142 write(cfile3(1:80),'(80a)') ' '
143 if ( lxxadxx ) then
144 write(cfile2(1:80),'(a,I2.2,a,I4.4,a)')
145 & 'diag_unpack_nondim_ctrl_',
146 & ivartype, '_', optimcycle, '.bin'
147 write(cfile3(1:80),'(a,I2.2,a,I4.4,a)')
148 & 'diag_unpack_dimens_ctrl_',
149 & ivartype, '_', optimcycle, '.bin'
150 else
151 write(cfile2(1:80),'(a,I2.2,a,I4.4,a)')
152 & 'diag_unpack_nondim_grad_',
153 & ivartype, '_', optimcycle, '.bin'
154 write(cfile3(1:80),'(a,I2.2,a,I4.4,a)')
155 & 'diag_unpack_dimens_grad_',
156 & ivartype, '_', optimcycle, '.bin'
157 endif
158
159 reclen = snx*nsx*npx*nsy*npy*4
160 call mdsfindunit( cunit2, mythid )
161 open( cunit2, file=cfile2, status='unknown',
162 & access='direct', recl=reclen )
163 call mdsfindunit( cunit3, mythid )
164 open( cunit3, file=cfile3, status='unknown',
165 & access='direct', recl=reclen )
166 endif
167
168 do iobcs=1,nobcs
169 call MDSREADFIELD_XZ_GL(
170 & masktype, ctrlprec, 'RL',
171 & Nr, globmskxz(1,1,1,1,1,1,iobcs), iobcs,mythid)
172 #ifdef CTRL_UNPACK_PRECISE
173 il=ilnblnk( weighttype)
174 write(weightname(1:80),'(80a)') ' '
175 write(weightname(1:80),'(a)') weighttype(1:il)
176 call MDSREADFIELD_XZ_GL(
177 & weightname, ctrlprec, 'RL',
178 & Nr, weightfldxz(1,1,1,1,1,1,iobcs), iobcs, mythid)
179 #endif /* CTRL_UNPACK_PRECISE */
180 enddo
181
182 if ( useSingleCPUio ) then
183 C MDSWRITEFIELD_XZ_GL does not know about useSingleCPUio, so the faster
184 C method that works for .not.useSingleCPUio cannot be used
185 nrec_nl = 0
186 else
187 nrec_nl = int(ncvarrecs(ivartype)/Ny)
188 endif
189 do irec = 1, nrec_nl
190 c And now back-calculate what iobcs should be.
191 do j=1,sny
192 iobcs= mod((irec-1)*sny+j-1,nobcs)+1
193
194 read(cunit) filencvarindex(ivartype)
195 if (filencvarindex(ivartype) .NE. ncvarindex(ivartype))
196 & then
197 print *, 'ctrl-set_unpack:xz:WARNING: wrong ncvarindex ',
198 & filencvarindex(ivartype), ncvarindex(ivartype)
199 STOP 'in S/R ctrl_set_unpack_xz'
200 endif
201 read(cunit) filej
202 read(cunit) filei
203 do k = 1, Nr
204 irectrue = (irec-1)*nobcs*nr + (iobcs-1)*nr + k
205 cbuffindex = nwetglobal(k,iobcs)
206 if ( cbuffindex .gt. 0 ) then
207 read(cunit) filencbuffindex
208 if (filencbuffindex .NE. cbuffindex) then
209 print *, 'WARNING: wrong cbuffindex ',
210 & filencbuffindex, cbuffindex
211 STOP 'in S/R ctrl_set_unpack_xz'
212 endif
213 read(cunit) filek
214 if (filek .NE. k) then
215 print *, 'WARNING: wrong k ',
216 & filek, k
217 STOP 'in S/R ctrl_set_unpack_xz'
218 endif
219 read(cunit) (cbuff(ii), ii=1,cbuffindex)
220 endif
221
222 cbuffindex = 0
223 jj=mod((j-1)*nr+k-1,sny)+1
224 kk=int((j-1)*nr+k-1)/sny+1
225 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 globfld3d(i,bi,ip,jj,bj,jp,kk) =
233 & cbuff(cbuffindex)
234 cph(
235 globfldtmp2(i,bi,ip,bj,jp) =
236 & cbuff(cbuffindex)
237 cph)
238 #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
239 globfld3d(i,bi,ip,jj,bj,jp,kk) =
240 & globfld3d(i,bi,ip,jj,bj,jp,kk)/
241 # ifdef CTRL_UNPACK_PRECISE
242 & sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
243 # else
244 & sqrt(weightfld(k,iobcs))
245 # endif
246 #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
247 else
248 globfld3d(i,bi,ip,jj,bj,jp,kk) = 0. _d 0
249 endif
250 cph(
251 globfldtmp3(i,bi,ip,bj,jp) =
252 & globfld3d(i,bi,ip,jj,bj,jp,kk)
253 cph)
254 enddo
255 enddo
256 enddo
257 enddo
258 enddo
259 c
260 if ( doPackDiag ) then
261 write(cunit2,rec=irectrue) globfldtmp2
262 write(cunit3,rec=irectrue) globfldtmp3
263 endif
264 c
265 c -- end of k loop --
266 enddo
267 c -- end of j loop --
268 enddo
269
270 call MDSWRITEFIELD_3D_GL( fname, ctrlprec, 'RL',
271 & Nr, globfld3d, irec,
272 & optimcycle, mythid)
273
274 c -- end of irec loop --
275 enddo
276
277 do irec = nrec_nl*ny+1, ncvarrecs(ivartype)
278 c And now back-calculate what iobcs should be.
279 iobcs= mod(irec-1,nobcs)+1
280
281 read(cunit) filencvarindex(ivartype)
282 if (filencvarindex(ivartype) .NE. ncvarindex(ivartype))
283 & then
284 print *, 'ctrl-set_unpack:xz:WARNING: wrong ncvarindex ',
285 & filencvarindex(ivartype), ncvarindex(ivartype)
286 STOP 'in S/R ctrl_set_unpack_xz'
287 endif
288 read(cunit) filej
289 read(cunit) filei
290 do k = 1, Nr
291 irectrue = (irec-1)*nobcs*nr + (iobcs-1)*nr + k
292 cbuffindex = nwetglobal(k,iobcs)
293 if ( cbuffindex .gt. 0 ) then
294 read(cunit) filencbuffindex
295 if (filencbuffindex .NE. cbuffindex) then
296 print *, 'WARNING: wrong cbuffindex ',
297 & filencbuffindex, cbuffindex
298 STOP 'in S/R ctrl_set_unpack_xz'
299 endif
300 read(cunit) filek
301 if (filek .NE. k) then
302 print *, 'WARNING: wrong k ',
303 & filek, k
304 STOP 'in S/R ctrl_set_unpack_xz'
305 endif
306 read(cunit) (cbuff(ii), ii=1,cbuffindex)
307 endif
308
309 cbuffindex = 0
310 do jp = 1,nPy
311 do bj = jtlo,jthi
312 do ip = 1,nPx
313 do bi = itlo,ithi
314 do i = imin,imax
315 if ( globmskxz(i,bi,ip,bj,jp,k,iobcs) .ne. 0. ) then
316 cbuffindex = cbuffindex + 1
317 globfldxz(i,bi,ip,bj,jp,k) = cbuff(cbuffindex)
318 cph(
319 globfldtmp2(i,bi,ip,bj,jp) = cbuff(cbuffindex)
320 cph)
321 #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
322 globfldxz(i,bi,ip,bj,jp,k) =
323 & globfldxz(i,bi,ip,bj,jp,k)/
324 # ifdef CTRL_UNPACK_PRECISE
325 & sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
326 # else
327 & sqrt(weightfld(k,iobcs))
328 # endif
329 #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
330 else
331 globfldxz(i,bi,ip,bj,jp,k) = 0. _d 0
332 endif
333 cph(
334 globfldtmp3(i,bi,ip,bj,jp) =
335 & globfldxz(i,bi,ip,bj,jp,k)
336 cph)
337 enddo
338 enddo
339 enddo
340 enddo
341 enddo
342 c
343 if ( doPackDiag ) then
344 write(cunit2,rec=irectrue) globfldtmp2
345 write(cunit3,rec=irectrue) globfldtmp3
346 endif
347 c
348 c -- end of k loop --
349 enddo
350
351 call MDSWRITEFIELD_XZ_GL( fname, ctrlprec, 'RL',
352 & Nr, globfldxz, irec,
353 & optimcycle, mythid)
354
355 c -- end of irec loop --
356 enddo
357
358 _END_MASTER( mythid )
359
360 #endif
361
362 return
363 end
364
365
366
367
368

  ViewVC Help
Powered by ViewVC 1.1.22