/[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.28 - (hide annotations) (download)
Wed Jun 7 01:55:14 2006 UTC (17 years, 11 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint58m_post, checkpoint58o_post, checkpoint58p_post, checkpoint58q_post, checkpoint58n_post, checkpoint58k_post, checkpoint58l_post, checkpoint58h_post, checkpoint58j_post, checkpoint58i_post
Changes since 1.27: +6 -12 lines
Modifications for bottom topography control
o replace hFacC by _hFacC at various places
o replace ALLOW_HFACC_CONTROL by ALLOW_DEPTH_CONTROL
o add non-self-adjoint cg2d_nsa
o update autodiff support routines
o re-initialise hfac after ctrl_depth_ini
o works for 5x5 box, doesnt work for global_ocean.90x40x15

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