/[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.44 - (show annotations) (download)
Tue Sep 11 01:32:02 2012 UTC (11 years, 9 months ago) by heimbach
Branch: MAIN
Changes since 1.43: +5 -1 lines
Merge OpenAD-specific code into main branch.

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

  ViewVC Help
Powered by ViewVC 1.1.22