/[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.24 - (show annotations) (download)
Wed Aug 31 00:03:30 2005 UTC (18 years, 8 months ago) by heimbach
Branch: MAIN
Changes since 1.23: +7 -11 lines
Adding time-dependent SST, SSS control.

1 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_pack.F,v 1.23 2005/08/06 11:02:01 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) || 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) || defined (ALLOW_SFLUX0_CONTROL))
237 ivartype = 4
238 write(weighttype(1:80),'(80a)') ' '
239 write(weighttype(1:80),'(a)') "wsflux"
240 call ctrl_set_pack_xy(
241 & cunit, ivartype, fname_sflux(ictrlgrad), "maskCtrlC",
242 & weighttype, lxxadxx, mythid)
243 #endif
244
245 #if (defined (ALLOW_USTRESS_CONTROL) || 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) || defined (ALLOW_TAUV0_CONTROL))
255 ivartype = 6
256 write(weighttype(1:80),'(80a)') ' '
257 write(weighttype(1:80),'(a)') "wtauv"
258 call ctrl_set_pack_xy(
259 & cunit, ivartype, fname_tauv(ictrlgrad), "maskCtrlS",
260 & weighttype, lxxadxx, mythid)
261 #endif
262
263 #ifdef ALLOW_ATEMP_CONTROL
264 ivartype = 7
265 write(weighttype(1:80),'(80a)') ' '
266 write(weighttype(1:80),'(a)') "watemp"
267 call ctrl_set_pack_xy(
268 & cunit, ivartype, fname_atemp(ictrlgrad), "maskCtrlC",
269 & weighttype, lxxadxx, mythid)
270 #endif
271
272 #ifdef ALLOW_AQH_CONTROL
273 ivartype = 8
274 write(weighttype(1:80),'(80a)') ' '
275 write(weighttype(1:80),'(a)') "waqh"
276 call ctrl_set_pack_xy(
277 & cunit, ivartype, fname_aqh(ictrlgrad), "maskCtrlC",
278 & weighttype, lxxadxx, mythid)
279 #endif
280
281 #ifdef ALLOW_UWIND_CONTROL
282 ivartype = 9
283 write(weighttype(1:80),'(80a)') ' '
284 write(weighttype(1:80),'(a)') "wuwind"
285 call ctrl_set_pack_xy(
286 & cunit, ivartype, fname_uwind(ictrlgrad), "maskCtrlC",
287 & weighttype, lxxadxx, mythid)
288 #endif
289
290 #ifdef ALLOW_VWIND_CONTROL
291 ivartype = 10
292 write(weighttype(1:80),'(80a)') ' '
293 write(weighttype(1:80),'(a)') "wvwind"
294 call ctrl_set_pack_xy(
295 & cunit, ivartype, fname_vwind(ictrlgrad), "maskCtrlC",
296 & weighttype, lxxadxx, mythid)
297 #endif
298
299 #ifdef ALLOW_OBCSN_CONTROL
300 ivartype = 11
301 write(weighttype(1:80),'(80a)') ' '
302 write(weighttype(1:80),'(a)') "wobcsn"
303 call ctrl_set_pack_xz(
304 & cunit, ivartype, fname_obcsn(ictrlgrad), "maskobcsn",
305 & weighttype, wobcsn, lxxadxx, mythid)
306 #endif
307
308 #ifdef ALLOW_OBCSS_CONTROL
309 ivartype = 12
310 write(weighttype(1:80),'(80a)') ' '
311 write(weighttype(1:80),'(a)') "wobcss"
312 call ctrl_set_pack_xz(
313 & cunit, ivartype, fname_obcss(ictrlgrad), "maskobcss",
314 & weighttype, wobcss, lxxadxx, mythid)
315 #endif
316
317 #ifdef ALLOW_OBCSW_CONTROL
318 ivartype = 13
319 write(weighttype(1:80),'(80a)') ' '
320 write(weighttype(1:80),'(a)') "wobcsw"
321 call ctrl_set_pack_yz(
322 & cunit, ivartype, fname_obcsw(ictrlgrad), "maskobcsw",
323 & weighttype, wobcsw, lxxadxx, mythid)
324 #endif
325
326 #ifdef ALLOW_OBCSE_CONTROL
327 ivartype = 14
328 write(weighttype(1:80),'(80a)') ' '
329 write(weighttype(1:80),'(a)') "wobcse"
330 call ctrl_set_pack_yz(
331 & cunit, ivartype, fname_obcse(ictrlgrad), "maskobcse",
332 & weighttype, wobcse, lxxadxx, mythid)
333 #endif
334
335 #ifdef ALLOW_DIFFKR_CONTROL
336 ivartype = 15
337 write(weighttype(1:80),'(80a)') ' '
338 write(weighttype(1:80),'(a)') "wdiffkr"
339 call ctrl_set_pack_xyz(
340 & cunit, ivartype, fname_diffkr(ictrlgrad), "maskCtrlC",
341 & weighttype, wunit, lxxadxx, mythid)
342 #endif
343
344 #ifdef ALLOW_KAPGM_CONTROL
345 ivartype = 16
346 write(weighttype(1:80),'(80a)') ' '
347 write(weighttype(1:80),'(a)') "wkapgm"
348 call ctrl_set_pack_xyz(
349 & cunit, ivartype, fname_kapgm(ictrlgrad), "maskCtrlC",
350 & weighttype, wunit, lxxadxx, mythid)
351 #endif
352
353 #ifdef ALLOW_TR10_CONTROL
354 ivartype = 17
355 write(weighttype(1:80),'(80a)') ' '
356 write(weighttype(1:80),'(a)') "wtr1"
357 call ctrl_set_pack_xyz(
358 & cunit, ivartype, fname_tr1(ictrlgrad), "maskCtrlC",
359 & weighttype, wunit, lxxadxx, mythid)
360 #endif
361
362 #if (defined (ALLOW_SST_CONTROL) || defined (ALLOW_SST0_CONTROL))
363 ivartype = 18
364 write(weighttype(1:80),'(80a)') ' '
365 write(weighttype(1:80),'(a)') "wsst0"
366 call ctrl_set_pack_xy(
367 & cunit, ivartype, fname_sst(ictrlgrad), "maskCtrlC",
368 & weighttype, lxxadxx, mythid)
369 #endif
370
371 #if (defined (ALLOW_SSS_CONTROL) || defined (ALLOW_SSS0_CONTROL))
372 ivartype = 19
373 write(weighttype(1:80),'(80a)') ' '
374 write(weighttype(1:80),'(a)') "wsss0"
375 call ctrl_set_pack_xy(
376 & cunit, ivartype, fname_sss(ictrlgrad), "maskCtrlC",
377 & weighttype, lxxadxx, mythid)
378 #endif
379
380 #ifdef ALLOW_HFACC_CONTROL
381 ivartype = 20
382 write(weighttype(1:80),'(80a)') ' '
383 write(weighttype(1:80),'(a)') "whfacc"
384 # ifdef ALLOW_HFACC3D_CONTROL
385 call ctrl_set_pack_xyz(
386 & cunit, ivartype, fname_hfacc(ictrlgrad), "maskCtrlC",
387 & weighttype, wunit, lxxadxx, mythid)
388 # else
389 call ctrl_set_pack_xy(
390 & cunit, ivartype, fname_hfacc(ictrlgrad), "maskCtrlC",
391 & weighttype, lxxadxx, mythid)
392 # endif
393 #endif
394
395 #ifdef ALLOW_EFLUXY0_CONTROL
396 ivartype = 21
397 write(weighttype(1:80),'(80a)') ' '
398 write(weighttype(1:80),'(a)') "wefluxy0"
399 call ctrl_set_pack_xyz(
400 & cunit, ivartype, fname_efluxy(ictrlgrad), "maskCtrlS",
401 & weighttype, wunit, lxxadxx, mythid)
402 #endif
403
404 #ifdef ALLOW_EFLUXP0_CONTROL
405 ivartype = 22
406 write(weighttype(1:80),'(80a)') ' '
407 write(weighttype(1:80),'(a)') "wefluxp0"
408 call ctrl_set_pack_xyz(
409 & cunit, ivartype, fname_efluxp(ictrlgrad), "maskhFacV",
410 & weighttype, wunit, lxxadxx, mythid)
411 #endif
412
413 #ifdef ALLOW_BOTTOMDRAG_CONTROL
414 ivartype = 23
415 write(weighttype(1:80),'(80a)') ' '
416 write(weighttype(1:80),'(a)') "wbottomdrag"
417 call ctrl_set_pack_xy(
418 & cunit, ivartype, fname_bottomdrag(ictrlgrad), "maskCtrlC",
419 & weighttype, lxxadxx, mythid)
420 #endif
421
422 #ifdef ALLOW_EDTAUX_CONTROL
423 ivartype = 25
424 write(weighttype(1:80),'(80a)') ' '
425 write(weighttype(1:80),'(a)') "wedtaux"
426 call ctrl_set_pack_xyz(
427 & cunit, ivartype, fname_edtaux(ictrlgrad), "maskCtrlW",
428 & weighttype, wunit, lxxadxx, mythid)
429 #endif
430
431 #ifdef ALLOW_EDTAUY_CONTROL
432 ivartype = 26
433 write(weighttype(1:80),'(80a)') ' '
434 write(weighttype(1:80),'(a)') "wedtauy"
435 call ctrl_set_pack_xyz(
436 & cunit, ivartype, fname_edtauy(ictrlgrad), "maskCtrlS",
437 & weighttype, wunit, lxxadxx, mythid)
438 #endif
439
440 #ifdef ALLOW_UVEL0_CONTROL
441 ivartype = 27
442 write(weighttype(1:80),'(80a)') ' '
443 write(weighttype(1:80),'(a)') "wuvel"
444 call ctrl_set_pack_xyz(
445 & cunit, ivartype, fname_uvel(ictrlgrad), "maskCtrlW",
446 & weighttype, wunit, lxxadxx, mythid)
447 #endif
448
449 #ifdef ALLOW_VVEL0_CONTROL
450 ivartype = 28
451 write(weighttype(1:80),'(80a)') ' '
452 write(weighttype(1:80),'(a)') "wvvel"
453 call ctrl_set_pack_xyz(
454 & cunit, ivartype, fname_vvel(ictrlgrad), "maskCtrlS",
455 & weighttype, wunit, lxxadxx, mythid)
456 #endif
457
458 #ifdef ALLOW_ETAN0_CONTROL
459 ivartype = 29
460 write(weighttype(1:80),'(80a)') ' '
461 write(weighttype(1:80),'(a)') "wetan"
462 call ctrl_set_pack_xy(
463 & cunit, ivartype, fname_etan(ictrlgrad), "maskCtrlC",
464 & weighttype, lxxadxx, mythid)
465 #endif
466
467 #ifdef ALLOW_RELAXSST_CONTROL
468 ivartype = 30
469 write(weighttype(1:80),'(80a)') ' '
470 write(weighttype(1:80),'(a)') "wrelaxsst"
471 call ctrl_set_pack_xy(
472 & cunit, ivartype, fname_relaxsst(ictrlgrad), "maskCtrlC",
473 & weighttype, lxxadxx, mythid)
474 #endif
475
476 #ifdef ALLOW_RELAXSSS_CONTROL
477 ivartype = 31
478 write(weighttype(1:80),'(80a)') ' '
479 write(weighttype(1:80),'(a)') "wrelaxsss"
480 call ctrl_set_pack_xy(
481 & cunit, ivartype, fname_relaxsss(ictrlgrad), "maskCtrlC",
482 & weighttype, lxxadxx, mythid)
483 #endif
484
485 #ifdef ALLOW_PRECIP_CONTROL
486 ivartype = 32
487 write(weighttype(1:80),'(80a)') ' '
488 write(weighttype(1:80),'(a)') "wprecip"
489 call ctrl_set_pack_xy(
490 & cunit, ivartype, fname_precip(ictrlgrad), "maskCtrlC",
491 & weighttype, lxxadxx, mythid)
492 #endif
493
494 #ifdef ALLOW_SWFLUX_CONTROL
495 ivartype = 33
496 write(weighttype(1:80),'(80a)') ' '
497 write(weighttype(1:80),'(a)') "wswflux"
498 call ctrl_set_pack_xy(
499 & cunit, ivartype, fname_swflux(ictrlgrad), "maskCtrlC",
500 & weighttype, lxxadxx, mythid)
501 #endif
502
503 #ifdef ALLOW_SWDOWN_CONTROL
504 ivartype = 34
505 write(weighttype(1:80),'(80a)') ' '
506 write(weighttype(1:80),'(a)') "wswdown"
507 call ctrl_set_pack_xy(
508 & cunit, ivartype, fname_swdown(ictrlgrad), "maskCtrlC",
509 & weighttype, lxxadxx, mythid)
510 #endif
511
512 close ( cunit )
513
514 _END_MASTER( mythid )
515
516 #endif /* EXCLUDE_CTRL_PACK */
517
518 return
519 end
520

  ViewVC Help
Powered by ViewVC 1.1.22