/[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.22 - (hide annotations) (download)
Thu Jul 28 19:51:22 2005 UTC (18 years, 10 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57n_post, checkpoint57p_post, checkpoint57o_post
Changes since 1.21: +21 -11 lines
Adding swflux control

1 heimbach 1.22 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_pack.F,v 1.21 2005/07/28 13:47:49 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.11 call ctrl_set_fname(xx_uwind_file, fname_uwind, mythid)
118     call ctrl_set_fname(xx_vwind_file, fname_vwind, mythid)
119     call ctrl_set_fname(xx_obcsn_file, fname_obcsn, mythid)
120     call ctrl_set_fname(xx_obcss_file, fname_obcss, mythid)
121     call ctrl_set_fname(xx_obcsw_file, fname_obcsw, mythid)
122     call ctrl_set_fname(xx_obcse_file, fname_obcse, mythid)
123     call ctrl_set_fname(xx_diffkr_file, fname_diffkr, mythid)
124     call ctrl_set_fname(xx_kapgm_file, fname_kapgm, mythid)
125     call ctrl_set_fname(xx_tr1_file, fname_tr1, mythid)
126     call ctrl_set_fname(xx_sst_file, fname_sst, mythid)
127     call ctrl_set_fname(xx_sss_file, fname_sss, mythid)
128     call ctrl_set_fname(xx_hfacc_file, fname_hfacc, mythid)
129     call ctrl_set_fname(xx_efluxy_file, fname_efluxy, mythid)
130     call ctrl_set_fname(xx_efluxp_file, fname_efluxp, mythid)
131     call ctrl_set_fname(xx_bottomdrag_file, fname_bottomdrag, mythid)
132 heimbach 1.18 call ctrl_set_fname(xx_edtaux_file, fname_edtaux, mythid)
133     call ctrl_set_fname(xx_edtauy_file, fname_edtauy, mythid)
134 heimbach 1.19 call ctrl_set_fname(xx_uvel_file, fname_uvel, mythid)
135     call ctrl_set_fname(xx_vvel_file, fname_vvel, mythid)
136     call ctrl_set_fname(xx_etan_file, fname_etan, mythid)
137     call ctrl_set_fname(xx_relaxsst_file, fname_relaxsst, mythid)
138     call ctrl_set_fname(xx_relaxsss_file, fname_relaxsss, mythid)
139 heimbach 1.5
140 heimbach 1.19 c-- Only the master thread will do I/O.
141 heimbach 1.1 _BEGIN_MASTER( mythid )
142    
143 heimbach 1.13 if ( first ) then
144 heimbach 1.11 c >>> Initialise control vector for optimcycle=0 <<<
145     lxxadxx = .TRUE.
146     ictrlgrad = 1
147     fcloc = fmin
148     write(cfile(1:128),'(4a,i4.4)')
149 heimbach 1.13 & ctrlname(1:9),'_',yctrlid(1:10),
150     & yctrlpospack, optimcycle
151 heimbach 1.17 print *, 'ph-pack: packing ', ctrlname(1:9)
152 heimbach 1.11 else
153 heimbach 1.1 c >>> Write gradient vector <<<
154 heimbach 1.11 lxxadxx = .FALSE.
155     ictrlgrad = 2
156     fcloc = fc
157 heimbach 1.5 write(cfile(1:128),'(4a,i4.4)')
158 heimbach 1.13 & costname(1:9),'_',yctrlid(1:10),
159     & yctrlpospack, optimcycle
160 heimbach 1.17 print *, 'ph-pack: packing ', costname(1:9)
161 heimbach 1.11 endif
162 heimbach 1.1
163 heimbach 1.11 call mdsfindunit( cunit, mythid )
164     open( cunit, file = cfile,
165     & status = 'unknown',
166     & form = 'unformatted',
167     & access = 'sequential' )
168 heimbach 1.1
169     c-- Header information.
170 mlosch 1.15 write(cunit) nvartype
171     write(cunit) nvarlength
172     write(cunit) yctrlid
173     write(cunit) optimCycle
174     write(cunit) fc
175     C place holder of obsolete variable iG
176     write(cunit) 1
177     C place holder of obsolete variable jG
178     write(cunit) 1
179     write(cunit) nsx
180     write(cunit) nsy
181     write(cunit) (nWetcGlobal(k), k=1,nr)
182     write(cunit) (nWetsGlobal(k), k=1,nr)
183     write(cunit) (nWetwGlobal(k), k=1,nr)
184 heimbach 1.7 #ifdef ALLOW_CTRL_WETV
185 mlosch 1.15 write(cunit) (nWetvGlobal(k), k=1,nr)
186 heimbach 1.7 #endif
187 heimbach 1.11
188 heimbach 1.5 #ifdef ALLOW_OBCSN_CONTROL
189     write(cunit) ((nWetobcsnGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
190     #endif
191     #ifdef ALLOW_OBCSS_CONTROL
192     write(cunit) ((nWetobcssGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
193     #endif
194     #ifdef ALLOW_OBCSW_CONTROL
195     write(cunit) ((nWetobcswGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
196     #endif
197     #ifdef ALLOW_OBCSE_CONTROL
198     write(cunit) ((nWetobcseGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
199     #endif
200 mlosch 1.15 write(cunit) (ncvarindex(i), i=1,maxcvars)
201     write(cunit) (ncvarrecs(i), i=1,maxcvars)
202     write(cunit) (ncvarxmax(i), i=1,maxcvars)
203     write(cunit) (ncvarymax(i), i=1,maxcvars)
204     write(cunit) (ncvarnrmax(i), i=1,maxcvars)
205     write(cunit) (ncvargrd(i), i=1,maxcvars)
206 heimbach 1.1 write(cunit)
207    
208     #ifdef ALLOW_THETA0_CONTROL
209 heimbach 1.5 ivartype = 1
210 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
211     write(weighttype(1:80),'(a)') "wtheta"
212 heimbach 1.5 call ctrl_set_pack_xyz(
213 heimbach 1.19 & cunit, ivartype, fname_theta(ictrlgrad), "maskCtrlC",
214 heimbach 1.8 & weighttype, wtheta, lxxadxx, mythid)
215 heimbach 1.1 #endif
216    
217     #ifdef ALLOW_SALT0_CONTROL
218 heimbach 1.5 ivartype = 2
219 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
220     write(weighttype(1:80),'(a)') "wsalt"
221 heimbach 1.5 call ctrl_set_pack_xyz(
222 heimbach 1.19 & cunit, ivartype, fname_salt(ictrlgrad), "maskCtrlC",
223 heimbach 1.8 & weighttype, wsalt, lxxadxx, mythid)
224 heimbach 1.5 #endif
225    
226     #if (defined (ALLOW_HFLUX_CONTROL) || \
227     defined (ALLOW_HFLUX0_CONTROL))
228     ivartype = 3
229     write(weighttype(1:80),'(80a)') ' '
230     write(weighttype(1:80),'(a)') "whflux"
231     call ctrl_set_pack_xy(
232 heimbach 1.19 & cunit, ivartype, fname_hflux(ictrlgrad), "maskCtrlC",
233 heimbach 1.11 & weighttype, lxxadxx, mythid)
234 heimbach 1.5 #endif
235    
236     #if (defined (ALLOW_SFLUX_CONTROL) || \
237     defined (ALLOW_SFLUX0_CONTROL))
238     ivartype = 4
239     write(weighttype(1:80),'(80a)') ' '
240     write(weighttype(1:80),'(a)') "wsflux"
241     call ctrl_set_pack_xy(
242 heimbach 1.19 & cunit, ivartype, fname_sflux(ictrlgrad), "maskCtrlC",
243 heimbach 1.11 & weighttype, lxxadxx, mythid)
244 heimbach 1.5 #endif
245    
246     #if (defined (ALLOW_USTRESS_CONTROL) || \
247     defined (ALLOW_TAUU0_CONTROL))
248     ivartype = 5
249     write(weighttype(1:80),'(80a)') ' '
250     write(weighttype(1:80),'(a)') "wtauu"
251     call ctrl_set_pack_xy(
252 heimbach 1.19 & cunit, ivartype, fname_tauu(ictrlgrad), "maskCtrlW",
253 heimbach 1.11 & weighttype, lxxadxx, mythid)
254 heimbach 1.5 #endif
255    
256     #if (defined (ALLOW_VSTRESS_CONTROL) || \
257     defined (ALLOW_TAUV0_CONTROL))
258     ivartype = 6
259     write(weighttype(1:80),'(80a)') ' '
260     write(weighttype(1:80),'(a)') "wtauv"
261     call ctrl_set_pack_xy(
262 heimbach 1.19 & cunit, ivartype, fname_tauv(ictrlgrad), "maskCtrlS",
263 heimbach 1.11 & weighttype, lxxadxx, mythid)
264 heimbach 1.5 #endif
265    
266     #ifdef ALLOW_ATEMP_CONTROL
267     ivartype = 7
268     write(weighttype(1:80),'(80a)') ' '
269     write(weighttype(1:80),'(a)') "watemp"
270     call ctrl_set_pack_xy(
271 heimbach 1.19 & cunit, ivartype, fname_atemp(ictrlgrad), "maskCtrlC",
272 heimbach 1.11 & weighttype, lxxadxx, mythid)
273 heimbach 1.5 #endif
274    
275     #ifdef ALLOW_AQH_CONTROL
276     ivartype = 8
277     write(weighttype(1:80),'(80a)') ' '
278     write(weighttype(1:80),'(a)') "waqh"
279     call ctrl_set_pack_xy(
280 heimbach 1.19 & cunit, ivartype, fname_aqh(ictrlgrad), "maskCtrlC",
281 heimbach 1.11 & weighttype, lxxadxx, mythid)
282 heimbach 1.5 #endif
283    
284     #ifdef ALLOW_UWIND_CONTROL
285     ivartype = 9
286     write(weighttype(1:80),'(80a)') ' '
287     write(weighttype(1:80),'(a)') "wuwind"
288     call ctrl_set_pack_xy(
289 heimbach 1.20 & cunit, ivartype, fname_uwind(ictrlgrad), "maskCtrlC",
290 heimbach 1.11 & weighttype, lxxadxx, mythid)
291 heimbach 1.5 #endif
292    
293     #ifdef ALLOW_VWIND_CONTROL
294     ivartype = 10
295     write(weighttype(1:80),'(80a)') ' '
296     write(weighttype(1:80),'(a)') "wvwind"
297     call ctrl_set_pack_xy(
298 heimbach 1.20 & cunit, ivartype, fname_vwind(ictrlgrad), "maskCtrlC",
299 heimbach 1.11 & weighttype, lxxadxx, mythid)
300 heimbach 1.5 #endif
301    
302     #ifdef ALLOW_OBCSN_CONTROL
303     ivartype = 11
304 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
305     write(weighttype(1:80),'(a)') "wobcsn"
306 heimbach 1.5 call ctrl_set_pack_xz(
307 heimbach 1.11 & cunit, ivartype, fname_obcsn(ictrlgrad), "maskobcsn",
308 heimbach 1.8 & weighttype, wobcsn, lxxadxx, mythid)
309 heimbach 1.5 #endif
310    
311     #ifdef ALLOW_OBCSS_CONTROL
312     ivartype = 12
313 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
314     write(weighttype(1:80),'(a)') "wobcss"
315 heimbach 1.5 call ctrl_set_pack_xz(
316 heimbach 1.11 & cunit, ivartype, fname_obcss(ictrlgrad), "maskobcss",
317 heimbach 1.8 & weighttype, wobcss, lxxadxx, mythid)
318 heimbach 1.5 #endif
319    
320     #ifdef ALLOW_OBCSW_CONTROL
321     ivartype = 13
322 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
323     write(weighttype(1:80),'(a)') "wobcsw"
324 heimbach 1.5 call ctrl_set_pack_yz(
325 heimbach 1.11 & cunit, ivartype, fname_obcsw(ictrlgrad), "maskobcsw",
326 heimbach 1.8 & weighttype, wobcsw, lxxadxx, mythid)
327 heimbach 1.5 #endif
328    
329     #ifdef ALLOW_OBCSE_CONTROL
330     ivartype = 14
331 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
332     write(weighttype(1:80),'(a)') "wobcse"
333 heimbach 1.5 call ctrl_set_pack_yz(
334 heimbach 1.11 & cunit, ivartype, fname_obcse(ictrlgrad), "maskobcse",
335 heimbach 1.8 & weighttype, wobcse, lxxadxx, mythid)
336 heimbach 1.1 #endif
337 heimbach 1.3
338     #ifdef ALLOW_DIFFKR_CONTROL
339 heimbach 1.5 ivartype = 15
340 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
341     write(weighttype(1:80),'(a)') "wdiffkr"
342 heimbach 1.5 call ctrl_set_pack_xyz(
343 heimbach 1.19 & cunit, ivartype, fname_diffkr(ictrlgrad), "maskCtrlC",
344 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
345 heimbach 1.3 #endif
346    
347     #ifdef ALLOW_KAPGM_CONTROL
348 heimbach 1.5 ivartype = 16
349 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
350     write(weighttype(1:80),'(a)') "wkapgm"
351 heimbach 1.5 call ctrl_set_pack_xyz(
352 heimbach 1.19 & cunit, ivartype, fname_kapgm(ictrlgrad), "maskCtrlC",
353 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
354 heimbach 1.3 #endif
355    
356 heimbach 1.22 #ifdef ALLOW_TR10_CONTROL
357 heimbach 1.5 ivartype = 17
358 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
359 heimbach 1.22 write(weighttype(1:80),'(a)') "wtr1"
360     call ctrl_set_pack_xyz(
361     & cunit, ivartype, fname_tr1(ictrlgrad), "maskCtrlC",
362     & weighttype, wunit, lxxadxx, mythid)
363 heimbach 1.5 #endif
364    
365 heimbach 1.6 #ifdef ALLOW_SST0_CONTROL
366     ivartype = 18
367     write(weighttype(1:80),'(80a)') ' '
368     write(weighttype(1:80),'(a)') "wsst0"
369     call ctrl_set_pack_xy(
370 heimbach 1.19 & cunit, ivartype, fname_sst(ictrlgrad), "maskCtrlC",
371 heimbach 1.11 & weighttype, lxxadxx, mythid)
372 heimbach 1.6 #endif
373    
374     #ifdef ALLOW_SSS0_CONTROL
375     ivartype = 19
376     write(weighttype(1:80),'(80a)') ' '
377     write(weighttype(1:80),'(a)') "wsss0"
378     call ctrl_set_pack_xy(
379 heimbach 1.19 & cunit, ivartype, fname_sss(ictrlgrad), "maskCtrlC",
380 heimbach 1.11 & weighttype, lxxadxx, mythid)
381 heimbach 1.6 #endif
382    
383     #ifdef ALLOW_HFACC_CONTROL
384     ivartype = 20
385 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
386     write(weighttype(1:80),'(a)') "whfacc"
387     # ifdef ALLOW_HFACC3D_CONTROL
388 heimbach 1.6 call ctrl_set_pack_xyz(
389 heimbach 1.19 & cunit, ivartype, fname_hfacc(ictrlgrad), "maskCtrlC",
390 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
391     # else
392 heimbach 1.6 call ctrl_set_pack_xy(
393 heimbach 1.19 & cunit, ivartype, fname_hfacc(ictrlgrad), "maskCtrlC",
394 heimbach 1.11 & weighttype, lxxadxx, mythid)
395 heimbach 1.8 # endif
396 heimbach 1.6 #endif
397    
398 heimbach 1.5 #ifdef ALLOW_EFLUXY0_CONTROL
399     ivartype = 21
400 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
401     write(weighttype(1:80),'(a)') "wefluxy0"
402 heimbach 1.5 call ctrl_set_pack_xyz(
403 heimbach 1.19 & cunit, ivartype, fname_efluxy(ictrlgrad), "maskCtrlS",
404 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
405 heimbach 1.5 #endif
406    
407     #ifdef ALLOW_EFLUXP0_CONTROL
408     ivartype = 22
409 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
410     write(weighttype(1:80),'(a)') "wefluxp0"
411 heimbach 1.5 call ctrl_set_pack_xyz(
412 heimbach 1.19 & cunit, ivartype, fname_efluxp(ictrlgrad), "maskhFacV",
413 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
414 heimbach 1.6 #endif
415    
416     #ifdef ALLOW_BOTTOMDRAG_CONTROL
417     ivartype = 23
418     write(weighttype(1:80),'(80a)') ' '
419     write(weighttype(1:80),'(a)') "wbottomdrag"
420     call ctrl_set_pack_xy(
421 heimbach 1.19 & cunit, ivartype, fname_bottomdrag(ictrlgrad), "maskCtrlC",
422     & weighttype, lxxadxx, mythid)
423 heimbach 1.5 #endif
424    
425 heimbach 1.18 #ifdef ALLOW_EDTAUX_CONTROL
426     ivartype = 25
427     write(weighttype(1:80),'(80a)') ' '
428     write(weighttype(1:80),'(a)') "wedtaux"
429     call ctrl_set_pack_xyz(
430 heimbach 1.19 & cunit, ivartype, fname_edtaux(ictrlgrad), "maskCtrlW",
431 heimbach 1.18 & weighttype, wunit, lxxadxx, mythid)
432     #endif
433    
434     #ifdef ALLOW_EDTAUY_CONTROL
435     ivartype = 26
436     write(weighttype(1:80),'(80a)') ' '
437     write(weighttype(1:80),'(a)') "wedtauy"
438     call ctrl_set_pack_xyz(
439 heimbach 1.19 & cunit, ivartype, fname_edtauy(ictrlgrad), "maskCtrlS",
440     & weighttype, wunit, lxxadxx, mythid)
441     #endif
442    
443     #ifdef ALLOW_UVEL0_CONTROL
444     ivartype = 27
445     write(weighttype(1:80),'(80a)') ' '
446     write(weighttype(1:80),'(a)') "wuvel"
447     call ctrl_set_pack_xyz(
448     & cunit, ivartype, fname_uvel(ictrlgrad), "maskCtrlW",
449 heimbach 1.18 & weighttype, wunit, lxxadxx, mythid)
450     #endif
451    
452 heimbach 1.19 #ifdef ALLOW_VVEL0_CONTROL
453     ivartype = 28
454     write(weighttype(1:80),'(80a)') ' '
455     write(weighttype(1:80),'(a)') "wvvel"
456     call ctrl_set_pack_xyz(
457     & cunit, ivartype, fname_vvel(ictrlgrad), "maskCtrlS",
458     & weighttype, wunit, lxxadxx, mythid)
459     #endif
460    
461     #ifdef ALLOW_ETAN0_CONTROL
462     ivartype = 29
463     write(weighttype(1:80),'(80a)') ' '
464     write(weighttype(1:80),'(a)') "wetan"
465     call ctrl_set_pack_xy(
466     & cunit, ivartype, fname_etan(ictrlgrad), "maskCtrlC",
467     & weighttype, lxxadxx, mythid)
468     #endif
469    
470     #ifdef ALLOW_RELAXSST_CONTROL
471     ivartype = 30
472     write(weighttype(1:80),'(80a)') ' '
473     write(weighttype(1:80),'(a)') "wrelaxsst"
474     call ctrl_set_pack_xy(
475     & cunit, ivartype, fname_relaxsst(ictrlgrad), "maskCtrlC",
476     & weighttype, lxxadxx, mythid)
477     #endif
478    
479     #ifdef ALLOW_RELAXSSS_CONTROL
480     ivartype = 31
481     write(weighttype(1:80),'(80a)') ' '
482     write(weighttype(1:80),'(a)') "wrelaxsss"
483     call ctrl_set_pack_xy(
484     & cunit, ivartype, fname_relaxsss(ictrlgrad), "maskCtrlC",
485     & weighttype, lxxadxx, mythid)
486     #endif
487 heimbach 1.18
488 heimbach 1.22 #ifdef ALLOW_PRECIP_CONTROL
489 heimbach 1.21 ivartype = 32
490     write(weighttype(1:80),'(80a)') ' '
491 heimbach 1.22 write(weighttype(1:80),'(a)') "wprecip"
492     call ctrl_set_pack_xy(
493     & cunit, ivartype, fname_precip(ictrlgrad), "maskCtrlC",
494     & weighttype, lxxadxx, mythid)
495     #endif
496    
497     #ifdef ALLOW_SWFLUX_CONTROL
498     ivartype = 33
499     write(weighttype(1:80),'(80a)') ' '
500     write(weighttype(1:80),'(a)') "wswflux"
501     call ctrl_set_pack_xy(
502     & cunit, ivartype, fname_swflux(ictrlgrad), "maskCtrlC",
503     & weighttype, lxxadxx, mythid)
504 heimbach 1.21 #endif
505    
506 heimbach 1.5 close ( cunit )
507 heimbach 1.1
508 heimbach 1.5 _END_MASTER( mythid )
509 heimbach 1.11
510     #endif /* EXCLUDE_CTRL_PACK */
511 heimbach 1.1
512     return
513     end
514    

  ViewVC Help
Powered by ViewVC 1.1.22