/[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.38 - (hide annotations) (download)
Tue May 10 07:30:15 2011 UTC (13 years ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint62x
Changes since 1.37: +15 -1 lines
add new control variable xx_shifwflx (fresh water flux underneath ice
shelves). This is almost as tedious as obcs-ctrl, because the
variables needs its own mask.

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

  ViewVC Help
Powered by ViewVC 1.1.22