/[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.23 - (show annotations) (download)
Sat Aug 6 11:02:01 2005 UTC (18 years, 10 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57q_post
Changes since 1.22: +11 -1 lines
Adding swdown control.

1 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_pack.F,v 1.22 2005/07/28 19:51:22 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_precip_file, fname_precip, mythid)
116 call ctrl_set_fname(xx_swflux_file, fname_swflux, mythid)
117 call ctrl_set_fname(xx_swdown_file, fname_swdown, mythid)
118 call ctrl_set_fname(xx_uwind_file, fname_uwind, mythid)
119 call ctrl_set_fname(xx_vwind_file, fname_vwind, mythid)
120 call ctrl_set_fname(xx_obcsn_file, fname_obcsn, mythid)
121 call ctrl_set_fname(xx_obcss_file, fname_obcss, mythid)
122 call ctrl_set_fname(xx_obcsw_file, fname_obcsw, mythid)
123 call ctrl_set_fname(xx_obcse_file, fname_obcse, mythid)
124 call ctrl_set_fname(xx_diffkr_file, fname_diffkr, mythid)
125 call ctrl_set_fname(xx_kapgm_file, fname_kapgm, mythid)
126 call ctrl_set_fname(xx_tr1_file, fname_tr1, mythid)
127 call ctrl_set_fname(xx_sst_file, fname_sst, mythid)
128 call ctrl_set_fname(xx_sss_file, fname_sss, mythid)
129 call ctrl_set_fname(xx_hfacc_file, fname_hfacc, mythid)
130 call ctrl_set_fname(xx_efluxy_file, fname_efluxy, mythid)
131 call ctrl_set_fname(xx_efluxp_file, fname_efluxp, mythid)
132 call ctrl_set_fname(xx_bottomdrag_file, fname_bottomdrag, mythid)
133 call ctrl_set_fname(xx_edtaux_file, fname_edtaux, mythid)
134 call ctrl_set_fname(xx_edtauy_file, fname_edtauy, mythid)
135 call ctrl_set_fname(xx_uvel_file, fname_uvel, mythid)
136 call ctrl_set_fname(xx_vvel_file, fname_vvel, mythid)
137 call ctrl_set_fname(xx_etan_file, fname_etan, mythid)
138 call ctrl_set_fname(xx_relaxsst_file, fname_relaxsst, mythid)
139 call ctrl_set_fname(xx_relaxsss_file, fname_relaxsss, mythid)
140
141 c-- Only the master thread will do I/O.
142 _BEGIN_MASTER( mythid )
143
144 if ( first ) then
145 c >>> Initialise control vector for optimcycle=0 <<<
146 lxxadxx = .TRUE.
147 ictrlgrad = 1
148 fcloc = fmin
149 write(cfile(1:128),'(4a,i4.4)')
150 & ctrlname(1:9),'_',yctrlid(1:10),
151 & yctrlpospack, optimcycle
152 print *, 'ph-pack: packing ', ctrlname(1:9)
153 else
154 c >>> Write gradient vector <<<
155 lxxadxx = .FALSE.
156 ictrlgrad = 2
157 fcloc = fc
158 write(cfile(1:128),'(4a,i4.4)')
159 & costname(1:9),'_',yctrlid(1:10),
160 & yctrlpospack, optimcycle
161 print *, 'ph-pack: packing ', costname(1:9)
162 endif
163
164 call mdsfindunit( cunit, mythid )
165 open( cunit, file = cfile,
166 & status = 'unknown',
167 & form = 'unformatted',
168 & access = 'sequential' )
169
170 c-- Header information.
171 write(cunit) nvartype
172 write(cunit) nvarlength
173 write(cunit) yctrlid
174 write(cunit) optimCycle
175 write(cunit) fc
176 C place holder of obsolete variable iG
177 write(cunit) 1
178 C place holder of obsolete variable jG
179 write(cunit) 1
180 write(cunit) nsx
181 write(cunit) nsy
182 write(cunit) (nWetcGlobal(k), k=1,nr)
183 write(cunit) (nWetsGlobal(k), k=1,nr)
184 write(cunit) (nWetwGlobal(k), k=1,nr)
185 #ifdef ALLOW_CTRL_WETV
186 write(cunit) (nWetvGlobal(k), k=1,nr)
187 #endif
188
189 #ifdef ALLOW_OBCSN_CONTROL
190 write(cunit) ((nWetobcsnGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
191 #endif
192 #ifdef ALLOW_OBCSS_CONTROL
193 write(cunit) ((nWetobcssGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
194 #endif
195 #ifdef ALLOW_OBCSW_CONTROL
196 write(cunit) ((nWetobcswGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
197 #endif
198 #ifdef ALLOW_OBCSE_CONTROL
199 write(cunit) ((nWetobcseGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
200 #endif
201 write(cunit) (ncvarindex(i), i=1,maxcvars)
202 write(cunit) (ncvarrecs(i), i=1,maxcvars)
203 write(cunit) (ncvarxmax(i), i=1,maxcvars)
204 write(cunit) (ncvarymax(i), i=1,maxcvars)
205 write(cunit) (ncvarnrmax(i), i=1,maxcvars)
206 write(cunit) (ncvargrd(i), i=1,maxcvars)
207 write(cunit)
208
209 #ifdef ALLOW_THETA0_CONTROL
210 ivartype = 1
211 write(weighttype(1:80),'(80a)') ' '
212 write(weighttype(1:80),'(a)') "wtheta"
213 call ctrl_set_pack_xyz(
214 & cunit, ivartype, fname_theta(ictrlgrad), "maskCtrlC",
215 & weighttype, wtheta, lxxadxx, mythid)
216 #endif
217
218 #ifdef ALLOW_SALT0_CONTROL
219 ivartype = 2
220 write(weighttype(1:80),'(80a)') ' '
221 write(weighttype(1:80),'(a)') "wsalt"
222 call ctrl_set_pack_xyz(
223 & cunit, ivartype, fname_salt(ictrlgrad), "maskCtrlC",
224 & weighttype, wsalt, lxxadxx, mythid)
225 #endif
226
227 #if (defined (ALLOW_HFLUX_CONTROL) || \
228 defined (ALLOW_HFLUX0_CONTROL))
229 ivartype = 3
230 write(weighttype(1:80),'(80a)') ' '
231 write(weighttype(1:80),'(a)') "whflux"
232 call ctrl_set_pack_xy(
233 & cunit, ivartype, fname_hflux(ictrlgrad), "maskCtrlC",
234 & weighttype, lxxadxx, mythid)
235 #endif
236
237 #if (defined (ALLOW_SFLUX_CONTROL) || \
238 defined (ALLOW_SFLUX0_CONTROL))
239 ivartype = 4
240 write(weighttype(1:80),'(80a)') ' '
241 write(weighttype(1:80),'(a)') "wsflux"
242 call ctrl_set_pack_xy(
243 & cunit, ivartype, fname_sflux(ictrlgrad), "maskCtrlC",
244 & weighttype, lxxadxx, mythid)
245 #endif
246
247 #if (defined (ALLOW_USTRESS_CONTROL) || \
248 defined (ALLOW_TAUU0_CONTROL))
249 ivartype = 5
250 write(weighttype(1:80),'(80a)') ' '
251 write(weighttype(1:80),'(a)') "wtauu"
252 call ctrl_set_pack_xy(
253 & cunit, ivartype, fname_tauu(ictrlgrad), "maskCtrlW",
254 & weighttype, lxxadxx, mythid)
255 #endif
256
257 #if (defined (ALLOW_VSTRESS_CONTROL) || \
258 defined (ALLOW_TAUV0_CONTROL))
259 ivartype = 6
260 write(weighttype(1:80),'(80a)') ' '
261 write(weighttype(1:80),'(a)') "wtauv"
262 call ctrl_set_pack_xy(
263 & cunit, ivartype, fname_tauv(ictrlgrad), "maskCtrlS",
264 & weighttype, lxxadxx, mythid)
265 #endif
266
267 #ifdef ALLOW_ATEMP_CONTROL
268 ivartype = 7
269 write(weighttype(1:80),'(80a)') ' '
270 write(weighttype(1:80),'(a)') "watemp"
271 call ctrl_set_pack_xy(
272 & cunit, ivartype, fname_atemp(ictrlgrad), "maskCtrlC",
273 & weighttype, lxxadxx, mythid)
274 #endif
275
276 #ifdef ALLOW_AQH_CONTROL
277 ivartype = 8
278 write(weighttype(1:80),'(80a)') ' '
279 write(weighttype(1:80),'(a)') "waqh"
280 call ctrl_set_pack_xy(
281 & cunit, ivartype, fname_aqh(ictrlgrad), "maskCtrlC",
282 & weighttype, lxxadxx, mythid)
283 #endif
284
285 #ifdef ALLOW_UWIND_CONTROL
286 ivartype = 9
287 write(weighttype(1:80),'(80a)') ' '
288 write(weighttype(1:80),'(a)') "wuwind"
289 call ctrl_set_pack_xy(
290 & cunit, ivartype, fname_uwind(ictrlgrad), "maskCtrlC",
291 & weighttype, lxxadxx, mythid)
292 #endif
293
294 #ifdef ALLOW_VWIND_CONTROL
295 ivartype = 10
296 write(weighttype(1:80),'(80a)') ' '
297 write(weighttype(1:80),'(a)') "wvwind"
298 call ctrl_set_pack_xy(
299 & cunit, ivartype, fname_vwind(ictrlgrad), "maskCtrlC",
300 & weighttype, lxxadxx, mythid)
301 #endif
302
303 #ifdef ALLOW_OBCSN_CONTROL
304 ivartype = 11
305 write(weighttype(1:80),'(80a)') ' '
306 write(weighttype(1:80),'(a)') "wobcsn"
307 call ctrl_set_pack_xz(
308 & cunit, ivartype, fname_obcsn(ictrlgrad), "maskobcsn",
309 & weighttype, wobcsn, lxxadxx, mythid)
310 #endif
311
312 #ifdef ALLOW_OBCSS_CONTROL
313 ivartype = 12
314 write(weighttype(1:80),'(80a)') ' '
315 write(weighttype(1:80),'(a)') "wobcss"
316 call ctrl_set_pack_xz(
317 & cunit, ivartype, fname_obcss(ictrlgrad), "maskobcss",
318 & weighttype, wobcss, lxxadxx, mythid)
319 #endif
320
321 #ifdef ALLOW_OBCSW_CONTROL
322 ivartype = 13
323 write(weighttype(1:80),'(80a)') ' '
324 write(weighttype(1:80),'(a)') "wobcsw"
325 call ctrl_set_pack_yz(
326 & cunit, ivartype, fname_obcsw(ictrlgrad), "maskobcsw",
327 & weighttype, wobcsw, lxxadxx, mythid)
328 #endif
329
330 #ifdef ALLOW_OBCSE_CONTROL
331 ivartype = 14
332 write(weighttype(1:80),'(80a)') ' '
333 write(weighttype(1:80),'(a)') "wobcse"
334 call ctrl_set_pack_yz(
335 & cunit, ivartype, fname_obcse(ictrlgrad), "maskobcse",
336 & weighttype, wobcse, lxxadxx, mythid)
337 #endif
338
339 #ifdef ALLOW_DIFFKR_CONTROL
340 ivartype = 15
341 write(weighttype(1:80),'(80a)') ' '
342 write(weighttype(1:80),'(a)') "wdiffkr"
343 call ctrl_set_pack_xyz(
344 & cunit, ivartype, fname_diffkr(ictrlgrad), "maskCtrlC",
345 & weighttype, wunit, lxxadxx, mythid)
346 #endif
347
348 #ifdef ALLOW_KAPGM_CONTROL
349 ivartype = 16
350 write(weighttype(1:80),'(80a)') ' '
351 write(weighttype(1:80),'(a)') "wkapgm"
352 call ctrl_set_pack_xyz(
353 & cunit, ivartype, fname_kapgm(ictrlgrad), "maskCtrlC",
354 & weighttype, wunit, lxxadxx, mythid)
355 #endif
356
357 #ifdef ALLOW_TR10_CONTROL
358 ivartype = 17
359 write(weighttype(1:80),'(80a)') ' '
360 write(weighttype(1:80),'(a)') "wtr1"
361 call ctrl_set_pack_xyz(
362 & cunit, ivartype, fname_tr1(ictrlgrad), "maskCtrlC",
363 & weighttype, wunit, lxxadxx, mythid)
364 #endif
365
366 #ifdef ALLOW_SST0_CONTROL
367 ivartype = 18
368 write(weighttype(1:80),'(80a)') ' '
369 write(weighttype(1:80),'(a)') "wsst0"
370 call ctrl_set_pack_xy(
371 & cunit, ivartype, fname_sst(ictrlgrad), "maskCtrlC",
372 & weighttype, lxxadxx, mythid)
373 #endif
374
375 #ifdef ALLOW_SSS0_CONTROL
376 ivartype = 19
377 write(weighttype(1:80),'(80a)') ' '
378 write(weighttype(1:80),'(a)') "wsss0"
379 call ctrl_set_pack_xy(
380 & cunit, ivartype, fname_sss(ictrlgrad), "maskCtrlC",
381 & weighttype, lxxadxx, mythid)
382 #endif
383
384 #ifdef ALLOW_HFACC_CONTROL
385 ivartype = 20
386 write(weighttype(1:80),'(80a)') ' '
387 write(weighttype(1:80),'(a)') "whfacc"
388 # ifdef ALLOW_HFACC3D_CONTROL
389 call ctrl_set_pack_xyz(
390 & cunit, ivartype, fname_hfacc(ictrlgrad), "maskCtrlC",
391 & weighttype, wunit, lxxadxx, mythid)
392 # else
393 call ctrl_set_pack_xy(
394 & cunit, ivartype, fname_hfacc(ictrlgrad), "maskCtrlC",
395 & weighttype, lxxadxx, mythid)
396 # endif
397 #endif
398
399 #ifdef ALLOW_EFLUXY0_CONTROL
400 ivartype = 21
401 write(weighttype(1:80),'(80a)') ' '
402 write(weighttype(1:80),'(a)') "wefluxy0"
403 call ctrl_set_pack_xyz(
404 & cunit, ivartype, fname_efluxy(ictrlgrad), "maskCtrlS",
405 & weighttype, wunit, lxxadxx, mythid)
406 #endif
407
408 #ifdef ALLOW_EFLUXP0_CONTROL
409 ivartype = 22
410 write(weighttype(1:80),'(80a)') ' '
411 write(weighttype(1:80),'(a)') "wefluxp0"
412 call ctrl_set_pack_xyz(
413 & cunit, ivartype, fname_efluxp(ictrlgrad), "maskhFacV",
414 & weighttype, wunit, lxxadxx, mythid)
415 #endif
416
417 #ifdef ALLOW_BOTTOMDRAG_CONTROL
418 ivartype = 23
419 write(weighttype(1:80),'(80a)') ' '
420 write(weighttype(1:80),'(a)') "wbottomdrag"
421 call ctrl_set_pack_xy(
422 & cunit, ivartype, fname_bottomdrag(ictrlgrad), "maskCtrlC",
423 & weighttype, lxxadxx, mythid)
424 #endif
425
426 #ifdef ALLOW_EDTAUX_CONTROL
427 ivartype = 25
428 write(weighttype(1:80),'(80a)') ' '
429 write(weighttype(1:80),'(a)') "wedtaux"
430 call ctrl_set_pack_xyz(
431 & cunit, ivartype, fname_edtaux(ictrlgrad), "maskCtrlW",
432 & weighttype, wunit, lxxadxx, mythid)
433 #endif
434
435 #ifdef ALLOW_EDTAUY_CONTROL
436 ivartype = 26
437 write(weighttype(1:80),'(80a)') ' '
438 write(weighttype(1:80),'(a)') "wedtauy"
439 call ctrl_set_pack_xyz(
440 & cunit, ivartype, fname_edtauy(ictrlgrad), "maskCtrlS",
441 & weighttype, wunit, lxxadxx, mythid)
442 #endif
443
444 #ifdef ALLOW_UVEL0_CONTROL
445 ivartype = 27
446 write(weighttype(1:80),'(80a)') ' '
447 write(weighttype(1:80),'(a)') "wuvel"
448 call ctrl_set_pack_xyz(
449 & cunit, ivartype, fname_uvel(ictrlgrad), "maskCtrlW",
450 & weighttype, wunit, lxxadxx, mythid)
451 #endif
452
453 #ifdef ALLOW_VVEL0_CONTROL
454 ivartype = 28
455 write(weighttype(1:80),'(80a)') ' '
456 write(weighttype(1:80),'(a)') "wvvel"
457 call ctrl_set_pack_xyz(
458 & cunit, ivartype, fname_vvel(ictrlgrad), "maskCtrlS",
459 & weighttype, wunit, lxxadxx, mythid)
460 #endif
461
462 #ifdef ALLOW_ETAN0_CONTROL
463 ivartype = 29
464 write(weighttype(1:80),'(80a)') ' '
465 write(weighttype(1:80),'(a)') "wetan"
466 call ctrl_set_pack_xy(
467 & cunit, ivartype, fname_etan(ictrlgrad), "maskCtrlC",
468 & weighttype, lxxadxx, mythid)
469 #endif
470
471 #ifdef ALLOW_RELAXSST_CONTROL
472 ivartype = 30
473 write(weighttype(1:80),'(80a)') ' '
474 write(weighttype(1:80),'(a)') "wrelaxsst"
475 call ctrl_set_pack_xy(
476 & cunit, ivartype, fname_relaxsst(ictrlgrad), "maskCtrlC",
477 & weighttype, lxxadxx, mythid)
478 #endif
479
480 #ifdef ALLOW_RELAXSSS_CONTROL
481 ivartype = 31
482 write(weighttype(1:80),'(80a)') ' '
483 write(weighttype(1:80),'(a)') "wrelaxsss"
484 call ctrl_set_pack_xy(
485 & cunit, ivartype, fname_relaxsss(ictrlgrad), "maskCtrlC",
486 & weighttype, lxxadxx, mythid)
487 #endif
488
489 #ifdef ALLOW_PRECIP_CONTROL
490 ivartype = 32
491 write(weighttype(1:80),'(80a)') ' '
492 write(weighttype(1:80),'(a)') "wprecip"
493 call ctrl_set_pack_xy(
494 & cunit, ivartype, fname_precip(ictrlgrad), "maskCtrlC",
495 & weighttype, lxxadxx, mythid)
496 #endif
497
498 #ifdef ALLOW_SWFLUX_CONTROL
499 ivartype = 33
500 write(weighttype(1:80),'(80a)') ' '
501 write(weighttype(1:80),'(a)') "wswflux"
502 call ctrl_set_pack_xy(
503 & cunit, ivartype, fname_swflux(ictrlgrad), "maskCtrlC",
504 & weighttype, lxxadxx, mythid)
505 #endif
506
507 #ifdef ALLOW_SWDOWN_CONTROL
508 ivartype = 34
509 write(weighttype(1:80),'(80a)') ' '
510 write(weighttype(1:80),'(a)') "wswdown"
511 call ctrl_set_pack_xy(
512 & cunit, ivartype, fname_swdown(ictrlgrad), "maskCtrlC",
513 & weighttype, lxxadxx, mythid)
514 #endif
515
516 close ( cunit )
517
518 _END_MASTER( mythid )
519
520 #endif /* EXCLUDE_CTRL_PACK */
521
522 return
523 end
524

  ViewVC Help
Powered by ViewVC 1.1.22