/[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.33 - (show annotations) (download)
Fri May 30 02:48:28 2008 UTC (16 years ago) by gforget
Branch: MAIN
CVS Tags: checkpoint62c, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62h, checkpoint60, checkpoint61, checkpoint62, checkpoint62b, checkpoint61f, checkpoint61n, checkpoint61q, checkpoint61e, checkpoint61g, checkpoint61d, checkpoint61b, checkpoint61c, checkpoint61a, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61j, checkpoint61k, checkpoint61h, checkpoint61i, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61z, checkpoint61x, checkpoint61y
Changes since 1.32: +2 -4 lines
o bridging the gap between eddy stress and GM.
  -> eddyTau is replaced with eddyPsi (eddyTau = f x rho0 x eddyPsi)
      along with a change in CPP option (now ALLOW_EDDYPSI).
  -> when using GM w/ GM_AdvForm:
      The total eddy streamfunction (Psi = eddyPsi + K x Slope)
      is applied either in the tracer Eq. or in momentum Eq.
      depending on data.gmredi (intro. GM_InMomAsStress).
  -> ALLOW_EDDYPSI_CONTROL for estimation purpose.
  The key modifications are in model/src/taueddy_external_forcing.F
  pkg/gmredi/gmredi_calc_*F pkg/gmredi/gmredi_*transport.F

1 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_pack.F,v 1.32 2008/02/02 02:34:50 gforget 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_lwflux_file, fname_lwflux, mythid)
119 call ctrl_set_fname(xx_lwdown_file, fname_lwdown, mythid)
120 call ctrl_set_fname(xx_evap_file, fname_evap, mythid)
121 call ctrl_set_fname(xx_snowprecip_file, fname_snowprecip, mythid)
122 call ctrl_set_fname(xx_apressure_file, fname_apressure, mythid)
123 call ctrl_set_fname(xx_runoff_file, fname_runoff, mythid)
124
125 call ctrl_set_fname(xx_uwind_file, fname_uwind, mythid)
126 call ctrl_set_fname(xx_vwind_file, fname_vwind, mythid)
127 call ctrl_set_fname(xx_obcsn_file, fname_obcsn, mythid)
128 call ctrl_set_fname(xx_obcss_file, fname_obcss, mythid)
129 call ctrl_set_fname(xx_obcsw_file, fname_obcsw, mythid)
130 call ctrl_set_fname(xx_obcse_file, fname_obcse, mythid)
131 call ctrl_set_fname(xx_diffkr_file, fname_diffkr, mythid)
132 call ctrl_set_fname(xx_kapgm_file, fname_kapgm, mythid)
133 call ctrl_set_fname(xx_kapredi_file, fname_kapredi, mythid)
134 call ctrl_set_fname(xx_tr1_file, fname_tr1, mythid)
135 call ctrl_set_fname(xx_sst_file, fname_sst, mythid)
136 call ctrl_set_fname(xx_sss_file, fname_sss, mythid)
137 call ctrl_set_fname(xx_depth_file, fname_depth, mythid)
138 call ctrl_set_fname(xx_efluxy_file, fname_efluxy, mythid)
139 call ctrl_set_fname(xx_efluxp_file, fname_efluxp, mythid)
140 call ctrl_set_fname(xx_bottomdrag_file, fname_bottomdrag, mythid)
141 call ctrl_set_fname(xx_edtaux_file, fname_edtaux, mythid)
142 call ctrl_set_fname(xx_edtauy_file, fname_edtauy, mythid)
143 call ctrl_set_fname(xx_uvel_file, fname_uvel, mythid)
144 call ctrl_set_fname(xx_vvel_file, fname_vvel, mythid)
145 call ctrl_set_fname(xx_etan_file, fname_etan, mythid)
146 call ctrl_set_fname(xx_relaxsst_file, fname_relaxsst, mythid)
147 call ctrl_set_fname(xx_relaxsss_file, fname_relaxsss, mythid)
148 call ctrl_set_fname(xx_siarea_file, fname_siarea, mythid)
149 call ctrl_set_fname(xx_siheff_file, fname_siheff, mythid)
150 call ctrl_set_fname(xx_sihsnow_file, fname_sihsnow, mythid)
151 cHFLUXM_CONTROL
152 call ctrl_set_fname(xx_hfluxm_file, fname_hfluxm, mythid)
153 cHFLUXM_CONTROL
154
155 c-- Only the master thread will do I/O.
156 _BEGIN_MASTER( mythid )
157
158 if ( first ) then
159 c >>> Initialise control vector for optimcycle=0 <<<
160 lxxadxx = .TRUE.
161 ictrlgrad = 1
162 fcloc = fmin
163 write(cfile(1:128),'(4a,i4.4)')
164 & ctrlname(1:9),'_',yctrlid(1:10),
165 & yctrlpospack, optimcycle
166 print *, 'ph-pack: packing ', ctrlname(1:9)
167 else
168 c >>> Write gradient vector <<<
169 lxxadxx = .FALSE.
170 ictrlgrad = 2
171 fcloc = fc
172 write(cfile(1:128),'(4a,i4.4)')
173 & costname(1:9),'_',yctrlid(1:10),
174 & yctrlpospack, optimcycle
175 print *, 'ph-pack: packing ', costname(1:9)
176 endif
177
178 call mdsfindunit( cunit, mythid )
179 open( cunit, file = cfile,
180 & status = 'unknown',
181 & form = 'unformatted',
182 & access = 'sequential' )
183
184 c-- Header information.
185 write(cunit) nvartype
186 write(cunit) nvarlength
187 write(cunit) yctrlid
188 write(cunit) optimCycle
189 write(cunit) fc
190 C place holder of obsolete variable iG
191 write(cunit) 1
192 C place holder of obsolete variable jG
193 write(cunit) 1
194 write(cunit) nsx
195 write(cunit) nsy
196 write(cunit) (nWetcGlobal(k), k=1,nr)
197 write(cunit) (nWetsGlobal(k), k=1,nr)
198 write(cunit) (nWetwGlobal(k), k=1,nr)
199 #ifdef ALLOW_CTRL_WETV
200 write(cunit) (nWetvGlobal(k), k=1,nr)
201 #endif
202
203 #ifdef ALLOW_OBCSN_CONTROL
204 write(cunit) ((nWetobcsnGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
205 #endif
206 #ifdef ALLOW_OBCSS_CONTROL
207 write(cunit) ((nWetobcssGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
208 #endif
209 #ifdef ALLOW_OBCSW_CONTROL
210 write(cunit) ((nWetobcswGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
211 #endif
212 #ifdef ALLOW_OBCSE_CONTROL
213 write(cunit) ((nWetobcseGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
214 #endif
215 write(cunit) (ncvarindex(i), i=1,maxcvars)
216 write(cunit) (ncvarrecs(i), i=1,maxcvars)
217 write(cunit) (ncvarxmax(i), i=1,maxcvars)
218 write(cunit) (ncvarymax(i), i=1,maxcvars)
219 write(cunit) (ncvarnrmax(i), i=1,maxcvars)
220 write(cunit) (ncvargrd(i), i=1,maxcvars)
221 write(cunit)
222
223 #ifdef ALLOW_THETA0_CONTROL
224 ivartype = 1
225 write(weighttype(1:80),'(80a)') ' '
226 write(weighttype(1:80),'(a)') "wthetaLev"
227 call ctrl_set_pack_xyz(
228 & cunit, ivartype, fname_theta(ictrlgrad), "maskCtrlC",
229 & weighttype, wtheta, lxxadxx, mythid)
230 #endif
231
232 #ifdef ALLOW_SALT0_CONTROL
233 ivartype = 2
234 write(weighttype(1:80),'(80a)') ' '
235 write(weighttype(1:80),'(a)') "wsaltLev"
236 call ctrl_set_pack_xyz(
237 & cunit, ivartype, fname_salt(ictrlgrad), "maskCtrlC",
238 & weighttype, wsalt, lxxadxx, mythid)
239 #endif
240
241 #if (defined (ALLOW_HFLUX_CONTROL) || defined (ALLOW_HFLUX0_CONTROL))
242 ivartype = 3
243 write(weighttype(1:80),'(80a)') ' '
244 write(weighttype(1:80),'(a)') "whflux"
245 call ctrl_set_pack_xy(
246 & cunit, ivartype, fname_hflux(ictrlgrad), "maskCtrlC",
247 & weighttype, lxxadxx, mythid)
248 #endif
249
250 #if (defined (ALLOW_SFLUX_CONTROL) || defined (ALLOW_SFLUX0_CONTROL))
251 ivartype = 4
252 write(weighttype(1:80),'(80a)') ' '
253 write(weighttype(1:80),'(a)') "wsflux"
254 call ctrl_set_pack_xy(
255 & cunit, ivartype, fname_sflux(ictrlgrad), "maskCtrlC",
256 & weighttype, lxxadxx, mythid)
257 #endif
258
259 #if (defined (ALLOW_USTRESS_CONTROL) || defined (ALLOW_TAUU0_CONTROL))
260 ivartype = 5
261 write(weighttype(1:80),'(80a)') ' '
262 write(weighttype(1:80),'(a)') "wtauu"
263 call ctrl_set_pack_xy(
264 & cunit, ivartype, fname_tauu(ictrlgrad), "maskCtrlW",
265 & weighttype, lxxadxx, mythid)
266 #endif
267
268 #if (defined (ALLOW_VSTRESS_CONTROL) || defined (ALLOW_TAUV0_CONTROL))
269 ivartype = 6
270 write(weighttype(1:80),'(80a)') ' '
271 write(weighttype(1:80),'(a)') "wtauv"
272 call ctrl_set_pack_xy(
273 & cunit, ivartype, fname_tauv(ictrlgrad), "maskCtrlS",
274 & weighttype, lxxadxx, mythid)
275 #endif
276
277 #ifdef ALLOW_ATEMP_CONTROL
278 ivartype = 7
279 write(weighttype(1:80),'(80a)') ' '
280 write(weighttype(1:80),'(a)') "watemp"
281 call ctrl_set_pack_xy(
282 & cunit, ivartype, fname_atemp(ictrlgrad), "maskCtrlC",
283 & weighttype, lxxadxx, mythid)
284 #endif
285
286 #ifdef ALLOW_AQH_CONTROL
287 ivartype = 8
288 write(weighttype(1:80),'(80a)') ' '
289 write(weighttype(1:80),'(a)') "waqh"
290 call ctrl_set_pack_xy(
291 & cunit, ivartype, fname_aqh(ictrlgrad), "maskCtrlC",
292 & weighttype, lxxadxx, mythid)
293 #endif
294
295 #ifdef ALLOW_UWIND_CONTROL
296 ivartype = 9
297 write(weighttype(1:80),'(80a)') ' '
298 write(weighttype(1:80),'(a)') "wuwind"
299 call ctrl_set_pack_xy(
300 & cunit, ivartype, fname_uwind(ictrlgrad), "maskCtrlC",
301 & weighttype, lxxadxx, mythid)
302 #endif
303
304 #ifdef ALLOW_VWIND_CONTROL
305 ivartype = 10
306 write(weighttype(1:80),'(80a)') ' '
307 write(weighttype(1:80),'(a)') "wvwind"
308 call ctrl_set_pack_xy(
309 & cunit, ivartype, fname_vwind(ictrlgrad), "maskCtrlC",
310 & weighttype, lxxadxx, mythid)
311 #endif
312
313 #ifdef ALLOW_OBCSN_CONTROL
314 ivartype = 11
315 write(weighttype(1:80),'(80a)') ' '
316 write(weighttype(1:80),'(a)') "wobcsn"
317 call ctrl_set_pack_xz(
318 & cunit, ivartype, fname_obcsn(ictrlgrad), "maskobcsn",
319 & weighttype, wobcsn, lxxadxx, mythid)
320 #endif
321
322 #ifdef ALLOW_OBCSS_CONTROL
323 ivartype = 12
324 write(weighttype(1:80),'(80a)') ' '
325 write(weighttype(1:80),'(a)') "wobcss"
326 call ctrl_set_pack_xz(
327 & cunit, ivartype, fname_obcss(ictrlgrad), "maskobcss",
328 & weighttype, wobcss, lxxadxx, mythid)
329 #endif
330
331 #ifdef ALLOW_OBCSW_CONTROL
332 ivartype = 13
333 write(weighttype(1:80),'(80a)') ' '
334 write(weighttype(1:80),'(a)') "wobcsw"
335 call ctrl_set_pack_yz(
336 & cunit, ivartype, fname_obcsw(ictrlgrad), "maskobcsw",
337 & weighttype, wobcsw, lxxadxx, mythid)
338 #endif
339
340 #ifdef ALLOW_OBCSE_CONTROL
341 ivartype = 14
342 write(weighttype(1:80),'(80a)') ' '
343 write(weighttype(1:80),'(a)') "wobcse"
344 call ctrl_set_pack_yz(
345 & cunit, ivartype, fname_obcse(ictrlgrad), "maskobcse",
346 & weighttype, wobcse, lxxadxx, mythid)
347 #endif
348
349 #ifdef ALLOW_DIFFKR_CONTROL
350 ivartype = 15
351 write(weighttype(1:80),'(80a)') ' '
352 write(weighttype(1:80),'(a)') "wdiffkr"
353 call ctrl_set_pack_xyz(
354 & cunit, ivartype, fname_diffkr(ictrlgrad), "maskCtrlC",
355 & weighttype, wdiffkr, lxxadxx, mythid)
356 #endif
357
358 #ifdef ALLOW_KAPGM_CONTROL
359 ivartype = 16
360 write(weighttype(1:80),'(80a)') ' '
361 write(weighttype(1:80),'(a)') "wkapgm"
362 call ctrl_set_pack_xyz(
363 & cunit, ivartype, fname_kapgm(ictrlgrad), "maskCtrlC",
364 & weighttype, wkapgm, lxxadxx, mythid)
365 #endif
366
367 #ifdef ALLOW_TR10_CONTROL
368 ivartype = 17
369 write(weighttype(1:80),'(80a)') ' '
370 write(weighttype(1:80),'(a)') "wtr1"
371 call ctrl_set_pack_xyz(
372 & cunit, ivartype, fname_tr1(ictrlgrad), "maskCtrlC",
373 & weighttype, wunit, lxxadxx, mythid)
374 #endif
375
376 #if (defined (ALLOW_SST_CONTROL) || defined (ALLOW_SST0_CONTROL))
377 ivartype = 18
378 write(weighttype(1:80),'(80a)') ' '
379 write(weighttype(1:80),'(a)') "wsst"
380 call ctrl_set_pack_xy(
381 & cunit, ivartype, fname_sst(ictrlgrad), "maskCtrlC",
382 & weighttype, lxxadxx, mythid)
383 #endif
384
385 #if (defined (ALLOW_SSS_CONTROL) || defined (ALLOW_SSS0_CONTROL))
386 ivartype = 19
387 write(weighttype(1:80),'(80a)') ' '
388 write(weighttype(1:80),'(a)') "wsss"
389 call ctrl_set_pack_xy(
390 & cunit, ivartype, fname_sss(ictrlgrad),
391 & "maskCtrlC", weighttype, lxxadxx, mythid)
392 #endif
393
394 #ifdef ALLOW_DEPTH_CONTROL
395 ivartype = 20
396 write(weighttype(1:80),'(80a)') ' '
397 write(weighttype(1:80),'(a)') "wdepth"
398 call ctrl_set_pack_xy(
399 & cunit, ivartype, fname_depth(ictrlgrad),
400 & "maskCtrlC", weighttype, lxxadxx, mythid)
401 #endif /* ALLOW_DEPTH_CONTROL */
402
403 #ifdef ALLOW_EFLUXY0_CONTROL
404 ivartype = 21
405 write(weighttype(1:80),'(80a)') ' '
406 write(weighttype(1:80),'(a)') "wefluxy0"
407 call ctrl_set_pack_xyz(
408 & cunit, ivartype, fname_efluxy(ictrlgrad), "maskCtrlS",
409 & weighttype, wunit, lxxadxx, mythid)
410 #endif
411
412 #ifdef ALLOW_EFLUXP0_CONTROL
413 ivartype = 22
414 write(weighttype(1:80),'(80a)') ' '
415 write(weighttype(1:80),'(a)') "wefluxp0"
416 call ctrl_set_pack_xyz(
417 & cunit, ivartype, fname_efluxp(ictrlgrad), "maskhFacV",
418 & weighttype, wunit, lxxadxx, mythid)
419 #endif
420
421 #ifdef ALLOW_BOTTOMDRAG_CONTROL
422 ivartype = 23
423 write(weighttype(1:80),'(80a)') ' '
424 write(weighttype(1:80),'(a)') "wbottomdrag"
425 call ctrl_set_pack_xy(
426 & cunit, ivartype, fname_bottomdrag(ictrlgrad), "maskCtrlC",
427 & weighttype, lxxadxx, mythid)
428 #endif
429
430 #ifdef ALLOW_HFLUXM_CONTROL
431 ivartype = 24
432 write(weighttype(1:80),'(80a)') ' '
433 write(weighttype(1:80),'(a)') "whfluxm"
434 call ctrl_set_pack_xy(
435 & cunit, ivartype, fname_hfluxm(ictrlgrad), "maskCtrlC",
436 & weighttype, lxxadxx, mythid)
437 #endif
438
439 #ifdef ALLOW_EDDYPSI_CONTROL
440 ivartype = 25
441 write(weighttype(1:80),'(80a)') ' '
442 write(weighttype(1:80),'(a)') "wedtaux"
443 call ctrl_set_pack_xyz(
444 & cunit, ivartype, fname_edtaux(ictrlgrad), "maskCtrlW",
445 & weighttype, wedtaux, lxxadxx, mythid)
446
447 ivartype = 26
448 write(weighttype(1:80),'(80a)') ' '
449 write(weighttype(1:80),'(a)') "wedtauy"
450 call ctrl_set_pack_xyz(
451 & cunit, ivartype, fname_edtauy(ictrlgrad), "maskCtrlS",
452 & weighttype, wedtauy, lxxadxx, mythid)
453 #endif
454
455 #ifdef ALLOW_UVEL0_CONTROL
456 ivartype = 27
457 write(weighttype(1:80),'(80a)') ' '
458 write(weighttype(1:80),'(a)') "wuvel"
459 call ctrl_set_pack_xyz(
460 & cunit, ivartype, fname_uvel(ictrlgrad), "maskCtrlW",
461 & weighttype, wunit, lxxadxx, mythid)
462 #endif
463
464 #ifdef ALLOW_VVEL0_CONTROL
465 ivartype = 28
466 write(weighttype(1:80),'(80a)') ' '
467 write(weighttype(1:80),'(a)') "wvvel"
468 call ctrl_set_pack_xyz(
469 & cunit, ivartype, fname_vvel(ictrlgrad), "maskCtrlS",
470 & weighttype, wunit, lxxadxx, mythid)
471 #endif
472
473 #ifdef ALLOW_ETAN0_CONTROL
474 ivartype = 29
475 write(weighttype(1:80),'(80a)') ' '
476 write(weighttype(1:80),'(a)') "wetan"
477 call ctrl_set_pack_xy(
478 & cunit, ivartype, fname_etan(ictrlgrad),
479 & "maskCtrlC", weighttype, lxxadxx, mythid)
480 #endif
481
482 #ifdef ALLOW_RELAXSST_CONTROL
483 ivartype = 30
484 write(weighttype(1:80),'(80a)') ' '
485 write(weighttype(1:80),'(a)') "wrelaxsst"
486 call ctrl_set_pack_xy(
487 & cunit, ivartype, fname_relaxsst(ictrlgrad),
488 & "maskCtrlC", weighttype, lxxadxx, mythid)
489 #endif
490
491 #ifdef ALLOW_RELAXSSS_CONTROL
492 ivartype = 31
493 write(weighttype(1:80),'(80a)') ' '
494 write(weighttype(1:80),'(a)') "wrelaxsss"
495 call ctrl_set_pack_xy(
496 & cunit, ivartype, fname_relaxsss(ictrlgrad),
497 & "maskCtrlC", weighttype, lxxadxx, mythid)
498 #endif
499
500 #ifdef ALLOW_PRECIP_CONTROL
501 ivartype = 32
502 write(weighttype(1:80),'(80a)') ' '
503 write(weighttype(1:80),'(a)') "wprecip"
504 call ctrl_set_pack_xy(
505 & cunit, ivartype, fname_precip(ictrlgrad),
506 & "maskCtrlC", weighttype, lxxadxx, mythid)
507 #endif
508
509 #ifdef ALLOW_SWFLUX_CONTROL
510 ivartype = 33
511 write(weighttype(1:80),'(80a)') ' '
512 write(weighttype(1:80),'(a)') "wswflux"
513 call ctrl_set_pack_xy(
514 & cunit, ivartype, fname_swflux(ictrlgrad),
515 & "maskCtrlC", weighttype, lxxadxx, mythid)
516 #endif
517
518 #ifdef ALLOW_SWDOWN_CONTROL
519 ivartype = 34
520 write(weighttype(1:80),'(80a)') ' '
521 write(weighttype(1:80),'(a)') "wswdown"
522 call ctrl_set_pack_xy(
523 & cunit, ivartype, fname_swdown(ictrlgrad),
524 & "maskCtrlC", weighttype, lxxadxx, mythid)
525 #endif
526
527 #ifdef ALLOW_LWFLUX_CONTROL
528 ivartype = 35
529 write(weighttype(1:80),'(80a)') ' '
530 write(weighttype(1:80),'(a)') "wlwflux"
531 call ctrl_set_pack_xy(
532 & cunit, ivartype, fname_lwflux(ictrlgrad),
533 & "maskCtrlC", weighttype, lxxadxx, mythid)
534 #endif
535
536 #ifdef ALLOW_LWDOWN_CONTROL
537 ivartype = 36
538 write(weighttype(1:80),'(80a)') ' '
539 write(weighttype(1:80),'(a)') "wlwdown"
540 call ctrl_set_pack_xy(
541 & cunit, ivartype, fname_lwdown(ictrlgrad),
542 & "maskCtrlC", weighttype, lxxadxx, mythid)
543 #endif
544
545 #ifdef ALLOW_EVAP_CONTROL
546 ivartype = 37
547 write(weighttype(1:80),'(80a)') ' '
548 write(weighttype(1:80),'(a)') "wevap"
549 call ctrl_set_pack_xy(
550 & cunit, ivartype, fname_evap(ictrlgrad),
551 & "maskCtrlC", weighttype, lxxadxx, mythid)
552 #endif
553
554 #ifdef ALLOW_SNOWPRECIP_CONTROL
555 ivartype = 38
556 write(weighttype(1:80),'(80a)') ' '
557 write(weighttype(1:80),'(a)') "wsnowprecip"
558 call ctrl_set_pack_xy(
559 & cunit, ivartype, fname_snowprecip(ictrlgrad),
560 & "maskCtrlC", weighttype, lxxadxx, mythid)
561 #endif
562
563 #ifdef ALLOW_APRESSURE_CONTROL
564 ivartype = 39
565 write(weighttype(1:80),'(80a)') ' '
566 write(weighttype(1:80),'(a)') "wapressure"
567 call ctrl_set_pack_xy(
568 & cunit, ivartype, fname_apressure(ictrlgrad),
569 & "maskCtrlC", weighttype, lxxadxx, mythid)
570 #endif
571
572 #ifdef ALLOW_RUNOFF_CONTROL
573 ivartype = 40
574 write(weighttype(1:80),'(80a)') ' '
575 write(weighttype(1:80),'(a)') "wrunoff"
576 call ctrl_set_pack_xy(
577 & cunit, ivartype, fname_runoff(ictrlgrad),
578 & "maskCtrlC", weighttype, lxxadxx, mythid)
579 #endif
580
581 #ifdef ALLOW_SIAREA_CONTROL
582 ivartype = 41
583 write(weighttype(1:80),'(80a)') ' '
584 write(weighttype(1:80),'(a)') "wunit"
585 call ctrl_set_pack_xy(
586 & cunit, ivartype, fname_siarea(ictrlgrad),
587 & "maskCtrlC", weighttype, lxxadxx, mythid)
588 #endif
589
590 #ifdef ALLOW_SIHEFF_CONTROL
591 ivartype = 42
592 write(weighttype(1:80),'(80a)') ' '
593 write(weighttype(1:80),'(a)') "wunit"
594 call ctrl_set_pack_xy(
595 & cunit, ivartype, fname_siheff(ictrlgrad),
596 & "maskCtrlC", weighttype, lxxadxx, mythid)
597 #endif
598
599 #ifdef ALLOW_SIHSNOW_CONTROL
600 ivartype = 43
601 write(weighttype(1:80),'(80a)') ' '
602 write(weighttype(1:80),'(a)') "wunit"
603 call ctrl_set_pack_xy(
604 & cunit, ivartype, fname_sihsnow(ictrlgrad),
605 & "maskCtrlC", weighttype, lxxadxx, mythid)
606 #endif
607
608 #ifdef ALLOW_KAPREDI_CONTROL
609 ivartype = 44
610 write(weighttype(1:80),'(80a)') ' '
611 write(weighttype(1:80),'(a)') "wkapredi"
612 call ctrl_set_pack_xyz(
613 & cunit, ivartype, fname_kapredi(ictrlgrad), "maskCtrlC",
614 & weighttype, wkapredi, lxxadxx, mythid)
615 #endif
616
617 close ( cunit )
618
619 _END_MASTER( mythid )
620
621 #endif /* EXCLUDE_CTRL_PACK */
622
623 return
624 end
625

  ViewVC Help
Powered by ViewVC 1.1.22