/[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.20 - (hide annotations) (download)
Tue Apr 12 22:22:21 2005 UTC (19 years, 1 month ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57m_post, checkpoint57g_pre, checkpoint57g_post, checkpoint57i_post, checkpoint57l_post, checkpoint57h_pre, checkpoint57h_post, checkpoint57h_done, checkpoint57j_post, checkpoint57k_post
Changes since 1.19: +3 -3 lines
nwetcglobal was changed, but not mask.

1 heimbach 1.20 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_pack.F,v 1.19 2005/04/07 23:38:43 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     call ctrl_set_fname(xx_uwind_file, fname_uwind, mythid)
116     call ctrl_set_fname(xx_vwind_file, fname_vwind, mythid)
117     call ctrl_set_fname(xx_obcsn_file, fname_obcsn, mythid)
118     call ctrl_set_fname(xx_obcss_file, fname_obcss, mythid)
119     call ctrl_set_fname(xx_obcsw_file, fname_obcsw, mythid)
120     call ctrl_set_fname(xx_obcse_file, fname_obcse, mythid)
121     call ctrl_set_fname(xx_diffkr_file, fname_diffkr, mythid)
122     call ctrl_set_fname(xx_kapgm_file, fname_kapgm, mythid)
123     call ctrl_set_fname(xx_tr1_file, fname_tr1, mythid)
124     call ctrl_set_fname(xx_sst_file, fname_sst, mythid)
125     call ctrl_set_fname(xx_sss_file, fname_sss, mythid)
126     call ctrl_set_fname(xx_hfacc_file, fname_hfacc, mythid)
127     call ctrl_set_fname(xx_efluxy_file, fname_efluxy, mythid)
128     call ctrl_set_fname(xx_efluxp_file, fname_efluxp, mythid)
129     call ctrl_set_fname(xx_bottomdrag_file, fname_bottomdrag, mythid)
130 heimbach 1.18 call ctrl_set_fname(xx_edtaux_file, fname_edtaux, mythid)
131     call ctrl_set_fname(xx_edtauy_file, fname_edtauy, mythid)
132 heimbach 1.19 call ctrl_set_fname(xx_uvel_file, fname_uvel, mythid)
133     call ctrl_set_fname(xx_vvel_file, fname_vvel, mythid)
134     call ctrl_set_fname(xx_etan_file, fname_etan, mythid)
135     call ctrl_set_fname(xx_relaxsst_file, fname_relaxsst, mythid)
136     call ctrl_set_fname(xx_relaxsss_file, fname_relaxsss, mythid)
137 heimbach 1.5
138 heimbach 1.19 c-- Only the master thread will do I/O.
139 heimbach 1.1 _BEGIN_MASTER( mythid )
140    
141 heimbach 1.13 if ( first ) then
142 heimbach 1.11 c >>> Initialise control vector for optimcycle=0 <<<
143     lxxadxx = .TRUE.
144     ictrlgrad = 1
145     fcloc = fmin
146     write(cfile(1:128),'(4a,i4.4)')
147 heimbach 1.13 & ctrlname(1:9),'_',yctrlid(1:10),
148     & yctrlpospack, optimcycle
149 heimbach 1.17 print *, 'ph-pack: packing ', ctrlname(1:9)
150 heimbach 1.11 else
151 heimbach 1.1 c >>> Write gradient vector <<<
152 heimbach 1.11 lxxadxx = .FALSE.
153     ictrlgrad = 2
154     fcloc = fc
155 heimbach 1.5 write(cfile(1:128),'(4a,i4.4)')
156 heimbach 1.13 & costname(1:9),'_',yctrlid(1:10),
157     & yctrlpospack, optimcycle
158 heimbach 1.17 print *, 'ph-pack: packing ', costname(1:9)
159 heimbach 1.11 endif
160 heimbach 1.1
161 heimbach 1.11 call mdsfindunit( cunit, mythid )
162     open( cunit, file = cfile,
163     & status = 'unknown',
164     & form = 'unformatted',
165     & access = 'sequential' )
166 heimbach 1.1
167     c-- Header information.
168 mlosch 1.15 write(cunit) nvartype
169     write(cunit) nvarlength
170     write(cunit) yctrlid
171     write(cunit) optimCycle
172     write(cunit) fc
173     C place holder of obsolete variable iG
174     write(cunit) 1
175     C place holder of obsolete variable jG
176     write(cunit) 1
177     write(cunit) nsx
178     write(cunit) nsy
179     write(cunit) (nWetcGlobal(k), k=1,nr)
180     write(cunit) (nWetsGlobal(k), k=1,nr)
181     write(cunit) (nWetwGlobal(k), k=1,nr)
182 heimbach 1.7 #ifdef ALLOW_CTRL_WETV
183 mlosch 1.15 write(cunit) (nWetvGlobal(k), k=1,nr)
184 heimbach 1.7 #endif
185 heimbach 1.11
186 heimbach 1.5 #ifdef ALLOW_OBCSN_CONTROL
187     write(cunit) ((nWetobcsnGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
188     #endif
189     #ifdef ALLOW_OBCSS_CONTROL
190     write(cunit) ((nWetobcssGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
191     #endif
192     #ifdef ALLOW_OBCSW_CONTROL
193     write(cunit) ((nWetobcswGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
194     #endif
195     #ifdef ALLOW_OBCSE_CONTROL
196     write(cunit) ((nWetobcseGlo(k,iobcs), k=1,nr),iobcs= 1,nobcs)
197     #endif
198 mlosch 1.15 write(cunit) (ncvarindex(i), i=1,maxcvars)
199     write(cunit) (ncvarrecs(i), i=1,maxcvars)
200     write(cunit) (ncvarxmax(i), i=1,maxcvars)
201     write(cunit) (ncvarymax(i), i=1,maxcvars)
202     write(cunit) (ncvarnrmax(i), i=1,maxcvars)
203     write(cunit) (ncvargrd(i), i=1,maxcvars)
204 heimbach 1.1 write(cunit)
205    
206     #ifdef ALLOW_THETA0_CONTROL
207 heimbach 1.5 ivartype = 1
208 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
209     write(weighttype(1:80),'(a)') "wtheta"
210 heimbach 1.5 call ctrl_set_pack_xyz(
211 heimbach 1.19 & cunit, ivartype, fname_theta(ictrlgrad), "maskCtrlC",
212 heimbach 1.8 & weighttype, wtheta, lxxadxx, mythid)
213 heimbach 1.1 #endif
214    
215     #ifdef ALLOW_SALT0_CONTROL
216 heimbach 1.5 ivartype = 2
217 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
218     write(weighttype(1:80),'(a)') "wsalt"
219 heimbach 1.5 call ctrl_set_pack_xyz(
220 heimbach 1.19 & cunit, ivartype, fname_salt(ictrlgrad), "maskCtrlC",
221 heimbach 1.8 & weighttype, wsalt, lxxadxx, mythid)
222 heimbach 1.5 #endif
223    
224     #if (defined (ALLOW_HFLUX_CONTROL) || \
225     defined (ALLOW_HFLUX0_CONTROL))
226     ivartype = 3
227     write(weighttype(1:80),'(80a)') ' '
228     write(weighttype(1:80),'(a)') "whflux"
229     call ctrl_set_pack_xy(
230 heimbach 1.19 & cunit, ivartype, fname_hflux(ictrlgrad), "maskCtrlC",
231 heimbach 1.11 & weighttype, lxxadxx, mythid)
232 heimbach 1.5 #endif
233    
234     #if (defined (ALLOW_SFLUX_CONTROL) || \
235     defined (ALLOW_SFLUX0_CONTROL))
236     ivartype = 4
237     write(weighttype(1:80),'(80a)') ' '
238     write(weighttype(1:80),'(a)') "wsflux"
239     call ctrl_set_pack_xy(
240 heimbach 1.19 & cunit, ivartype, fname_sflux(ictrlgrad), "maskCtrlC",
241 heimbach 1.11 & weighttype, lxxadxx, mythid)
242 heimbach 1.5 #endif
243    
244     #if (defined (ALLOW_USTRESS_CONTROL) || \
245     defined (ALLOW_TAUU0_CONTROL))
246     ivartype = 5
247     write(weighttype(1:80),'(80a)') ' '
248     write(weighttype(1:80),'(a)') "wtauu"
249     call ctrl_set_pack_xy(
250 heimbach 1.19 & cunit, ivartype, fname_tauu(ictrlgrad), "maskCtrlW",
251 heimbach 1.11 & weighttype, lxxadxx, mythid)
252 heimbach 1.5 #endif
253    
254     #if (defined (ALLOW_VSTRESS_CONTROL) || \
255     defined (ALLOW_TAUV0_CONTROL))
256     ivartype = 6
257     write(weighttype(1:80),'(80a)') ' '
258     write(weighttype(1:80),'(a)') "wtauv"
259     call ctrl_set_pack_xy(
260 heimbach 1.19 & cunit, ivartype, fname_tauv(ictrlgrad), "maskCtrlS",
261 heimbach 1.11 & weighttype, lxxadxx, mythid)
262 heimbach 1.5 #endif
263    
264     #ifdef ALLOW_ATEMP_CONTROL
265     ivartype = 7
266     write(weighttype(1:80),'(80a)') ' '
267     write(weighttype(1:80),'(a)') "watemp"
268     call ctrl_set_pack_xy(
269 heimbach 1.19 & cunit, ivartype, fname_atemp(ictrlgrad), "maskCtrlC",
270 heimbach 1.11 & weighttype, lxxadxx, mythid)
271 heimbach 1.5 #endif
272    
273     #ifdef ALLOW_AQH_CONTROL
274     ivartype = 8
275     write(weighttype(1:80),'(80a)') ' '
276     write(weighttype(1:80),'(a)') "waqh"
277     call ctrl_set_pack_xy(
278 heimbach 1.19 & cunit, ivartype, fname_aqh(ictrlgrad), "maskCtrlC",
279 heimbach 1.11 & weighttype, lxxadxx, mythid)
280 heimbach 1.5 #endif
281    
282     #ifdef ALLOW_UWIND_CONTROL
283     ivartype = 9
284     write(weighttype(1:80),'(80a)') ' '
285     write(weighttype(1:80),'(a)') "wuwind"
286     call ctrl_set_pack_xy(
287 heimbach 1.20 & cunit, ivartype, fname_uwind(ictrlgrad), "maskCtrlC",
288 heimbach 1.11 & weighttype, lxxadxx, mythid)
289 heimbach 1.5 #endif
290    
291     #ifdef ALLOW_VWIND_CONTROL
292     ivartype = 10
293     write(weighttype(1:80),'(80a)') ' '
294     write(weighttype(1:80),'(a)') "wvwind"
295     call ctrl_set_pack_xy(
296 heimbach 1.20 & cunit, ivartype, fname_vwind(ictrlgrad), "maskCtrlC",
297 heimbach 1.11 & weighttype, lxxadxx, mythid)
298 heimbach 1.5 #endif
299    
300     #ifdef ALLOW_OBCSN_CONTROL
301     ivartype = 11
302 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
303     write(weighttype(1:80),'(a)') "wobcsn"
304 heimbach 1.5 call ctrl_set_pack_xz(
305 heimbach 1.11 & cunit, ivartype, fname_obcsn(ictrlgrad), "maskobcsn",
306 heimbach 1.8 & weighttype, wobcsn, lxxadxx, mythid)
307 heimbach 1.5 #endif
308    
309     #ifdef ALLOW_OBCSS_CONTROL
310     ivartype = 12
311 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
312     write(weighttype(1:80),'(a)') "wobcss"
313 heimbach 1.5 call ctrl_set_pack_xz(
314 heimbach 1.11 & cunit, ivartype, fname_obcss(ictrlgrad), "maskobcss",
315 heimbach 1.8 & weighttype, wobcss, lxxadxx, mythid)
316 heimbach 1.5 #endif
317    
318     #ifdef ALLOW_OBCSW_CONTROL
319     ivartype = 13
320 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
321     write(weighttype(1:80),'(a)') "wobcsw"
322 heimbach 1.5 call ctrl_set_pack_yz(
323 heimbach 1.11 & cunit, ivartype, fname_obcsw(ictrlgrad), "maskobcsw",
324 heimbach 1.8 & weighttype, wobcsw, lxxadxx, mythid)
325 heimbach 1.5 #endif
326    
327     #ifdef ALLOW_OBCSE_CONTROL
328     ivartype = 14
329 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
330     write(weighttype(1:80),'(a)') "wobcse"
331 heimbach 1.5 call ctrl_set_pack_yz(
332 heimbach 1.11 & cunit, ivartype, fname_obcse(ictrlgrad), "maskobcse",
333 heimbach 1.8 & weighttype, wobcse, lxxadxx, mythid)
334 heimbach 1.1 #endif
335 heimbach 1.3
336     #ifdef ALLOW_DIFFKR_CONTROL
337 heimbach 1.5 ivartype = 15
338 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
339     write(weighttype(1:80),'(a)') "wdiffkr"
340 heimbach 1.5 call ctrl_set_pack_xyz(
341 heimbach 1.19 & cunit, ivartype, fname_diffkr(ictrlgrad), "maskCtrlC",
342 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
343 heimbach 1.3 #endif
344    
345     #ifdef ALLOW_KAPGM_CONTROL
346 heimbach 1.5 ivartype = 16
347 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
348     write(weighttype(1:80),'(a)') "wkapgm"
349 heimbach 1.5 call ctrl_set_pack_xyz(
350 heimbach 1.19 & cunit, ivartype, fname_kapgm(ictrlgrad), "maskCtrlC",
351 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
352 heimbach 1.3 #endif
353    
354 heimbach 1.5 #ifdef ALLOW_TR10_CONTROL
355     ivartype = 17
356 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
357     write(weighttype(1:80),'(a)') "wtr1"
358 heimbach 1.5 call ctrl_set_pack_xyz(
359 heimbach 1.19 & cunit, ivartype, fname_tr1(ictrlgrad), "maskCtrlC",
360 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
361 heimbach 1.5 #endif
362    
363 heimbach 1.6 #ifdef ALLOW_SST0_CONTROL
364     ivartype = 18
365     write(weighttype(1:80),'(80a)') ' '
366     write(weighttype(1:80),'(a)') "wsst0"
367     call ctrl_set_pack_xy(
368 heimbach 1.19 & cunit, ivartype, fname_sst(ictrlgrad), "maskCtrlC",
369 heimbach 1.11 & weighttype, lxxadxx, mythid)
370 heimbach 1.6 #endif
371    
372     #ifdef ALLOW_SSS0_CONTROL
373     ivartype = 19
374     write(weighttype(1:80),'(80a)') ' '
375     write(weighttype(1:80),'(a)') "wsss0"
376     call ctrl_set_pack_xy(
377 heimbach 1.19 & cunit, ivartype, fname_sss(ictrlgrad), "maskCtrlC",
378 heimbach 1.11 & weighttype, lxxadxx, mythid)
379 heimbach 1.6 #endif
380    
381     #ifdef ALLOW_HFACC_CONTROL
382     ivartype = 20
383 heimbach 1.8 write(weighttype(1:80),'(80a)') ' '
384     write(weighttype(1:80),'(a)') "whfacc"
385     # ifdef ALLOW_HFACC3D_CONTROL
386 heimbach 1.6 call ctrl_set_pack_xyz(
387 heimbach 1.19 & cunit, ivartype, fname_hfacc(ictrlgrad), "maskCtrlC",
388 heimbach 1.8 & weighttype, wunit, lxxadxx, mythid)
389     # else
390 heimbach 1.6 call ctrl_set_pack_xy(
391 heimbach 1.19 & cunit, ivartype, fname_hfacc(ictrlgrad), "maskCtrlC",
392 heimbach 1.11 & weighttype, lxxadxx, mythid)
393 heimbach 1.8 # endif
394 heimbach 1.6 #endif
395    
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.18 & weighttype, wunit, lxxadxx, mythid)
430     #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     & weighttype, wunit, lxxadxx, mythid)
439     #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     & cunit, ivartype, fname_etan(ictrlgrad), "maskCtrlC",
465     & weighttype, lxxadxx, mythid)
466     #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     & cunit, ivartype, fname_relaxsst(ictrlgrad), "maskCtrlC",
474     & weighttype, lxxadxx, mythid)
475     #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     & cunit, ivartype, fname_relaxsss(ictrlgrad), "maskCtrlC",
483     & weighttype, lxxadxx, mythid)
484     #endif
485 heimbach 1.18
486 heimbach 1.5 close ( cunit )
487 heimbach 1.1
488 heimbach 1.5 _END_MASTER( mythid )
489 heimbach 1.11
490     #endif /* EXCLUDE_CTRL_PACK */
491 heimbach 1.1
492     return
493     end
494    

  ViewVC Help
Powered by ViewVC 1.1.22