/[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.6 - (show annotations) (download)
Thu Oct 30 19:09:06 2003 UTC (20 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint51t_post, checkpoint51s_post, checkpoint51q_post, checkpoint51r_post
Branch point for: branch-nonh
Changes since 1.5: +2 -3 lines
ctrl package totally restructured
o pack/unpack now optional and decoupled from
  xx_/adxx_ I/O
o ctrl_pack/unpack cleaned
  (new routines ctrl_init_ctrlvar.F, pkg/ctrl/ctrl_init_wet.F)
o confined inclusion of AD_CONFIG.h to where necessary.

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

  ViewVC Help
Powered by ViewVC 1.1.22