/[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.29 - (hide annotations) (download)
Fri Oct 27 05:16:54 2006 UTC (17 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint59, checkpoint58y_post, checkpoint58t_post, checkpoint58w_post, checkpoint58r_post, checkpoint59d, checkpoint59a, checkpoint59c, checkpoint59b, checkpoint58v_post, checkpoint58x_post, checkpoint58u_post, checkpoint58s_post
Changes since 1.28: +78 -17 lines
Adding new control variables:
lwflux, lwdown, evap, snowprecip, apressure, runoff.

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

  ViewVC Help
Powered by ViewVC 1.1.22