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

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

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


Revision 1.19 - (show annotations) (download)
Thu Apr 7 23:38:43 2005 UTC (19 years, 1 month ago) by heimbach
Branch: MAIN
Changes since 1.18: +74 -26 lines
o separate masks used for ctrl_pack/unpack 'from write_grid' output
  (suggested by G. Forget)
o added new control variables
  * init. uVel, vVel, etanN
  * lambda[Theta,Salt]ClimRelax

1 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_pack.F,v 1.18 2005/02/28 17:29:38 heimbach Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CTRL_CPPOPTIONS.h"
6
7 subroutine ctrl_pack( first, mythid )
8
9 c ==================================================================
10 c SUBROUTINE ctrl_pack
11 c ==================================================================
12 c
13 c o Compress the control vector such that only ocean points are
14 c written to file.
15 c
16 c started: Christian Eckert eckert@mit.edu 10-Mar=2000
17 c
18 c changed: Patrick Heimbach heimbach@mit.edu 06-Jun-2000
19 c - Transferred some filename declarations
20 c from here to namelist in ctrl_init
21 c
22 c Patrick Heimbach heimbach@mit.edu 16-Jun-2000
23 c - single file name convention with or without
24 c ALLOW_ECCO_OPTIMIZATION
25 c
26 c G. Gebbie, added open boundary control packing,
27 c gebbie@mit.edu 18 -Mar- 2003
28 c
29 c heimbach@mit.edu totally restructured 28-Oct-2003
30 c
31 c ==================================================================
32 c SUBROUTINE ctrl_pack
33 c ==================================================================
34
35 implicit none
36
37 c == global variables ==
38
39 #include "EEPARAMS.h"
40 #include "SIZE.h"
41 #include "PARAMS.h"
42 #include "GRID.h"
43
44 #include "ctrl.h"
45 #include "optim.h"
46
47 #ifdef ALLOW_COST
48 # include "cost.h"
49 #endif
50 #ifdef ALLOW_ECCO
51 # include "ecco_cost.h"
52 #else
53 # include "ctrl_weights.h"
54 #endif
55
56 c == routine arguments ==
57
58 logical first
59 integer mythid
60
61 #ifndef EXCLUDE_CTRL_PACK
62 c == local variables ==
63
64 _RL fcloc
65
66 integer i, j, k
67 integer ii
68 integer il
69 integer irec
70 integer ig,jg
71 integer ivartype
72 integer iobcs
73
74 logical doglobalread
75 logical ladinit
76 integer cbuffindex
77 logical lxxadxx
78
79 integer cunit
80 integer ictrlgrad
81
82 character*(128) cfile
83 character*( 80) weighttype
84
85 c == external ==
86
87 integer ilnblnk
88 external ilnblnk
89
90 c == end of interface ==
91
92 #ifndef ALLOW_ECCO_OPTIMIZATION
93 fmin = 0. _d 0
94 #endif
95
96 c-- Tiled files are used.
97 doglobalread = .false.
98
99 c-- Initialise adjoint variables on active files.
100 ladinit = .false.
101
102 c-- Initialise global buffer index
103 nbuffglobal = 0
104
105 c-- Assign file names.
106
107 call ctrl_set_fname(xx_theta_file, fname_theta, mythid)
108 call ctrl_set_fname(xx_salt_file, fname_salt, mythid)
109 call ctrl_set_fname(xx_hflux_file, fname_hflux, mythid)
110 call ctrl_set_fname(xx_sflux_file, fname_sflux, mythid)
111 call ctrl_set_fname(xx_tauu_file, fname_tauu, mythid)
112 call ctrl_set_fname(xx_tauv_file, fname_tauv, mythid)
113 call ctrl_set_fname(xx_atemp_file, fname_atemp, mythid)
114 call ctrl_set_fname(xx_aqh_file, fname_aqh, mythid)
115 call ctrl_set_fname(xx_uwind_file, fname_uwind, mythid)
116 call ctrl_set_fname(xx_vwind_file, fname_vwind, mythid)
117 call ctrl_set_fname(xx_obcsn_file, fname_obcsn, mythid)
118 call ctrl_set_fname(xx_obcss_file, fname_obcss, mythid)
119 call ctrl_set_fname(xx_obcsw_file, fname_obcsw, mythid)
120 call ctrl_set_fname(xx_obcse_file, fname_obcse, mythid)
121 call ctrl_set_fname(xx_diffkr_file, fname_diffkr, mythid)
122 call ctrl_set_fname(xx_kapgm_file, fname_kapgm, mythid)
123 call ctrl_set_fname(xx_tr1_file, fname_tr1, mythid)
124 call ctrl_set_fname(xx_sst_file, fname_sst, mythid)
125 call ctrl_set_fname(xx_sss_file, fname_sss, mythid)
126 call ctrl_set_fname(xx_hfacc_file, fname_hfacc, mythid)
127 call ctrl_set_fname(xx_efluxy_file, fname_efluxy, mythid)
128 call ctrl_set_fname(xx_efluxp_file, fname_efluxp, mythid)
129 call ctrl_set_fname(xx_bottomdrag_file, fname_bottomdrag, mythid)
130 call ctrl_set_fname(xx_edtaux_file, fname_edtaux, mythid)
131 call ctrl_set_fname(xx_edtauy_file, fname_edtauy, mythid)
132 call ctrl_set_fname(xx_uvel_file, fname_uvel, mythid)
133 call ctrl_set_fname(xx_vvel_file, fname_vvel, mythid)
134 call ctrl_set_fname(xx_etan_file, fname_etan, mythid)
135 call ctrl_set_fname(xx_relaxsst_file, fname_relaxsst, mythid)
136 call ctrl_set_fname(xx_relaxsss_file, fname_relaxsss, mythid)
137
138 c-- Only the master thread will do I/O.
139 _BEGIN_MASTER( mythid )
140
141 if ( first ) then
142 c >>> Initialise control vector for optimcycle=0 <<<
143 lxxadxx = .TRUE.
144 ictrlgrad = 1
145 fcloc = fmin
146 write(cfile(1:128),'(4a,i4.4)')
147 & ctrlname(1:9),'_',yctrlid(1:10),
148 & yctrlpospack, optimcycle
149 print *, 'ph-pack: packing ', ctrlname(1:9)
150 else
151 c >>> Write gradient vector <<<
152 lxxadxx = .FALSE.
153 ictrlgrad = 2
154 fcloc = fc
155 write(cfile(1:128),'(4a,i4.4)')
156 & costname(1:9),'_',yctrlid(1:10),
157 & yctrlpospack, optimcycle
158 print *, 'ph-pack: packing ', costname(1:9)
159 endif
160
161 call mdsfindunit( cunit, mythid )
162 open( cunit, file = cfile,
163 & status = 'unknown',
164 & form = 'unformatted',
165 & access = 'sequential' )
166
167 c-- Header information.
168 write(cunit) nvartype
169 write(cunit) nvarlength
170 write(cunit) yctrlid
171 write(cunit) optimCycle
172 write(cunit) fc
173 C place holder of obsolete variable iG
174 write(cunit) 1
175 C place holder of obsolete variable jG
176 write(cunit) 1
177 write(cunit) nsx
178 write(cunit) nsy
179 write(cunit) (nWetcGlobal(k), k=1,nr)
180 write(cunit) (nWetsGlobal(k), k=1,nr)
181 write(cunit) (nWetwGlobal(k), k=1,nr)
182 #ifdef ALLOW_CTRL_WETV
183 write(cunit) (nWetvGlobal(k), k=1,nr)
184 #endif
185
186 #ifdef ALLOW_OBCSN_CONTROL
187 write(cunit) ((nWetobcsnGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
188 #endif
189 #ifdef ALLOW_OBCSS_CONTROL
190 write(cunit) ((nWetobcssGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
191 #endif
192 #ifdef ALLOW_OBCSW_CONTROL
193 write(cunit) ((nWetobcswGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
194 #endif
195 #ifdef ALLOW_OBCSE_CONTROL
196 write(cunit) ((nWetobcseGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
197 #endif
198 write(cunit) (ncvarindex(i), i=1,maxcvars)
199 write(cunit) (ncvarrecs(i), i=1,maxcvars)
200 write(cunit) (ncvarxmax(i), i=1,maxcvars)
201 write(cunit) (ncvarymax(i), i=1,maxcvars)
202 write(cunit) (ncvarnrmax(i), i=1,maxcvars)
203 write(cunit) (ncvargrd(i), i=1,maxcvars)
204 write(cunit)
205
206 #ifdef ALLOW_THETA0_CONTROL
207 ivartype = 1
208 write(weighttype(1:80),'(80a)') ' '
209 write(weighttype(1:80),'(a)') "wtheta"
210 call ctrl_set_pack_xyz(
211 & cunit, ivartype, fname_theta(ictrlgrad), "maskCtrlC",
212 & weighttype, wtheta, lxxadxx, mythid)
213 #endif
214
215 #ifdef ALLOW_SALT0_CONTROL
216 ivartype = 2
217 write(weighttype(1:80),'(80a)') ' '
218 write(weighttype(1:80),'(a)') "wsalt"
219 call ctrl_set_pack_xyz(
220 & cunit, ivartype, fname_salt(ictrlgrad), "maskCtrlC",
221 & weighttype, wsalt, lxxadxx, mythid)
222 #endif
223
224 #if (defined (ALLOW_HFLUX_CONTROL) || \
225 defined (ALLOW_HFLUX0_CONTROL))
226 ivartype = 3
227 write(weighttype(1:80),'(80a)') ' '
228 write(weighttype(1:80),'(a)') "whflux"
229 call ctrl_set_pack_xy(
230 & cunit, ivartype, fname_hflux(ictrlgrad), "maskCtrlC",
231 & weighttype, lxxadxx, mythid)
232 #endif
233
234 #if (defined (ALLOW_SFLUX_CONTROL) || \
235 defined (ALLOW_SFLUX0_CONTROL))
236 ivartype = 4
237 write(weighttype(1:80),'(80a)') ' '
238 write(weighttype(1:80),'(a)') "wsflux"
239 call ctrl_set_pack_xy(
240 & cunit, ivartype, fname_sflux(ictrlgrad), "maskCtrlC",
241 & weighttype, lxxadxx, mythid)
242 #endif
243
244 #if (defined (ALLOW_USTRESS_CONTROL) || \
245 defined (ALLOW_TAUU0_CONTROL))
246 ivartype = 5
247 write(weighttype(1:80),'(80a)') ' '
248 write(weighttype(1:80),'(a)') "wtauu"
249 call ctrl_set_pack_xy(
250 & cunit, ivartype, fname_tauu(ictrlgrad), "maskCtrlW",
251 & weighttype, lxxadxx, mythid)
252 #endif
253
254 #if (defined (ALLOW_VSTRESS_CONTROL) || \
255 defined (ALLOW_TAUV0_CONTROL))
256 ivartype = 6
257 write(weighttype(1:80),'(80a)') ' '
258 write(weighttype(1:80),'(a)') "wtauv"
259 call ctrl_set_pack_xy(
260 & cunit, ivartype, fname_tauv(ictrlgrad), "maskCtrlS",
261 & weighttype, lxxadxx, mythid)
262 #endif
263
264 #ifdef ALLOW_ATEMP_CONTROL
265 ivartype = 7
266 write(weighttype(1:80),'(80a)') ' '
267 write(weighttype(1:80),'(a)') "watemp"
268 call ctrl_set_pack_xy(
269 & cunit, ivartype, fname_atemp(ictrlgrad), "maskCtrlC",
270 & weighttype, lxxadxx, mythid)
271 #endif
272
273 #ifdef ALLOW_AQH_CONTROL
274 ivartype = 8
275 write(weighttype(1:80),'(80a)') ' '
276 write(weighttype(1:80),'(a)') "waqh"
277 call ctrl_set_pack_xy(
278 & cunit, ivartype, fname_aqh(ictrlgrad), "maskCtrlC",
279 & weighttype, lxxadxx, mythid)
280 #endif
281
282 #ifdef ALLOW_UWIND_CONTROL
283 ivartype = 9
284 write(weighttype(1:80),'(80a)') ' '
285 write(weighttype(1:80),'(a)') "wuwind"
286 call ctrl_set_pack_xy(
287 & cunit, ivartype, fname_uwind(ictrlgrad), "maskCtrlW",
288 & weighttype, lxxadxx, mythid)
289 #endif
290
291 #ifdef ALLOW_VWIND_CONTROL
292 ivartype = 10
293 write(weighttype(1:80),'(80a)') ' '
294 write(weighttype(1:80),'(a)') "wvwind"
295 call ctrl_set_pack_xy(
296 & cunit, ivartype, fname_vwind(ictrlgrad), "maskCtrlS",
297 & weighttype, lxxadxx, mythid)
298 #endif
299
300 #ifdef ALLOW_OBCSN_CONTROL
301 ivartype = 11
302 write(weighttype(1:80),'(80a)') ' '
303 write(weighttype(1:80),'(a)') "wobcsn"
304 call ctrl_set_pack_xz(
305 & cunit, ivartype, fname_obcsn(ictrlgrad), "maskobcsn",
306 & weighttype, wobcsn, lxxadxx, mythid)
307 #endif
308
309 #ifdef ALLOW_OBCSS_CONTROL
310 ivartype = 12
311 write(weighttype(1:80),'(80a)') ' '
312 write(weighttype(1:80),'(a)') "wobcss"
313 call ctrl_set_pack_xz(
314 & cunit, ivartype, fname_obcss(ictrlgrad), "maskobcss",
315 & weighttype, wobcss, lxxadxx, mythid)
316 #endif
317
318 #ifdef ALLOW_OBCSW_CONTROL
319 ivartype = 13
320 write(weighttype(1:80),'(80a)') ' '
321 write(weighttype(1:80),'(a)') "wobcsw"
322 call ctrl_set_pack_yz(
323 & cunit, ivartype, fname_obcsw(ictrlgrad), "maskobcsw",
324 & weighttype, wobcsw, lxxadxx, mythid)
325 #endif
326
327 #ifdef ALLOW_OBCSE_CONTROL
328 ivartype = 14
329 write(weighttype(1:80),'(80a)') ' '
330 write(weighttype(1:80),'(a)') "wobcse"
331 call ctrl_set_pack_yz(
332 & cunit, ivartype, fname_obcse(ictrlgrad), "maskobcse",
333 & weighttype, wobcse, lxxadxx, mythid)
334 #endif
335
336 #ifdef ALLOW_DIFFKR_CONTROL
337 ivartype = 15
338 write(weighttype(1:80),'(80a)') ' '
339 write(weighttype(1:80),'(a)') "wdiffkr"
340 call ctrl_set_pack_xyz(
341 & cunit, ivartype, fname_diffkr(ictrlgrad), "maskCtrlC",
342 & weighttype, wunit, lxxadxx, mythid)
343 #endif
344
345 #ifdef ALLOW_KAPGM_CONTROL
346 ivartype = 16
347 write(weighttype(1:80),'(80a)') ' '
348 write(weighttype(1:80),'(a)') "wkapgm"
349 call ctrl_set_pack_xyz(
350 & cunit, ivartype, fname_kapgm(ictrlgrad), "maskCtrlC",
351 & weighttype, wunit, lxxadxx, mythid)
352 #endif
353
354 #ifdef ALLOW_TR10_CONTROL
355 ivartype = 17
356 write(weighttype(1:80),'(80a)') ' '
357 write(weighttype(1:80),'(a)') "wtr1"
358 call ctrl_set_pack_xyz(
359 & cunit, ivartype, fname_tr1(ictrlgrad), "maskCtrlC",
360 & weighttype, wunit, lxxadxx, mythid)
361 #endif
362
363 #ifdef ALLOW_SST0_CONTROL
364 ivartype = 18
365 write(weighttype(1:80),'(80a)') ' '
366 write(weighttype(1:80),'(a)') "wsst0"
367 call ctrl_set_pack_xy(
368 & cunit, ivartype, fname_sst(ictrlgrad), "maskCtrlC",
369 & weighttype, lxxadxx, mythid)
370 #endif
371
372 #ifdef ALLOW_SSS0_CONTROL
373 ivartype = 19
374 write(weighttype(1:80),'(80a)') ' '
375 write(weighttype(1:80),'(a)') "wsss0"
376 call ctrl_set_pack_xy(
377 & cunit, ivartype, fname_sss(ictrlgrad), "maskCtrlC",
378 & weighttype, lxxadxx, mythid)
379 #endif
380
381 #ifdef ALLOW_HFACC_CONTROL
382 ivartype = 20
383 write(weighttype(1:80),'(80a)') ' '
384 write(weighttype(1:80),'(a)') "whfacc"
385 # ifdef ALLOW_HFACC3D_CONTROL
386 call ctrl_set_pack_xyz(
387 & cunit, ivartype, fname_hfacc(ictrlgrad), "maskCtrlC",
388 & weighttype, wunit, lxxadxx, mythid)
389 # else
390 call ctrl_set_pack_xy(
391 & cunit, ivartype, fname_hfacc(ictrlgrad), "maskCtrlC",
392 & weighttype, lxxadxx, mythid)
393 # endif
394 #endif
395
396 #ifdef ALLOW_EFLUXY0_CONTROL
397 ivartype = 21
398 write(weighttype(1:80),'(80a)') ' '
399 write(weighttype(1:80),'(a)') "wefluxy0"
400 call ctrl_set_pack_xyz(
401 & cunit, ivartype, fname_efluxy(ictrlgrad), "maskCtrlS",
402 & weighttype, wunit, lxxadxx, mythid)
403 #endif
404
405 #ifdef ALLOW_EFLUXP0_CONTROL
406 ivartype = 22
407 write(weighttype(1:80),'(80a)') ' '
408 write(weighttype(1:80),'(a)') "wefluxp0"
409 call ctrl_set_pack_xyz(
410 & cunit, ivartype, fname_efluxp(ictrlgrad), "maskhFacV",
411 & weighttype, wunit, lxxadxx, mythid)
412 #endif
413
414 #ifdef ALLOW_BOTTOMDRAG_CONTROL
415 ivartype = 23
416 write(weighttype(1:80),'(80a)') ' '
417 write(weighttype(1:80),'(a)') "wbottomdrag"
418 call ctrl_set_pack_xy(
419 & cunit, ivartype, fname_bottomdrag(ictrlgrad), "maskCtrlC",
420 & weighttype, lxxadxx, mythid)
421 #endif
422
423 #ifdef ALLOW_EDTAUX_CONTROL
424 ivartype = 25
425 write(weighttype(1:80),'(80a)') ' '
426 write(weighttype(1:80),'(a)') "wedtaux"
427 call ctrl_set_pack_xyz(
428 & cunit, ivartype, fname_edtaux(ictrlgrad), "maskCtrlW",
429 & weighttype, wunit, lxxadxx, mythid)
430 #endif
431
432 #ifdef ALLOW_EDTAUY_CONTROL
433 ivartype = 26
434 write(weighttype(1:80),'(80a)') ' '
435 write(weighttype(1:80),'(a)') "wedtauy"
436 call ctrl_set_pack_xyz(
437 & cunit, ivartype, fname_edtauy(ictrlgrad), "maskCtrlS",
438 & weighttype, wunit, lxxadxx, mythid)
439 #endif
440
441 #ifdef ALLOW_UVEL0_CONTROL
442 ivartype = 27
443 write(weighttype(1:80),'(80a)') ' '
444 write(weighttype(1:80),'(a)') "wuvel"
445 call ctrl_set_pack_xyz(
446 & cunit, ivartype, fname_uvel(ictrlgrad), "maskCtrlW",
447 & weighttype, wunit, lxxadxx, mythid)
448 #endif
449
450 #ifdef ALLOW_VVEL0_CONTROL
451 ivartype = 28
452 write(weighttype(1:80),'(80a)') ' '
453 write(weighttype(1:80),'(a)') "wvvel"
454 call ctrl_set_pack_xyz(
455 & cunit, ivartype, fname_vvel(ictrlgrad), "maskCtrlS",
456 & weighttype, wunit, lxxadxx, mythid)
457 #endif
458
459 #ifdef ALLOW_ETAN0_CONTROL
460 ivartype = 29
461 write(weighttype(1:80),'(80a)') ' '
462 write(weighttype(1:80),'(a)') "wetan"
463 call ctrl_set_pack_xy(
464 & cunit, ivartype, fname_etan(ictrlgrad), "maskCtrlC",
465 & weighttype, lxxadxx, mythid)
466 #endif
467
468 #ifdef ALLOW_RELAXSST_CONTROL
469 ivartype = 30
470 write(weighttype(1:80),'(80a)') ' '
471 write(weighttype(1:80),'(a)') "wrelaxsst"
472 call ctrl_set_pack_xy(
473 & cunit, ivartype, fname_relaxsst(ictrlgrad), "maskCtrlC",
474 & weighttype, lxxadxx, mythid)
475 #endif
476
477 #ifdef ALLOW_RELAXSSS_CONTROL
478 ivartype = 31
479 write(weighttype(1:80),'(80a)') ' '
480 write(weighttype(1:80),'(a)') "wrelaxsss"
481 call ctrl_set_pack_xy(
482 & cunit, ivartype, fname_relaxsss(ictrlgrad), "maskCtrlC",
483 & weighttype, lxxadxx, mythid)
484 #endif
485
486 close ( cunit )
487
488 _END_MASTER( mythid )
489
490 #endif /* EXCLUDE_CTRL_PACK */
491
492 return
493 end
494

  ViewVC Help
Powered by ViewVC 1.1.22