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

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