/[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.5 - (show annotations) (download)
Thu Oct 23 04:41:40 2003 UTC (20 years, 7 months ago) by edhill
Branch: MAIN
CVS Tags: checkpoint51o_pre, checkpoint51n_post, checkpoint51n_pre, checkpoint51o_post, checkpoint51p_post
Branch point for: checkpoint51n_branch
Changes since 1.4: +4 -0 lines
 o added the [#include "AD_CONFIG.h"] statement to all files that need
   it for adjoint/tl #defines
 o re-worked the build logic in genmake2 to support AD_CONFIG.h
 o removed tools/genmake since it no longer works

1 C
2 C $Header: $
3 C $Name: $
4
5 #include "AD_CONFIG.h"
6 #include "CTRL_CPPOPTIONS.h"
7
8
9 subroutine ctrl_set_pack_yz(
10 & cunit, ivartype, fname, masktype, weighttype,
11 & weightfld, lxxadxx, mythid)
12
13 c ==================================================================
14 c SUBROUTINE ctrl_set_pack_yz
15 c ==================================================================
16 c
17 c o Compress the control vector such that only ocean points are
18 c written to file.
19 c
20 c o Open boundary packing added :
21 c gebbie@mit.edu, 18-Mar-2003
22 c
23 c changed: heimbach@mit.edu 17-Jun-2003
24 c merged Armin's changes to replace write of
25 c nr * globfld2d by 1 * globfld3d
26 c (ad hoc fix to speed up global I/O)
27 c
28 c ==================================================================
29
30 implicit none
31
32 c == global variables ==
33
34 #include "EEPARAMS.h"
35 #include "SIZE.h"
36 #include "PARAMS.h"
37 #include "GRID.h"
38
39 #include "ctrl.h"
40 #include "cost.h"
41
42 #ifdef ALLOW_ECCO_OPTIMIZATION
43 #include "optim.h"
44 #endif
45
46 c == routine arguments ==
47
48 integer cunit
49 integer ivartype
50 character*( 80) fname
51 character*( 9) masktype
52 character*( 80) weighttype
53 _RL weightfld( nr,nobcs )
54 logical lxxadxx
55 integer mythid
56
57 c == local variables ==
58
59 #ifndef ALLOW_ECCO_OPTIMIZATION
60 integer optimcycle
61 #endif
62
63 integer bi,bj
64 integer ip,jp
65 integer i,j,k
66 integer ii,kk
67 integer il
68 integer irec,iobcs,nrec_nl
69 integer itlo,ithi
70 integer jtlo,jthi
71 integer jmin,jmax
72 integer imin,imax
73
74 integer cbuffindex
75 cgg(
76 integer igg
77 _RL gg
78 character*(80) weightname
79 cgg)
80 _RL cbuff ( nsx*npx*sny*nsy*npy )
81 _RL globfldyz ( nsx,npx,sny,nsy,npy,nr )
82 _RL globfld3d ( snx,nsx,npx,sny,nsy,npy,nr )
83 _RL globmskyz ( nsx,npx,sny,nsy,npy,nr,nobcs )
84 #ifdef CTRL_PACK_PRECISE
85 _RL weightfldyz( nsx,npx,sny,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 j = jmin,jmax
113 do ip = 1,nPx
114 do bi = itlo,ithi
115 globfldyz(bi,ip,j,bj,jp,k) = 0. _d 0
116 do iobcs=1,nobcs
117 globmskyz(bi,ip,j,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_YZ_GL(
147 & masktype, ctrlprec, 'RL',
148 & Nr, globmskyz(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_YZ_GL(
154 & weightname, ctrlprec, 'RL',
155 & Nr, weightfldyz(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 j = jmin,jmax
163 do ip = 1,nPx
164 do bi = itlo,ithi
165 weightfldyz(bi,ip,j,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)/snx)
176 do irec = 1, nrec_nl
177 cgg do iobcs = 1, nobcs
178 cgg Need to solve for what iobcs would have been.
179
180 call MDSREADFIELD_3D_GL( fname, ctrlprec, 'RL',
181 & nr, globfld3D, irec, mythid)
182
183 do i=1,snx
184 iobcs= mod((irec-1)*snx+i-1,nobcs)+1
185
186 CGG One special exception: barotropic velocity should be nondimensionalized
187 cgg differently. Probably introduce new variable.
188 if (iobcs .eq. 3 .or. iobcs .eq. 4) then
189 k = 1
190 do jp = 1,nPy
191 do bj = jtlo,jthi
192 do j = jmin,jmax
193 do ip = 1,nPx
194 do bi = itlo,ithi
195 #ifdef NO_CONTROL_BAROTROPIC_VELOCITY
196 if (.not. lxxadxx) then
197 cgg Get rid of any sensitivity to barotropic velocity.
198 globfld3d(i,bi,ip,j,bj,jp,k) = 0.
199 endif
200 #endif
201 enddo
202 enddo
203 enddo
204 enddo
205 enddo
206 endif
207
208 write(cunit) ncvarindex(ivartype)
209 write(cunit) 1
210 write(cunit) 1
211 do k = 1,nr
212 cbuffindex = 0
213 do jp = 1,nPy
214 do bj = jtlo,jthi
215 do ip = 1,nPx
216 do bi = itlo,ithi
217 do j = jmin,jmax
218 if (globmskyz(bi,ip,j,bj,jp,k,iobcs) .ne. 0. ) then
219 cbuffindex = cbuffindex + 1
220 ii=mod((i-1)*nr+k-1,snx)+1
221 kk=int((i-1)*nr+k-1)/snx+1
222 #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
223 if (lxxadxx) then
224 cbuff(cbuffindex) =
225 & globfld3d(ii,bi,ip,j,bj,jp,kk) *
226 #ifdef CTRL_PACK_PRECISE
227 & sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs))
228 #else
229 & sqrt(weightfld(k,iobcs))
230 #endif
231 else
232 cbuff(cbuffindex) =
233 & globfld3d(ii,bi,ip,j,bj,jp,kk) /
234 #ifdef CTRL_PACK_PRECISE
235 & sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs))
236 #else
237 & sqrt(weightfld(k,iobcs))
238 #endif
239 endif
240 #else /* ALLOW_NONDIMENSIONAL_CONTROL_IO undef */
241 cbuff(cbuffindex) = globfld3d(ii,bi,ip,j,bj,jp,kk)
242 #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
243 endif
244 enddo
245 enddo
246 enddo
247 enddo
248 enddo
249 c --> check cbuffindex.
250 if ( cbuffindex .gt. 0) then
251 write(cunit) cbuffindex
252 write(cunit) k
253 write(cunit) (cbuff(ii), ii=1,cbuffindex)
254 endif
255 c
256 c -- end of k loop --
257 enddo
258 c -- end of iobcs loop --
259 cgg enddo
260 c -- end of i loop --
261 enddo
262 c -- end of irec loop --
263 enddo
264
265 do irec = nrec_nl*snx+1, ncvarrecs(ivartype)
266 cgg do iobcs = 1, nobcs
267 cgg Need to solve for what iobcs would have been.
268 iobcs= mod(irec-1,nobcs)+1
269
270 call MDSREADFIELD_YZ_GL( fname, ctrlprec, 'RL',
271 & nr, globfldyz, irec, mythid)
272
273 CGG One special exception: barotropic velocity should be nondimensionalized
274 cgg differently. Probably introduce new variable.
275 if (iobcs .eq. 3 .or. iobcs .eq. 4) then
276 k = 1
277 do jp = 1,nPy
278 do bj = jtlo,jthi
279 do j = jmin,jmax
280 do ip = 1,nPx
281 do bi = itlo,ithi
282 #ifdef NO_CONTROL_BAROTROPIC_VELOCITY
283 if (.not. lxxadxx) then
284 cgg Get rid of any sensitivity to barotropic velocity.
285 globfldyz(bi,ip,j,bj,jp,k) = 0.
286 endif
287 #endif
288 enddo
289 enddo
290 enddo
291 enddo
292 enddo
293 endif
294
295 write(cunit) ncvarindex(ivartype)
296 write(cunit) 1
297 write(cunit) 1
298 do k = 1,nr
299 cbuffindex = 0
300 do jp = 1,nPy
301 do bj = jtlo,jthi
302 do ip = 1,nPx
303 do bi = itlo,ithi
304 do j = jmin,jmax
305 if (globmskyz(bi,ip,j,bj,jp,k,iobcs) .ne. 0. ) then
306 cbuffindex = cbuffindex + 1
307 #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
308 if (lxxadxx) then
309 cbuff(cbuffindex) =
310 & globfldyz(bi,ip,j,bj,jp,k) *
311 #ifdef CTRL_PACK_PRECISE
312 & sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs))
313 #else
314 & sqrt(weightfld(k,iobcs))
315 #endif
316 else
317 cbuff(cbuffindex) =
318 & globfldyz(bi,ip,j,bj,jp,k) /
319 #ifdef CTRL_PACK_PRECISE
320 & sqrt(weightfldyz(bi,ip,j,bj,jp,k,iobcs))
321 #else
322 & sqrt(weightfld(k,iobcs))
323 #endif
324 endif
325 #else /* ALLOW_NONDIMENSIONAL_CONTROL_IO undef */
326 cbuff(cbuffindex) = globfldyz(bi,ip,j,bj,jp,k)
327 #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
328 endif
329 enddo
330 enddo
331 enddo
332 enddo
333 enddo
334 c --> check cbuffindex.
335 if ( cbuffindex .gt. 0) then
336 write(cunit) cbuffindex
337 write(cunit) k
338 write(cunit) (cbuff(ii), ii=1,cbuffindex)
339 endif
340 c
341 c -- end of k loop --
342 enddo
343 c -- end of iobcs loop --
344 cgg enddo
345 c -- end of irec loop --
346 enddo
347
348 _END_MASTER( mythid )
349
350 return
351 end
352

  ViewVC Help
Powered by ViewVC 1.1.22