/[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.22 - (show annotations) (download)
Thu Jul 28 19:51:22 2005 UTC (18 years, 10 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57n_post, checkpoint57p_post, checkpoint57o_post
Changes since 1.21: +21 -11 lines
Adding swflux control

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

  ViewVC Help
Powered by ViewVC 1.1.22