/[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.45 - (show annotations) (download)
Tue Sep 11 23:46:15 2012 UTC (11 years, 9 months ago) by heimbach
Branch: MAIN
Changes since 1.44: +41 -3 lines
Enable packing/unpacking of generic control variables.
Compiles, but not yet tested.

1 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_pack.F,v 1.44 2012/09/11 01:32:02 heimbach Exp $
2 C $Name: $
3
4 #include "CTRL_OPTIONS.h"
5 #include "AD_CONFIG.h"
6 #ifdef ALLOW_EXF
7 # include "EXF_OPTIONS.h"
8 #endif
9
10 subroutine ctrl_pack( first, mythid )
11
12 c ==================================================================
13 c SUBROUTINE ctrl_pack
14 c ==================================================================
15 c
16 c o Compress the control vector such that only ocean points are
17 c written to file.
18 c
19 c started: Christian Eckert eckert@mit.edu 10-Mar=2000
20 c
21 c changed: Patrick Heimbach heimbach@mit.edu 06-Jun-2000
22 c - Transferred some filename declarations
23 c from here to namelist in ctrl_init
24 c
25 c Patrick Heimbach heimbach@mit.edu 16-Jun-2000
26 c - single file name convention with or without
27 c ALLOW_ECCO_OPTIMIZATION
28 c
29 c G. Gebbie, added open boundary control packing,
30 c gebbie@mit.edu 18 -Mar- 2003
31 c
32 c heimbach@mit.edu totally restructured 28-Oct-2003
33 c
34 c ==================================================================
35 c SUBROUTINE ctrl_pack
36 c ==================================================================
37
38 implicit none
39
40 c == global variables ==
41
42 #include "EEPARAMS.h"
43 #include "SIZE.h"
44 #include "PARAMS.h"
45 #include "GRID.h"
46
47 #include "ctrl.h"
48 #include "CTRL_SIZE.h"
49 #include "CTRL_GENARR.h"
50 #include "optim.h"
51
52 #ifdef ALLOW_COST
53 # include "cost.h"
54 #endif
55 #ifdef ALLOW_ECCO
56 # include "ecco_cost.h"
57 #else
58 # include "ctrl_weights.h"
59 #endif
60 #ifdef ALLOW_EXF
61 # include "EXF_PARAM.h"
62 #endif
63
64 c == routine arguments ==
65
66 logical first
67 integer mythid
68
69 #ifndef EXCLUDE_CTRL_PACK
70 #if (defined (ALLOW_ADJOINT_RUN)||defined (ALLOW_TANGENTLINEAR_RUN))
71 c == local variables ==
72
73 _RL fcloc
74
75 integer i, j, k
76 integer ii
77 integer il
78 integer irec
79 integer ig,jg
80 integer ivartype
81 integer iobcs
82 integer iarr
83
84 logical doglobalread
85 logical ladinit
86 integer cbuffindex
87 logical lxxadxx
88
89 integer cunit
90 integer ictrlgrad
91
92 character*(128) cfile
93 character*( 80) weighttype
94
95 c == external ==
96
97 integer ilnblnk
98 external ilnblnk
99
100 c == end of interface ==
101
102 #ifndef ALLOW_ECCO_OPTIMIZATION
103 fmin = 0. _d 0
104 #endif
105
106 c-- Tiled files are used.
107 doglobalread = .false.
108
109 c-- Initialise adjoint variables on active files.
110 ladinit = .false.
111
112 c-- Initialise global buffer index
113 nbuffglobal = 0
114
115 c-- Assign file names.
116
117 call ctrl_set_fname(xx_theta_file, fname_theta, mythid)
118 call ctrl_set_fname(xx_salt_file, fname_salt, mythid)
119 call ctrl_set_fname(xx_hflux_file, fname_hflux, mythid)
120 call ctrl_set_fname(xx_sflux_file, fname_sflux, mythid)
121 call ctrl_set_fname(xx_tauu_file, fname_tauu, mythid)
122 call ctrl_set_fname(xx_tauv_file, fname_tauv, mythid)
123 call ctrl_set_fname(xx_atemp_file, fname_atemp, mythid)
124 call ctrl_set_fname(xx_aqh_file, fname_aqh, mythid)
125 call ctrl_set_fname(xx_precip_file, fname_precip, mythid)
126 call ctrl_set_fname(xx_swflux_file, fname_swflux, mythid)
127 call ctrl_set_fname(xx_swdown_file, fname_swdown, mythid)
128 call ctrl_set_fname(xx_lwflux_file, fname_lwflux, mythid)
129 call ctrl_set_fname(xx_lwdown_file, fname_lwdown, mythid)
130 call ctrl_set_fname(xx_evap_file, fname_evap, mythid)
131 call ctrl_set_fname(xx_snowprecip_file, fname_snowprecip, mythid)
132 call ctrl_set_fname(xx_apressure_file, fname_apressure, mythid)
133 call ctrl_set_fname(xx_runoff_file, fname_runoff, mythid)
134 call ctrl_set_fname(xx_uwind_file, fname_uwind, mythid)
135 call ctrl_set_fname(xx_vwind_file, fname_vwind, mythid)
136 call ctrl_set_fname(xx_obcsn_file, fname_obcsn, mythid)
137 call ctrl_set_fname(xx_obcss_file, fname_obcss, mythid)
138 call ctrl_set_fname(xx_obcsw_file, fname_obcsw, mythid)
139 call ctrl_set_fname(xx_obcse_file, fname_obcse, mythid)
140 call ctrl_set_fname(xx_diffkr_file, fname_diffkr, mythid)
141 call ctrl_set_fname(xx_kapgm_file, fname_kapgm, mythid)
142 call ctrl_set_fname(xx_kapredi_file, fname_kapredi, mythid)
143 call ctrl_set_fname(xx_tr1_file, fname_tr1, mythid)
144 call ctrl_set_fname(xx_sst_file, fname_sst, mythid)
145 call ctrl_set_fname(xx_sss_file, fname_sss, mythid)
146 call ctrl_set_fname(xx_depth_file, fname_depth, mythid)
147 call ctrl_set_fname(xx_efluxy_file, fname_efluxy, mythid)
148 call ctrl_set_fname(xx_efluxp_file, fname_efluxp, mythid)
149 call ctrl_set_fname(xx_bottomdrag_file, fname_bottomdrag, mythid)
150 call ctrl_set_fname(xx_edtaux_file, fname_edtaux, mythid)
151 call ctrl_set_fname(xx_edtauy_file, fname_edtauy, mythid)
152 call ctrl_set_fname(xx_uvel_file, fname_uvel, mythid)
153 call ctrl_set_fname(xx_vvel_file, fname_vvel, mythid)
154 call ctrl_set_fname(xx_etan_file, fname_etan, mythid)
155 call ctrl_set_fname(xx_relaxsst_file, fname_relaxsst, mythid)
156 call ctrl_set_fname(xx_relaxsss_file, fname_relaxsss, mythid)
157 call ctrl_set_fname(xx_siarea_file, fname_siarea, mythid)
158 call ctrl_set_fname(xx_siheff_file, fname_siheff, mythid)
159 call ctrl_set_fname(xx_sihsnow_file, fname_sihsnow, mythid)
160 call ctrl_set_fname(xx_shifwflx_file, fname_shifwflx, mythid)
161 cHFLUXM_CONTROL
162 call ctrl_set_fname(xx_hfluxm_file, fname_hfluxm, mythid)
163 cHFLUXM_CONTROL
164 #ifdef ALLOW_GENARR2D_CONTROL
165 do iarr = 1, maxCtrlArr2D
166 call ctrl_set_fname(
167 & xx_genarr2d_file(iarr),
168 & fname_genarr2d(:,iarr), mythid)
169 enddo
170 #endif
171 #ifdef ALLOW_GENARR3D_CONTROL
172 do iarr = 1, maxCtrlArr3D
173 call ctrl_set_fname(
174 & xx_genarr3d_file(iarr),
175 & fname_genarr3d(:,iarr), mythid)
176 enddo
177 #endif
178
179 c-- Only the master thread will do I/O.
180 _BEGIN_MASTER( mythid )
181
182 if ( first ) then
183 c >>> Initialise control vector for optimcycle=0 <<<
184 lxxadxx = .TRUE.
185 ictrlgrad = 1
186 fcloc = fmin
187 write(cfile(1:128),'(4a,i4.4)')
188 & ctrlname(1:9),'_',yctrlid(1:10),
189 & yctrlpospack, optimcycle
190 print *, 'ph-pack: packing ', ctrlname(1:9)
191 else
192 c >>> Write gradient vector <<<
193 lxxadxx = .FALSE.
194 ictrlgrad = 2
195 #ifdef ALLOW_AUTODIFF_OPENAD
196 fcloc = fc%v
197 #else
198 fcloc = fc
199 #endif
200 write(cfile(1:128),'(4a,i4.4)')
201 & costname(1:9),'_',yctrlid(1:10),
202 & yctrlpospack, optimcycle
203 print *, 'ph-pack: packing ', costname(1:9)
204 endif
205
206 c-- Only Proc 0 will do I/O.
207 IF ( myProcId .eq. 0 ) THEN
208
209 call mdsfindunit( cunit, mythid )
210 open( cunit, file = cfile,
211 & status = 'unknown',
212 & form = 'unformatted',
213 & access = 'sequential' )
214
215 c-- Header information.
216 write(cunit) nvartype
217 write(cunit) nvarlength
218 write(cunit) yctrlid
219 write(cunit) optimCycle
220 write(cunit) fc
221 C place holder of obsolete variable iG
222 write(cunit) 1
223 C place holder of obsolete variable jG
224 write(cunit) 1
225 write(cunit) nsx
226 write(cunit) nsy
227 write(cunit) (nWetcGlobal(k), k=1,nr)
228 write(cunit) (nWetsGlobal(k), k=1,nr)
229 write(cunit) (nWetwGlobal(k), k=1,nr)
230 #ifdef ALLOW_CTRL_WETV
231 write(cunit) (nWetvGlobal(k), k=1,nr)
232 #endif
233 #ifdef ALLOW_SHIFWFLX_CONTROL
234 write(cunit) (nWetiGlobal(k), k=1,nr)
235 c write(cunit) nWetiGlobal(1)
236 #endif
237
238 #ifdef ALLOW_OBCSN_CONTROL
239 write(cunit) ((nWetobcsnGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
240 #endif
241 #ifdef ALLOW_OBCSS_CONTROL
242 write(cunit) ((nWetobcssGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
243 #endif
244 #ifdef ALLOW_OBCSW_CONTROL
245 write(cunit) ((nWetobcswGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
246 #endif
247 #ifdef ALLOW_OBCSE_CONTROL
248 write(cunit) ((nWetobcseGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
249 #endif
250 write(cunit) (ncvarindex(i), i=1,maxcvars)
251 write(cunit) (ncvarrecs(i), i=1,maxcvars)
252 write(cunit) (ncvarxmax(i), i=1,maxcvars)
253 write(cunit) (ncvarymax(i), i=1,maxcvars)
254 write(cunit) (ncvarnrmax(i), i=1,maxcvars)
255 write(cunit) (ncvargrd(i), i=1,maxcvars)
256 write(cunit)
257
258 #ifdef ALLOW_PACKUNPACK_METHOD2
259 ENDIF
260 _END_MASTER( mythid )
261 _BARRIER
262 #endif
263
264 #ifdef ALLOW_THETA0_CONTROL
265 ivartype = 1
266 write(weighttype(1:80),'(80a)') ' '
267 write(weighttype(1:80),'(a)') "wthetaLev"
268 call ctrl_set_pack_xyz(
269 & cunit, ivartype, fname_theta(ictrlgrad), "maskCtrlC",
270 & weighttype, wtheta, lxxadxx, mythid)
271 #endif
272
273 #ifdef ALLOW_SALT0_CONTROL
274 ivartype = 2
275 write(weighttype(1:80),'(80a)') ' '
276 write(weighttype(1:80),'(a)') "wsaltLev"
277 call ctrl_set_pack_xyz(
278 & cunit, ivartype, fname_salt(ictrlgrad), "maskCtrlC",
279 & weighttype, wsalt, lxxadxx, mythid)
280 #endif
281
282 #if (defined (ALLOW_HFLUX_CONTROL) || defined (ALLOW_HFLUX0_CONTROL))
283 ivartype = 3
284 write(weighttype(1:80),'(80a)') ' '
285 write(weighttype(1:80),'(a)') "whflux"
286 call ctrl_set_pack_xy(
287 & cunit, ivartype, fname_hflux(ictrlgrad), "maskCtrlC",
288 & weighttype, lxxadxx, mythid)
289 #endif
290
291 #if (defined (ALLOW_SFLUX_CONTROL) || defined (ALLOW_SFLUX0_CONTROL))
292 ivartype = 4
293 write(weighttype(1:80),'(80a)') ' '
294 write(weighttype(1:80),'(a)') "wsflux"
295 call ctrl_set_pack_xy(
296 & cunit, ivartype, fname_sflux(ictrlgrad), "maskCtrlC",
297 & weighttype, lxxadxx, mythid)
298 #endif
299
300 #if (defined (ALLOW_USTRESS_CONTROL) || defined (ALLOW_TAUU0_CONTROL))
301 #ifdef ALLOW_EXF
302 IF ( .NOT.useAtmWind ) THEN
303 #endif
304 ivartype = 5
305 write(weighttype(1:80),'(80a)') ' '
306 write(weighttype(1:80),'(a)') "wtauu"
307 call ctrl_set_pack_xy(
308 #ifndef ALLOW_ROTATE_UV_CONTROLS
309 & cunit, ivartype, fname_tauu(ictrlgrad), "maskCtrlW",
310 #else
311 & cunit, ivartype, fname_tauu(ictrlgrad), "maskCtrlC",
312 #endif
313 & weighttype, lxxadxx, mythid)
314 #ifdef ALLOW_EXF
315 ENDIF
316 #endif
317 #endif
318
319 #if (defined (ALLOW_VSTRESS_CONTROL) || defined (ALLOW_TAUV0_CONTROL))
320 #ifdef ALLOW_EXF
321 IF ( .NOT.useAtmWind ) THEN
322 #endif
323 ivartype = 6
324 write(weighttype(1:80),'(80a)') ' '
325 write(weighttype(1:80),'(a)') "wtauv"
326 call ctrl_set_pack_xy(
327 #ifndef ALLOW_ROTATE_UV_CONTROLS
328 & cunit, ivartype, fname_tauv(ictrlgrad), "maskCtrlS",
329 #else
330 & cunit, ivartype, fname_tauv(ictrlgrad), "maskCtrlC",
331 #endif
332 & weighttype, lxxadxx, mythid)
333 #ifdef ALLOW_EXF
334 ENDIF
335 #endif
336 #endif
337
338 #ifdef ALLOW_ATEMP_CONTROL
339 ivartype = 7
340 write(weighttype(1:80),'(80a)') ' '
341 write(weighttype(1:80),'(a)') "watemp"
342 call ctrl_set_pack_xy(
343 & cunit, ivartype, fname_atemp(ictrlgrad), "maskCtrlC",
344 & weighttype, lxxadxx, mythid)
345 #endif
346
347 #ifdef ALLOW_AQH_CONTROL
348 ivartype = 8
349 write(weighttype(1:80),'(80a)') ' '
350 write(weighttype(1:80),'(a)') "waqh"
351 call ctrl_set_pack_xy(
352 & cunit, ivartype, fname_aqh(ictrlgrad), "maskCtrlC",
353 & weighttype, lxxadxx, mythid)
354 #endif
355
356 #ifdef ALLOW_UWIND_CONTROL
357 #ifdef ALLOW_EXF
358 IF ( useAtmWind ) THEN
359 #endif
360 ivartype = 9
361 write(weighttype(1:80),'(80a)') ' '
362 write(weighttype(1:80),'(a)') "wuwind"
363 call ctrl_set_pack_xy(
364 & cunit, ivartype, fname_uwind(ictrlgrad), "maskCtrlC",
365 & weighttype, lxxadxx, mythid)
366 #ifdef ALLOW_EXF
367 ENDIF
368 #endif
369 #endif
370
371 #ifdef ALLOW_VWIND_CONTROL
372 #ifdef ALLOW_EXF
373 IF ( useAtmWind ) THEN
374 #endif
375 ivartype = 10
376 write(weighttype(1:80),'(80a)') ' '
377 write(weighttype(1:80),'(a)') "wvwind"
378 call ctrl_set_pack_xy(
379 & cunit, ivartype, fname_vwind(ictrlgrad), "maskCtrlC",
380 & weighttype, lxxadxx, mythid)
381 #ifdef ALLOW_EXF
382 ENDIF
383 #endif
384 #endif
385
386 #ifdef ALLOW_OBCSN_CONTROL
387 ivartype = 11
388 write(weighttype(1:80),'(80a)') ' '
389 write(weighttype(1:80),'(a)') "wobcsn"
390 call ctrl_set_pack_xz(
391 & cunit, ivartype, fname_obcsn(ictrlgrad), "maskobcsn",
392 & weighttype, wobcsn, lxxadxx, mythid)
393 #endif
394
395 #ifdef ALLOW_OBCSS_CONTROL
396 ivartype = 12
397 write(weighttype(1:80),'(80a)') ' '
398 write(weighttype(1:80),'(a)') "wobcss"
399 call ctrl_set_pack_xz(
400 & cunit, ivartype, fname_obcss(ictrlgrad), "maskobcss",
401 & weighttype, wobcss, lxxadxx, mythid)
402 #endif
403
404 #ifdef ALLOW_OBCSW_CONTROL
405 ivartype = 13
406 write(weighttype(1:80),'(80a)') ' '
407 write(weighttype(1:80),'(a)') "wobcsw"
408 call ctrl_set_pack_yz(
409 & cunit, ivartype, fname_obcsw(ictrlgrad), "maskobcsw",
410 & weighttype, wobcsw, lxxadxx, mythid)
411 #endif
412
413 #ifdef ALLOW_OBCSE_CONTROL
414 ivartype = 14
415 write(weighttype(1:80),'(80a)') ' '
416 write(weighttype(1:80),'(a)') "wobcse"
417 call ctrl_set_pack_yz(
418 & cunit, ivartype, fname_obcse(ictrlgrad), "maskobcse",
419 & weighttype, wobcse, lxxadxx, mythid)
420 #endif
421
422 #ifdef ALLOW_DIFFKR_CONTROL
423 ivartype = 15
424 write(weighttype(1:80),'(80a)') ' '
425 write(weighttype(1:80),'(a)') "wdiffkr"
426 call ctrl_set_pack_xyz(
427 & cunit, ivartype, fname_diffkr(ictrlgrad), "maskCtrlC",
428 & weighttype, wdiffkr, lxxadxx, mythid)
429 #endif
430
431 #ifdef ALLOW_KAPGM_CONTROL
432 ivartype = 16
433 write(weighttype(1:80),'(80a)') ' '
434 write(weighttype(1:80),'(a)') "wkapgm"
435 call ctrl_set_pack_xyz(
436 & cunit, ivartype, fname_kapgm(ictrlgrad), "maskCtrlC",
437 & weighttype, wkapgm, lxxadxx, mythid)
438 #endif
439
440 #ifdef ALLOW_TR10_CONTROL
441 ivartype = 17
442 write(weighttype(1:80),'(80a)') ' '
443 write(weighttype(1:80),'(a)') "wtr1"
444 call ctrl_set_pack_xyz(
445 & cunit, ivartype, fname_tr1(ictrlgrad), "maskCtrlC",
446 & weighttype, wunit, lxxadxx, mythid)
447 #endif
448
449 #if (defined (ALLOW_SST_CONTROL) || defined (ALLOW_SST0_CONTROL))
450 ivartype = 18
451 write(weighttype(1:80),'(80a)') ' '
452 write(weighttype(1:80),'(a)') "wsst"
453 call ctrl_set_pack_xy(
454 & cunit, ivartype, fname_sst(ictrlgrad), "maskCtrlC",
455 & weighttype, lxxadxx, mythid)
456 #endif
457
458 #if (defined (ALLOW_SSS_CONTROL) || defined (ALLOW_SSS0_CONTROL))
459 ivartype = 19
460 write(weighttype(1:80),'(80a)') ' '
461 write(weighttype(1:80),'(a)') "wsss"
462 call ctrl_set_pack_xy(
463 & cunit, ivartype, fname_sss(ictrlgrad),
464 & "maskCtrlC", weighttype, lxxadxx, mythid)
465 #endif
466
467 #ifdef ALLOW_DEPTH_CONTROL
468 ivartype = 20
469 write(weighttype(1:80),'(80a)') ' '
470 write(weighttype(1:80),'(a)') "wdepth"
471 call ctrl_set_pack_xy(
472 & cunit, ivartype, fname_depth(ictrlgrad),
473 & "maskCtrlC", weighttype, lxxadxx, mythid)
474 #endif /* ALLOW_DEPTH_CONTROL */
475
476 #ifdef ALLOW_EFLUXY0_CONTROL
477 ivartype = 21
478 write(weighttype(1:80),'(80a)') ' '
479 write(weighttype(1:80),'(a)') "wefluxy0"
480 call ctrl_set_pack_xyz(
481 & cunit, ivartype, fname_efluxy(ictrlgrad), "maskCtrlS",
482 & weighttype, wunit, lxxadxx, mythid)
483 #endif
484
485 #ifdef ALLOW_EFLUXP0_CONTROL
486 ivartype = 22
487 write(weighttype(1:80),'(80a)') ' '
488 write(weighttype(1:80),'(a)') "wefluxp0"
489 call ctrl_set_pack_xyz(
490 & cunit, ivartype, fname_efluxp(ictrlgrad), "maskhFacV",
491 & weighttype, wunit, lxxadxx, mythid)
492 #endif
493
494 #ifdef ALLOW_BOTTOMDRAG_CONTROL
495 ivartype = 23
496 write(weighttype(1:80),'(80a)') ' '
497 write(weighttype(1:80),'(a)') "wbottomdrag"
498 call ctrl_set_pack_xy(
499 & cunit, ivartype, fname_bottomdrag(ictrlgrad), "maskCtrlC",
500 & weighttype, lxxadxx, mythid)
501 #endif
502
503 #ifdef ALLOW_HFLUXM_CONTROL
504 ivartype = 24
505 write(weighttype(1:80),'(80a)') ' '
506 write(weighttype(1:80),'(a)') "whfluxm"
507 call ctrl_set_pack_xy(
508 & cunit, ivartype, fname_hfluxm(ictrlgrad), "maskCtrlC",
509 & weighttype, lxxadxx, mythid)
510 #endif
511
512 #ifdef ALLOW_EDDYPSI_CONTROL
513 ivartype = 25
514 write(weighttype(1:80),'(80a)') ' '
515 write(weighttype(1:80),'(a)') "wedtaux"
516 call ctrl_set_pack_xyz(
517 & cunit, ivartype, fname_edtaux(ictrlgrad), "maskCtrlW",
518 & weighttype, wedtaux, lxxadxx, mythid)
519
520 ivartype = 26
521 write(weighttype(1:80),'(80a)') ' '
522 write(weighttype(1:80),'(a)') "wedtauy"
523 call ctrl_set_pack_xyz(
524 & cunit, ivartype, fname_edtauy(ictrlgrad), "maskCtrlS",
525 & weighttype, wedtauy, lxxadxx, mythid)
526 #endif
527
528 #ifdef ALLOW_UVEL0_CONTROL
529 ivartype = 27
530 write(weighttype(1:80),'(80a)') ' '
531 write(weighttype(1:80),'(a)') "wuvel"
532 call ctrl_set_pack_xyz(
533 & cunit, ivartype, fname_uvel(ictrlgrad), "maskCtrlW",
534 & weighttype, wuvel, lxxadxx, mythid)
535 #endif
536
537 #ifdef ALLOW_VVEL0_CONTROL
538 ivartype = 28
539 write(weighttype(1:80),'(80a)') ' '
540 write(weighttype(1:80),'(a)') "wvvel"
541 call ctrl_set_pack_xyz(
542 & cunit, ivartype, fname_vvel(ictrlgrad), "maskCtrlS",
543 & weighttype, wvvel, lxxadxx, mythid)
544 #endif
545
546 #ifdef ALLOW_ETAN0_CONTROL
547 ivartype = 29
548 write(weighttype(1:80),'(80a)') ' '
549 write(weighttype(1:80),'(a)') "wetan"
550 call ctrl_set_pack_xy(
551 & cunit, ivartype, fname_etan(ictrlgrad),
552 & "maskCtrlC", weighttype, lxxadxx, mythid)
553 #endif
554
555 #ifdef ALLOW_RELAXSST_CONTROL
556 ivartype = 30
557 write(weighttype(1:80),'(80a)') ' '
558 write(weighttype(1:80),'(a)') "wrelaxsst"
559 call ctrl_set_pack_xy(
560 & cunit, ivartype, fname_relaxsst(ictrlgrad),
561 & "maskCtrlC", weighttype, lxxadxx, mythid)
562 #endif
563
564 #ifdef ALLOW_RELAXSSS_CONTROL
565 ivartype = 31
566 write(weighttype(1:80),'(80a)') ' '
567 write(weighttype(1:80),'(a)') "wrelaxsss"
568 call ctrl_set_pack_xy(
569 & cunit, ivartype, fname_relaxsss(ictrlgrad),
570 & "maskCtrlC", weighttype, lxxadxx, mythid)
571 #endif
572
573 #ifdef ALLOW_PRECIP_CONTROL
574 ivartype = 32
575 write(weighttype(1:80),'(80a)') ' '
576 write(weighttype(1:80),'(a)') "wprecip"
577 call ctrl_set_pack_xy(
578 & cunit, ivartype, fname_precip(ictrlgrad),
579 & "maskCtrlC", weighttype, lxxadxx, mythid)
580 #endif
581
582 #ifdef ALLOW_SWFLUX_CONTROL
583 ivartype = 33
584 write(weighttype(1:80),'(80a)') ' '
585 write(weighttype(1:80),'(a)') "wswflux"
586 call ctrl_set_pack_xy(
587 & cunit, ivartype, fname_swflux(ictrlgrad),
588 & "maskCtrlC", weighttype, lxxadxx, mythid)
589 #endif
590
591 #ifdef ALLOW_SWDOWN_CONTROL
592 ivartype = 34
593 write(weighttype(1:80),'(80a)') ' '
594 write(weighttype(1:80),'(a)') "wswdown"
595 call ctrl_set_pack_xy(
596 & cunit, ivartype, fname_swdown(ictrlgrad),
597 & "maskCtrlC", weighttype, lxxadxx, mythid)
598 #endif
599
600 #ifdef ALLOW_LWFLUX_CONTROL
601 ivartype = 35
602 write(weighttype(1:80),'(80a)') ' '
603 write(weighttype(1:80),'(a)') "wlwflux"
604 call ctrl_set_pack_xy(
605 & cunit, ivartype, fname_lwflux(ictrlgrad),
606 & "maskCtrlC", weighttype, lxxadxx, mythid)
607 #endif
608
609 #ifdef ALLOW_LWDOWN_CONTROL
610 ivartype = 36
611 write(weighttype(1:80),'(80a)') ' '
612 write(weighttype(1:80),'(a)') "wlwdown"
613 call ctrl_set_pack_xy(
614 & cunit, ivartype, fname_lwdown(ictrlgrad),
615 & "maskCtrlC", weighttype, lxxadxx, mythid)
616 #endif
617
618 #ifdef ALLOW_EVAP_CONTROL
619 ivartype = 37
620 write(weighttype(1:80),'(80a)') ' '
621 write(weighttype(1:80),'(a)') "wevap"
622 call ctrl_set_pack_xy(
623 & cunit, ivartype, fname_evap(ictrlgrad),
624 & "maskCtrlC", weighttype, lxxadxx, mythid)
625 #endif
626
627 #ifdef ALLOW_SNOWPRECIP_CONTROL
628 ivartype = 38
629 write(weighttype(1:80),'(80a)') ' '
630 write(weighttype(1:80),'(a)') "wsnowprecip"
631 call ctrl_set_pack_xy(
632 & cunit, ivartype, fname_snowprecip(ictrlgrad),
633 & "maskCtrlC", weighttype, lxxadxx, mythid)
634 #endif
635
636 #ifdef ALLOW_APRESSURE_CONTROL
637 ivartype = 39
638 write(weighttype(1:80),'(80a)') ' '
639 write(weighttype(1:80),'(a)') "wapressure"
640 call ctrl_set_pack_xy(
641 & cunit, ivartype, fname_apressure(ictrlgrad),
642 & "maskCtrlC", weighttype, lxxadxx, mythid)
643 #endif
644
645 #ifdef ALLOW_RUNOFF_CONTROL
646 ivartype = 40
647 write(weighttype(1:80),'(80a)') ' '
648 write(weighttype(1:80),'(a)') "wrunoff"
649 call ctrl_set_pack_xy(
650 & cunit, ivartype, fname_runoff(ictrlgrad),
651 & "maskCtrlC", weighttype, lxxadxx, mythid)
652 #endif
653
654 #ifdef ALLOW_SIAREA_CONTROL
655 ivartype = 41
656 write(weighttype(1:80),'(80a)') ' '
657 write(weighttype(1:80),'(a)') "wunit"
658 call ctrl_set_pack_xy(
659 & cunit, ivartype, fname_siarea(ictrlgrad),
660 & "maskCtrlC", weighttype, lxxadxx, mythid)
661 #endif
662
663 #ifdef ALLOW_SIHEFF_CONTROL
664 ivartype = 42
665 write(weighttype(1:80),'(80a)') ' '
666 write(weighttype(1:80),'(a)') "wunit"
667 call ctrl_set_pack_xy(
668 & cunit, ivartype, fname_siheff(ictrlgrad),
669 & "maskCtrlC", weighttype, lxxadxx, mythid)
670 #endif
671
672 #ifdef ALLOW_SIHSNOW_CONTROL
673 ivartype = 43
674 write(weighttype(1:80),'(80a)') ' '
675 write(weighttype(1:80),'(a)') "wunit"
676 call ctrl_set_pack_xy(
677 & cunit, ivartype, fname_sihsnow(ictrlgrad),
678 & "maskCtrlC", weighttype, lxxadxx, mythid)
679 #endif
680
681 #ifdef ALLOW_KAPREDI_CONTROL
682 ivartype = 44
683 write(weighttype(1:80),'(80a)') ' '
684 write(weighttype(1:80),'(a)') "wkapredi"
685 call ctrl_set_pack_xyz(
686 & cunit, ivartype, fname_kapredi(ictrlgrad), "maskCtrlC",
687 & weighttype, wkapredi, lxxadxx, mythid)
688 #endif
689
690 #ifdef ALLOW_SHIFWFLX_CONTROL
691 ivartype = 45
692 write(weighttype(1:80),'(80a)') ' '
693 write(weighttype(1:80),'(a)') "wshifwflx"
694 call ctrl_set_pack_xy(
695 & cunit, ivartype, fname_shifwflx(ictrlgrad),
696 & "maskCtrlI", weighttype, lxxadxx, mythid)
697 #endif
698
699 #ifdef ALLOW_GENARR2D_CONTROL
700 do iarr = 1, maxCtrlArr2D
701 ivartype = 100+iarr
702 write(weighttype(1:80),'(80a)') ' '
703 write(weighttype(1:80),'(a)') "wunit"
704 call ctrl_set_pack_xy(
705 & cunit, ivartype, fname_genarr2d(ictrlgrad,iarr),
706 & "maskCtrlC", weighttype, lxxadxx, mythid)
707 enddo
708 #endif
709
710 #ifdef ALLOW_GENARR3D_CONTROL
711 do iarr = 1, maxCtrlArr3D
712 ivartype = 200+iarr
713 write(weighttype(1:80),'(80a)') ' '
714 write(weighttype(1:80),'(a)') "wunit"
715 call ctrl_set_pack_xyz(
716 & cunit, ivartype, fname_genarr3d(ictrlgrad,iarr),
717 & "maskCtrlC", weighttype, wunit, lxxadxx, mythid)
718 enddo
719 #endif
720
721 #ifdef ALLOW_PACKUNPACK_METHOD2
722 _BEGIN_MASTER( mythid )
723 IF ( myProcId .eq. 0 ) THEN
724 #endif
725
726 close ( cunit )
727 ENDIF !IF ( myProcId .eq. 0 )
728 _END_MASTER( mythid )
729 _BARRIER
730 #endif /* (defined (ALLOW_ADJOINT_RUN)||defined (ALLOW_TANGENTLINEAR_RUN)) */
731 #endif /* EXCLUDE_CTRL_PACK */
732
733 return
734 end

  ViewVC Help
Powered by ViewVC 1.1.22