/[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.21 - (show annotations) (download)
Thu Jul 28 13:47:49 2005 UTC (18 years, 10 months ago) by heimbach
Branch: MAIN
Changes since 1.20: +16 -6 lines
Adding precip control

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

  ViewVC Help
Powered by ViewVC 1.1.22