/[MITgcm]/MITgcm/pkg/ctrl/ctrl_pack.F
ViewVC logotype

Annotation of /MITgcm/pkg/ctrl/ctrl_pack.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.34 - (hide annotations) (download)
Tue Jul 13 00:02:06 2010 UTC (13 years, 10 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint62k, checkpoint62j, checkpoint62i
Changes since 1.33: +9 -1 lines
ALLOW_ROTATE_UV_CONTROLS: when defined, we
rotate wind/stress controls adjustments
from Eastward/Northward to model grid directions.

1 gforget 1.34 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_pack.F,v 1.33 2008/05/30 02:48:28 gforget Exp $
2 heimbach 1.16 C $Name: $
3 heimbach 1.1
4 heimbach 1.12 #include "PACKAGES_CONFIG.h"
5 heimbach 1.1 #include "CTRL_CPPOPTIONS.h"
6    
7 heimbach 1.11 subroutine ctrl_pack( first, mythid )
8 heimbach 1.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 heimbach 1.12 c heimbach@mit.edu totally restructured 28-Oct-2003
30 heimbach 1.11 c
31 heimbach 1.8 c ==================================================================
32     c SUBROUTINE ctrl_pack
33     c ==================================================================
34    
35 heimbach 1.1 implicit none
36    
37     c == global variables ==
38 heimbach 1.5
39 heimbach 1.1 #include "EEPARAMS.h"
40     #include "SIZE.h"
41     #include "PARAMS.h"
42     #include "GRID.h"
43 heimbach 1.5
44 heimbach 1.1 #include "ctrl.h"
45 heimbach 1.14 #include "optim.h"
46 heimbach 1.5
47 heimbach 1.16 #ifdef ALLOW_COST
48     # include "cost.h"
49     #endif
50 heimbach 1.12 #ifdef ALLOW_ECCO
51     # include "ecco_cost.h"
52     #else
53     # include "ctrl_weights.h"
54     #endif
55    
56 heimbach 1.1 c == routine arguments ==
57 heimbach 1.5
58 heimbach 1.11 logical first
59 heimbach 1.1 integer mythid
60    
61 heimbach 1.11 #ifndef EXCLUDE_CTRL_PACK
62 heimbach 1.1 c == local variables ==
63    
64 heimbach 1.11 _RL fcloc
65    
66 heimbach 1.5 integer i, j, k
67 heimbach 1.1 integer ii
68     integer il
69     integer irec
70 heimbach 1.5 integer ig,jg
71     integer ivartype
72     integer iobcs
73 heimbach 1.1
74     logical doglobalread
75     logical ladinit
76 heimbach 1.5 integer cbuffindex
77 heimbach 1.11 logical lxxadxx
78    
79 heimbach 1.5 integer cunit
80 heimbach 1.11 integer ictrlgrad
81 heimbach 1.1
82     character*(128) cfile
83 heimbach 1.5 character*( 80) weighttype
84    
85 heimbach 1.1 c == external ==
86 heimbach 1.5
87 heimbach 1.1 integer ilnblnk
88     external ilnblnk
89    
90     c == end of interface ==
91    
92 heimbach 1.5 #ifndef ALLOW_ECCO_OPTIMIZATION
93 heimbach 1.11 fmin = 0. _d 0
94 heimbach 1.5 #endif
95 heimbach 1.1
96     c-- Tiled files are used.
97     doglobalread = .false.
98    
99     c-- Initialise adjoint variables on active files.
100     ladinit = .false.
101    
102 heimbach 1.14 c-- Initialise global buffer index
103     nbuffglobal = 0
104    
105 heimbach 1.5 c-- Assign file names.
106    
107 heimbach 1.11 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 heimbach 1.21 call ctrl_set_fname(xx_precip_file, fname_precip, mythid)
116 heimbach 1.22 call ctrl_set_fname(xx_swflux_file, fname_swflux, mythid)
117 heimbach 1.23 call ctrl_set_fname(xx_swdown_file, fname_swdown, mythid)
118 heimbach 1.29 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 heimbach 1.11 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 gforget 1.32 call ctrl_set_fname(xx_kapredi_file, fname_kapredi, mythid)
134 heimbach 1.11 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 heimbach 1.28 call ctrl_set_fname(xx_depth_file, fname_depth, mythid)
138 heimbach 1.11 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 heimbach 1.18 call ctrl_set_fname(xx_edtaux_file, fname_edtaux, mythid)
142     call ctrl_set_fname(xx_edtauy_file, fname_edtauy, mythid)
143 heimbach 1.19 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 heimbach 1.30 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 dfer 1.31 cHFLUXM_CONTROL
152     call ctrl_set_fname(xx_hfluxm_file, fname_hfluxm, mythid)
153     cHFLUXM_CONTROL
154 heimbach 1.5
155 heimbach 1.19 c-- Only the master thread will do I/O.
156 heimbach 1.1 _BEGIN_MASTER( mythid )
157    
158 heimbach 1.13 if ( first ) then
159 heimbach 1.11 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 heimbach 1.13 & ctrlname(1:9),'_',yctrlid(1:10),
165     & yctrlpospack, optimcycle
166 heimbach 1.17 print *, 'ph-pack: packing ', ctrlname(1:9)
167 heimbach 1.11 else
168 heimbach 1.1 c >>> Write gradient vector <<<
169 heimbach 1.11 lxxadxx = .FALSE.
170     ictrlgrad = 2
171     fcloc = fc
172 heimbach 1.5 write(cfile(1:128),'(4a,i4.4)')
173 heimbach 1.13 & costname(1:9),'_',yctrlid(1:10),
174     & yctrlpospack, optimcycle
175 heimbach 1.17 print *, 'ph-pack: packing ', costname(1:9)
176 heimbach 1.11 endif
177 heimbach 1.1
178 heimbach 1.11 call mdsfindunit( cunit, mythid )
179     open( cunit, file = cfile,
180     & status = 'unknown',
181     & form = 'unformatted',
182     & access = 'sequential' )
183 heimbach 1.1
184     c-- Header information.
185 mlosch 1.15 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 heimbach 1.7 #ifdef ALLOW_CTRL_WETV
200 mlosch 1.15 write(cunit) (nWetvGlobal(k), k=1,nr)
201 heimbach 1.7 #endif
202 heimbach 1.11
203 heimbach 1.5 #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 mlosch 1.15 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 heimbach 1.1 write(cunit)
222    
223     #ifdef ALLOW_THETA0_CONTROL
224 heimbach 1.5 ivartype = 1
225 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
226 heimbach 1.26 write(weighttype(1:80),'(a)') "wthetaLev"
227 heimbach 1.5 call ctrl_set_pack_xyz(
228 heimbach 1.19 & cunit, ivartype, fname_theta(ictrlgrad), "maskCtrlC",
229 heimbach 1.8 & weighttype, wtheta, lxxadxx, mythid)
230 heimbach 1.1 #endif
231    
232     #ifdef ALLOW_SALT0_CONTROL
233 heimbach 1.5 ivartype = 2
234 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
235 heimbach 1.26 write(weighttype(1:80),'(a)') "wsaltLev"
236 heimbach 1.5 call ctrl_set_pack_xyz(
237 heimbach 1.19 & cunit, ivartype, fname_salt(ictrlgrad), "maskCtrlC",
238 heimbach 1.8 & weighttype, wsalt, lxxadxx, mythid)
239 heimbach 1.5 #endif
240    
241 heimbach 1.24 #if (defined (ALLOW_HFLUX_CONTROL) || defined (ALLOW_HFLUX0_CONTROL))
242 heimbach 1.5 ivartype = 3
243     write(weighttype(1:80),'(80a)') ' '
244     write(weighttype(1:80),'(a)') "whflux"
245     call ctrl_set_pack_xy(
246 heimbach 1.19 & cunit, ivartype, fname_hflux(ictrlgrad), "maskCtrlC",
247 heimbach 1.11 & weighttype, lxxadxx, mythid)
248 heimbach 1.5 #endif
249    
250 heimbach 1.24 #if (defined (ALLOW_SFLUX_CONTROL) || defined (ALLOW_SFLUX0_CONTROL))
251 heimbach 1.5 ivartype = 4
252     write(weighttype(1:80),'(80a)') ' '
253     write(weighttype(1:80),'(a)') "wsflux"
254     call ctrl_set_pack_xy(
255 heimbach 1.19 & cunit, ivartype, fname_sflux(ictrlgrad), "maskCtrlC",
256 heimbach 1.11 & weighttype, lxxadxx, mythid)
257 heimbach 1.5 #endif
258    
259 heimbach 1.24 #if (defined (ALLOW_USTRESS_CONTROL) || defined (ALLOW_TAUU0_CONTROL))
260 heimbach 1.5 ivartype = 5
261     write(weighttype(1:80),'(80a)') ' '
262     write(weighttype(1:80),'(a)') "wtauu"
263     call ctrl_set_pack_xy(
264 gforget 1.34 #ifndef ALLOW_ROTATE_UV_CONTROLS
265 heimbach 1.19 & cunit, ivartype, fname_tauu(ictrlgrad), "maskCtrlW",
266 gforget 1.34 #else
267     & cunit, ivartype, fname_tauu(ictrlgrad), "maskCtrlC",
268     #endif
269 heimbach 1.11 & weighttype, lxxadxx, mythid)
270 heimbach 1.5 #endif
271    
272 heimbach 1.24 #if (defined (ALLOW_VSTRESS_CONTROL) || defined (ALLOW_TAUV0_CONTROL))
273 heimbach 1.5 ivartype = 6
274     write(weighttype(1:80),'(80a)') ' '
275     write(weighttype(1:80),'(a)') "wtauv"
276     call ctrl_set_pack_xy(
277 gforget 1.34 #ifndef ALLOW_ROTATE_UV_CONTROLS
278 heimbach 1.19 & cunit, ivartype, fname_tauv(ictrlgrad), "maskCtrlS",
279 gforget 1.34 #else
280     & cunit, ivartype, fname_tauv(ictrlgrad), "maskCtrlC",
281     #endif
282 heimbach 1.11 & weighttype, lxxadxx, mythid)
283 heimbach 1.5 #endif
284    
285     #ifdef ALLOW_ATEMP_CONTROL
286     ivartype = 7
287     write(weighttype(1:80),'(80a)') ' '
288     write(weighttype(1:80),'(a)') "watemp"
289     call ctrl_set_pack_xy(
290 heimbach 1.19 & cunit, ivartype, fname_atemp(ictrlgrad), "maskCtrlC",
291 heimbach 1.11 & weighttype, lxxadxx, mythid)
292 heimbach 1.5 #endif
293    
294     #ifdef ALLOW_AQH_CONTROL
295     ivartype = 8
296     write(weighttype(1:80),'(80a)') ' '
297     write(weighttype(1:80),'(a)') "waqh"
298     call ctrl_set_pack_xy(
299 heimbach 1.19 & cunit, ivartype, fname_aqh(ictrlgrad), "maskCtrlC",
300 heimbach 1.11 & weighttype, lxxadxx, mythid)
301 heimbach 1.5 #endif
302    
303     #ifdef ALLOW_UWIND_CONTROL
304     ivartype = 9
305     write(weighttype(1:80),'(80a)') ' '
306     write(weighttype(1:80),'(a)') "wuwind"
307     call ctrl_set_pack_xy(
308 heimbach 1.20 & cunit, ivartype, fname_uwind(ictrlgrad), "maskCtrlC",
309 heimbach 1.11 & weighttype, lxxadxx, mythid)
310 heimbach 1.5 #endif
311    
312     #ifdef ALLOW_VWIND_CONTROL
313     ivartype = 10
314     write(weighttype(1:80),'(80a)') ' '
315     write(weighttype(1:80),'(a)') "wvwind"
316     call ctrl_set_pack_xy(
317 heimbach 1.20 & cunit, ivartype, fname_vwind(ictrlgrad), "maskCtrlC",
318 heimbach 1.11 & weighttype, lxxadxx, mythid)
319 heimbach 1.5 #endif
320    
321     #ifdef ALLOW_OBCSN_CONTROL
322     ivartype = 11
323 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
324     write(weighttype(1:80),'(a)') "wobcsn"
325 heimbach 1.5 call ctrl_set_pack_xz(
326 heimbach 1.11 & cunit, ivartype, fname_obcsn(ictrlgrad), "maskobcsn",
327 heimbach 1.8 & weighttype, wobcsn, lxxadxx, mythid)
328 heimbach 1.5 #endif
329    
330     #ifdef ALLOW_OBCSS_CONTROL
331     ivartype = 12
332 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
333     write(weighttype(1:80),'(a)') "wobcss"
334 heimbach 1.5 call ctrl_set_pack_xz(
335 heimbach 1.11 & cunit, ivartype, fname_obcss(ictrlgrad), "maskobcss",
336 heimbach 1.8 & weighttype, wobcss, lxxadxx, mythid)
337 heimbach 1.5 #endif
338    
339     #ifdef ALLOW_OBCSW_CONTROL
340     ivartype = 13
341 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
342     write(weighttype(1:80),'(a)') "wobcsw"
343 heimbach 1.5 call ctrl_set_pack_yz(
344 heimbach 1.11 & cunit, ivartype, fname_obcsw(ictrlgrad), "maskobcsw",
345 heimbach 1.8 & weighttype, wobcsw, lxxadxx, mythid)
346 heimbach 1.5 #endif
347    
348     #ifdef ALLOW_OBCSE_CONTROL
349     ivartype = 14
350 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
351     write(weighttype(1:80),'(a)') "wobcse"
352 heimbach 1.5 call ctrl_set_pack_yz(
353 heimbach 1.11 & cunit, ivartype, fname_obcse(ictrlgrad), "maskobcse",
354 heimbach 1.8 & weighttype, wobcse, lxxadxx, mythid)
355 heimbach 1.1 #endif
356 heimbach 1.3
357     #ifdef ALLOW_DIFFKR_CONTROL
358 heimbach 1.5 ivartype = 15
359 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
360     write(weighttype(1:80),'(a)') "wdiffkr"
361 heimbach 1.5 call ctrl_set_pack_xyz(
362 heimbach 1.19 & cunit, ivartype, fname_diffkr(ictrlgrad), "maskCtrlC",
363 heimbach 1.27 & weighttype, wdiffkr, lxxadxx, mythid)
364 heimbach 1.3 #endif
365    
366     #ifdef ALLOW_KAPGM_CONTROL
367 heimbach 1.5 ivartype = 16
368 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
369     write(weighttype(1:80),'(a)') "wkapgm"
370 heimbach 1.5 call ctrl_set_pack_xyz(
371 heimbach 1.19 & cunit, ivartype, fname_kapgm(ictrlgrad), "maskCtrlC",
372 heimbach 1.27 & weighttype, wkapgm, lxxadxx, mythid)
373 heimbach 1.3 #endif
374    
375 heimbach 1.22 #ifdef ALLOW_TR10_CONTROL
376 heimbach 1.5 ivartype = 17
377 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
378 heimbach 1.22 write(weighttype(1:80),'(a)') "wtr1"
379     call ctrl_set_pack_xyz(
380     & cunit, ivartype, fname_tr1(ictrlgrad), "maskCtrlC",
381     & weighttype, wunit, lxxadxx, mythid)
382 heimbach 1.5 #endif
383    
384 heimbach 1.24 #if (defined (ALLOW_SST_CONTROL) || defined (ALLOW_SST0_CONTROL))
385 heimbach 1.6 ivartype = 18
386     write(weighttype(1:80),'(80a)') ' '
387 heimbach 1.25 write(weighttype(1:80),'(a)') "wsst"
388 heimbach 1.6 call ctrl_set_pack_xy(
389 heimbach 1.19 & cunit, ivartype, fname_sst(ictrlgrad), "maskCtrlC",
390 heimbach 1.11 & weighttype, lxxadxx, mythid)
391 heimbach 1.6 #endif
392    
393 heimbach 1.24 #if (defined (ALLOW_SSS_CONTROL) || defined (ALLOW_SSS0_CONTROL))
394 heimbach 1.6 ivartype = 19
395     write(weighttype(1:80),'(80a)') ' '
396 heimbach 1.25 write(weighttype(1:80),'(a)') "wsss"
397 heimbach 1.6 call ctrl_set_pack_xy(
398 heimbach 1.29 & cunit, ivartype, fname_sss(ictrlgrad),
399     & "maskCtrlC", weighttype, lxxadxx, mythid)
400 heimbach 1.6 #endif
401    
402 heimbach 1.28 #ifdef ALLOW_DEPTH_CONTROL
403 heimbach 1.6 ivartype = 20
404 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
405 heimbach 1.28 write(weighttype(1:80),'(a)') "wdepth"
406 heimbach 1.6 call ctrl_set_pack_xy(
407 heimbach 1.29 & cunit, ivartype, fname_depth(ictrlgrad),
408     & "maskCtrlC", weighttype, lxxadxx, mythid)
409 heimbach 1.28 #endif /* ALLOW_DEPTH_CONTROL */
410 heimbach 1.6
411 heimbach 1.5 #ifdef ALLOW_EFLUXY0_CONTROL
412     ivartype = 21
413 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
414     write(weighttype(1:80),'(a)') "wefluxy0"
415 heimbach 1.5 call ctrl_set_pack_xyz(
416 heimbach 1.19 & cunit, ivartype, fname_efluxy(ictrlgrad), "maskCtrlS",
417 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
418 heimbach 1.5 #endif
419    
420     #ifdef ALLOW_EFLUXP0_CONTROL
421     ivartype = 22
422 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
423     write(weighttype(1:80),'(a)') "wefluxp0"
424 heimbach 1.5 call ctrl_set_pack_xyz(
425 heimbach 1.19 & cunit, ivartype, fname_efluxp(ictrlgrad), "maskhFacV",
426 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
427 heimbach 1.6 #endif
428    
429     #ifdef ALLOW_BOTTOMDRAG_CONTROL
430     ivartype = 23
431     write(weighttype(1:80),'(80a)') ' '
432     write(weighttype(1:80),'(a)') "wbottomdrag"
433     call ctrl_set_pack_xy(
434 heimbach 1.19 & cunit, ivartype, fname_bottomdrag(ictrlgrad), "maskCtrlC",
435     & weighttype, lxxadxx, mythid)
436 heimbach 1.5 #endif
437    
438 dfer 1.31 #ifdef ALLOW_HFLUXM_CONTROL
439     ivartype = 24
440     write(weighttype(1:80),'(80a)') ' '
441     write(weighttype(1:80),'(a)') "whfluxm"
442     call ctrl_set_pack_xy(
443     & cunit, ivartype, fname_hfluxm(ictrlgrad), "maskCtrlC",
444     & weighttype, lxxadxx, mythid)
445     #endif
446    
447 gforget 1.33 #ifdef ALLOW_EDDYPSI_CONTROL
448 heimbach 1.18 ivartype = 25
449     write(weighttype(1:80),'(80a)') ' '
450     write(weighttype(1:80),'(a)') "wedtaux"
451     call ctrl_set_pack_xyz(
452 heimbach 1.19 & cunit, ivartype, fname_edtaux(ictrlgrad), "maskCtrlW",
453 heimbach 1.27 & weighttype, wedtaux, lxxadxx, mythid)
454 heimbach 1.18
455     ivartype = 26
456     write(weighttype(1:80),'(80a)') ' '
457     write(weighttype(1:80),'(a)') "wedtauy"
458     call ctrl_set_pack_xyz(
459 heimbach 1.19 & cunit, ivartype, fname_edtauy(ictrlgrad), "maskCtrlS",
460 heimbach 1.27 & weighttype, wedtauy, lxxadxx, mythid)
461 heimbach 1.19 #endif
462    
463     #ifdef ALLOW_UVEL0_CONTROL
464     ivartype = 27
465     write(weighttype(1:80),'(80a)') ' '
466     write(weighttype(1:80),'(a)') "wuvel"
467     call ctrl_set_pack_xyz(
468     & cunit, ivartype, fname_uvel(ictrlgrad), "maskCtrlW",
469 heimbach 1.18 & weighttype, wunit, lxxadxx, mythid)
470     #endif
471    
472 heimbach 1.19 #ifdef ALLOW_VVEL0_CONTROL
473     ivartype = 28
474     write(weighttype(1:80),'(80a)') ' '
475     write(weighttype(1:80),'(a)') "wvvel"
476     call ctrl_set_pack_xyz(
477     & cunit, ivartype, fname_vvel(ictrlgrad), "maskCtrlS",
478     & weighttype, wunit, lxxadxx, mythid)
479     #endif
480    
481     #ifdef ALLOW_ETAN0_CONTROL
482     ivartype = 29
483     write(weighttype(1:80),'(80a)') ' '
484     write(weighttype(1:80),'(a)') "wetan"
485     call ctrl_set_pack_xy(
486 heimbach 1.29 & cunit, ivartype, fname_etan(ictrlgrad),
487     & "maskCtrlC", weighttype, lxxadxx, mythid)
488 heimbach 1.19 #endif
489    
490     #ifdef ALLOW_RELAXSST_CONTROL
491     ivartype = 30
492     write(weighttype(1:80),'(80a)') ' '
493     write(weighttype(1:80),'(a)') "wrelaxsst"
494     call ctrl_set_pack_xy(
495 heimbach 1.29 & cunit, ivartype, fname_relaxsst(ictrlgrad),
496     & "maskCtrlC", weighttype, lxxadxx, mythid)
497 heimbach 1.19 #endif
498    
499     #ifdef ALLOW_RELAXSSS_CONTROL
500     ivartype = 31
501     write(weighttype(1:80),'(80a)') ' '
502     write(weighttype(1:80),'(a)') "wrelaxsss"
503     call ctrl_set_pack_xy(
504 heimbach 1.29 & cunit, ivartype, fname_relaxsss(ictrlgrad),
505     & "maskCtrlC", weighttype, lxxadxx, mythid)
506 heimbach 1.19 #endif
507 heimbach 1.18
508 heimbach 1.22 #ifdef ALLOW_PRECIP_CONTROL
509 heimbach 1.21 ivartype = 32
510     write(weighttype(1:80),'(80a)') ' '
511 heimbach 1.22 write(weighttype(1:80),'(a)') "wprecip"
512     call ctrl_set_pack_xy(
513 heimbach 1.29 & cunit, ivartype, fname_precip(ictrlgrad),
514     & "maskCtrlC", weighttype, lxxadxx, mythid)
515 heimbach 1.22 #endif
516    
517     #ifdef ALLOW_SWFLUX_CONTROL
518     ivartype = 33
519     write(weighttype(1:80),'(80a)') ' '
520     write(weighttype(1:80),'(a)') "wswflux"
521     call ctrl_set_pack_xy(
522 heimbach 1.29 & cunit, ivartype, fname_swflux(ictrlgrad),
523     & "maskCtrlC", weighttype, lxxadxx, mythid)
524 heimbach 1.21 #endif
525    
526 heimbach 1.23 #ifdef ALLOW_SWDOWN_CONTROL
527     ivartype = 34
528     write(weighttype(1:80),'(80a)') ' '
529     write(weighttype(1:80),'(a)') "wswdown"
530     call ctrl_set_pack_xy(
531 heimbach 1.29 & cunit, ivartype, fname_swdown(ictrlgrad),
532     & "maskCtrlC", weighttype, lxxadxx, mythid)
533     #endif
534    
535     #ifdef ALLOW_LWFLUX_CONTROL
536     ivartype = 35
537     write(weighttype(1:80),'(80a)') ' '
538     write(weighttype(1:80),'(a)') "wlwflux"
539     call ctrl_set_pack_xy(
540     & cunit, ivartype, fname_lwflux(ictrlgrad),
541     & "maskCtrlC", weighttype, lxxadxx, mythid)
542     #endif
543    
544     #ifdef ALLOW_LWDOWN_CONTROL
545     ivartype = 36
546     write(weighttype(1:80),'(80a)') ' '
547     write(weighttype(1:80),'(a)') "wlwdown"
548     call ctrl_set_pack_xy(
549     & cunit, ivartype, fname_lwdown(ictrlgrad),
550     & "maskCtrlC", weighttype, lxxadxx, mythid)
551     #endif
552    
553     #ifdef ALLOW_EVAP_CONTROL
554     ivartype = 37
555     write(weighttype(1:80),'(80a)') ' '
556     write(weighttype(1:80),'(a)') "wevap"
557     call ctrl_set_pack_xy(
558     & cunit, ivartype, fname_evap(ictrlgrad),
559     & "maskCtrlC", weighttype, lxxadxx, mythid)
560     #endif
561    
562     #ifdef ALLOW_SNOWPRECIP_CONTROL
563     ivartype = 38
564     write(weighttype(1:80),'(80a)') ' '
565     write(weighttype(1:80),'(a)') "wsnowprecip"
566     call ctrl_set_pack_xy(
567     & cunit, ivartype, fname_snowprecip(ictrlgrad),
568     & "maskCtrlC", weighttype, lxxadxx, mythid)
569     #endif
570    
571     #ifdef ALLOW_APRESSURE_CONTROL
572     ivartype = 39
573     write(weighttype(1:80),'(80a)') ' '
574     write(weighttype(1:80),'(a)') "wapressure"
575     call ctrl_set_pack_xy(
576     & cunit, ivartype, fname_apressure(ictrlgrad),
577     & "maskCtrlC", weighttype, lxxadxx, mythid)
578     #endif
579    
580     #ifdef ALLOW_RUNOFF_CONTROL
581     ivartype = 40
582     write(weighttype(1:80),'(80a)') ' '
583     write(weighttype(1:80),'(a)') "wrunoff"
584     call ctrl_set_pack_xy(
585     & cunit, ivartype, fname_runoff(ictrlgrad),
586     & "maskCtrlC", weighttype, lxxadxx, mythid)
587 heimbach 1.23 #endif
588    
589 heimbach 1.30 #ifdef ALLOW_SIAREA_CONTROL
590     ivartype = 41
591     write(weighttype(1:80),'(80a)') ' '
592     write(weighttype(1:80),'(a)') "wunit"
593     call ctrl_set_pack_xy(
594     & cunit, ivartype, fname_siarea(ictrlgrad),
595     & "maskCtrlC", weighttype, lxxadxx, mythid)
596     #endif
597    
598     #ifdef ALLOW_SIHEFF_CONTROL
599     ivartype = 42
600     write(weighttype(1:80),'(80a)') ' '
601     write(weighttype(1:80),'(a)') "wunit"
602     call ctrl_set_pack_xy(
603     & cunit, ivartype, fname_siheff(ictrlgrad),
604     & "maskCtrlC", weighttype, lxxadxx, mythid)
605     #endif
606    
607     #ifdef ALLOW_SIHSNOW_CONTROL
608     ivartype = 43
609     write(weighttype(1:80),'(80a)') ' '
610     write(weighttype(1:80),'(a)') "wunit"
611     call ctrl_set_pack_xy(
612     & cunit, ivartype, fname_sihsnow(ictrlgrad),
613     & "maskCtrlC", weighttype, lxxadxx, mythid)
614     #endif
615    
616 gforget 1.32 #ifdef ALLOW_KAPREDI_CONTROL
617     ivartype = 44
618     write(weighttype(1:80),'(80a)') ' '
619     write(weighttype(1:80),'(a)') "wkapredi"
620     call ctrl_set_pack_xyz(
621     & cunit, ivartype, fname_kapredi(ictrlgrad), "maskCtrlC",
622     & weighttype, wkapredi, lxxadxx, mythid)
623     #endif
624    
625 heimbach 1.5 close ( cunit )
626 heimbach 1.1
627 heimbach 1.5 _END_MASTER( mythid )
628 heimbach 1.11
629     #endif /* EXCLUDE_CTRL_PACK */
630 heimbach 1.1
631     return
632     end
633    

  ViewVC Help
Powered by ViewVC 1.1.22