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

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

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


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

  ViewVC Help
Powered by ViewVC 1.1.22