/[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.21 - (hide annotations) (download)
Thu Jul 28 13:47:49 2005 UTC (18 years, 9 months ago) by heimbach
Branch: MAIN
Changes since 1.20: +16 -6 lines
Adding precip control

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

  ViewVC Help
Powered by ViewVC 1.1.22