/[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.30 - (hide annotations) (download)
Thu Jun 21 04:06:21 2007 UTC (16 years, 11 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint59e, checkpoint59g, checkpoint59f, checkpoint59m, checkpoint59l, checkpoint59i, checkpoint59h, checkpoint59k, checkpoint59j
Changes since 1.29: +31 -1 lines
Adding AREA, HEFF, HSNOW as control variables.

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

  ViewVC Help
Powered by ViewVC 1.1.22