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

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

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


Revision 1.4 - (show annotations) (download)
Thu Jul 24 22:00:18 2003 UTC (20 years, 9 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint51k_post, checkpoint51l_post, checkpoint51f_post, checkpoint51j_post, checkpoint51l_pre, checkpoint51h_pre, branchpoint-genmake2, checkpoint51i_post, checkpoint51i_pre, checkpoint51e_post, checkpoint51f_pre, checkpoint51g_post, checkpoint51m_post
Branch point for: branch-genmake2, tg2-branch
Changes since 1.3: +8 -5 lines
bug fixes for 3d packing and I/O of sliced (xz/yz) fields
to increase I/O performance.

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

  ViewVC Help
Powered by ViewVC 1.1.22