/[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.31 - (hide annotations) (download)
Tue Jan 15 19:56:27 2008 UTC (16 years, 4 months ago) by dfer
Branch: MAIN
Changes since 1.30: +13 -1 lines
Bit of tutorial_global_oce_optim

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

  ViewVC Help
Powered by ViewVC 1.1.22