/[MITgcm]/MITgcm/pkg/ctrl/ctrl_init.F
ViewVC logotype

Annotation of /MITgcm/pkg/ctrl/ctrl_init.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.15 - (hide annotations) (download)
Mon Feb 28 17:29:38 2005 UTC (19 years, 2 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57f_post, checkpoint57e_post, checkpoint57f_pre
Changes since 1.14: +17 -1 lines
Adding eddy stress controls a la Ferreira et al.

1 edhill 1.10 C
2 heimbach 1.15 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_init.F,v 1.14 2005/02/24 05:08:34 heimbach Exp $
3 heimbach 1.14 C $Name: $
4 heimbach 1.1
5     #include "CTRL_CPPOPTIONS.h"
6    
7 heimbach 1.5 subroutine ctrl_init( mythid )
8 heimbach 1.1
9     c ==================================================================
10 heimbach 1.5 c SUBROUTINE ctrl_init
11 heimbach 1.1 c ==================================================================
12     c
13     c o Set parts of the vector of control variables and initialize the
14     c rest to zero.
15     c
16     c The vector of control variables is initialized here. The
17     c temperature and salinity contributions are read from file.
18     c Subsequently, the latter are dimensionalized and the tile
19     c edges are updated.
20     c
21     c started: Christian Eckert eckert@mit.edu 30-Jun-1999
22     c
23     c changed: Christian Eckert eckert@mit.edu 23-Feb-2000
24     c - Restructured the code in order to create a package
25     c for the MITgcmUV.
26     c
27     c Patrick Heimbach heimbach@mit.edu 30-May-2000
28     c - diffsec was falsely declared.
29     c
30     c Patrick Heimbach heimbach@mit.edu 06-Jun-2000
31     c - Transferred some filename declarations
32     c from ctrl_pack/ctrl_unpack to here
33     c - Transferred mask-per-tile to here
34     c - computation of control vector length here
35     c
36     c Patrick Heimbach heimbach@mit.edu 16-Jun-2000
37     c - Added call to ctrl_pack
38     c - Alternatively: transfer writing of scale files to
39     c ctrl_unpack
40     c
41 heimbach 1.9 c Dimitris Menemenlis menemenlis@mit.edu 7-Mar-2003
42     c - To be consistent with usage in ctrl_getrec.F,
43     c startrec and endrec need to be referenced to
44     c model time = 0, not to startTime.
45     c Also "- modelstep" -> "+ modelstep/2":
46     c old: startrec = int((modelstart - diffsecs)/
47     c old: & xx_???period) + 1
48     c old: endrec = int((modelend - diffsecs - modelstep)/
49     c old: & xx_???period) + 2
50     c new: startrec = int((modelstart + startTime - diffsecs)/
51     c new: & xx_???period) + 1
52     c new: endrec = int((modelend + startTime - diffsecs + modelstep/2)/
53     c new: & xx_???period) + 2
54     c
55 heimbach 1.11 c heimbach@mit.edu totally restructured 28-Oct-2003
56     c
57 heimbach 1.1 c ==================================================================
58 heimbach 1.5 c SUBROUTINE ctrl_init
59 heimbach 1.1 c ==================================================================
60    
61     implicit none
62    
63     c == global variables ==
64    
65     #include "EEPARAMS.h"
66     #include "SIZE.h"
67     #include "PARAMS.h"
68     #include "GRID.h"
69     #include "ctrl.h"
70 heimbach 1.12 #include "optim.h"
71 heimbach 1.1
72 heimbach 1.11 #ifdef ALLOW_CAL
73 heimbach 1.9 # include "cal.h"
74 heimbach 1.5 #endif
75     #ifdef ALLOW_OBCS_CONTROL
76     # include "OBCS.h"
77     #endif
78    
79 heimbach 1.1 c == routine arguments ==
80    
81     integer mythid
82    
83     c == local variables ==
84    
85     integer ntmp
86 heimbach 1.11 integer ivar
87 heimbach 1.1
88 heimbach 1.5 integer iobcs
89 heimbach 1.1 integer il
90     integer errio
91     integer startrec
92     integer endrec
93 heimbach 1.11 integer diffrec
94 heimbach 1.5 integer difftime(4)
95     _RL diffsecs
96 heimbach 1.1
97     character*(max_len_prec) record
98     character*(max_len_mbuf) msgbuf
99 heimbach 1.11 character*2 whichxyz
100 heimbach 1.5
101 heimbach 1.1 c == external ==
102    
103     integer ilnblnk
104     external ilnblnk
105    
106     c == end of interface ==
107    
108     c-- Set default values.
109 heimbach 1.11 do ivar = 1,maxcvars
110     ncvarindex(ivar) = -1
111     ncvarrecs(ivar) = 0
112     ncvarxmax(ivar) = 0
113     ncvarymax(ivar) = 0
114     ncvarnrmax(ivar) = 0
115     ncvargrd(ivar) = '?'
116     enddo
117 heimbach 1.1
118     _BARRIER
119    
120     c-- =====================
121     c-- Initial state fields.
122     c-- =====================
123    
124 heimbach 1.5 cph(
125     cph index 7-10 reserved for atmos. state,
126     cph index 11-14 reserved for open boundaries,
127     cph index 15-16 reserved for mixing coeff.
128     cph index 17 reserved for passive tracer TR1
129     cph index 18,19 reserved for sst, sss
130     cph index 20 for hFacC
131 heimbach 1.6 cph index 21-22 for efluxy, efluxp
132     cph index 23-24 for bottom drag
133 heimbach 1.5 cph)
134    
135 heimbach 1.9 c-------------------------------------------------------------------------------
136 heimbach 1.6 c--
137 heimbach 1.1 #ifdef ALLOW_THETA0_CONTROL
138 heimbach 1.5 c-- Initial state temperature contribution.
139 heimbach 1.11 call ctrl_init_ctrlvar (
140     & xx_theta_file, 1, 101, 1, 1, 1,
141     & snx, sny, nr, 'c', '3d', mythid )
142 heimbach 1.1 #endif /* ALLOW_THETA0_CONTROL */
143    
144 heimbach 1.9 c-------------------------------------------------------------------------------
145 heimbach 1.6 c--
146 heimbach 1.1 #ifdef ALLOW_SALT0_CONTROL
147 heimbach 1.5 c-- Initial state salinity contribution.
148 heimbach 1.11 call ctrl_init_ctrlvar (
149     & xx_salt_file, 2, 102, 1, 1, 1,
150     & snx, sny, nr, 'c', '3d', mythid )
151 heimbach 1.1 #endif /* ALLOW_SALT0_CONTROL */
152    
153 heimbach 1.5 c-- ===========================
154     c-- Surface flux contributions.
155     c-- ===========================
156    
157 heimbach 1.9 c-------------------------------------------------------------------------------
158 heimbach 1.6 c--
159 heimbach 1.5 #if (defined (ALLOW_HFLUX_CONTROL))
160     c-- Heat flux.
161    
162 heimbach 1.11 # ifdef ALLOW_CAL
163     call cal_FullDate( xx_hfluxstartdate1, xx_hfluxstartdate2,
164     & xx_hfluxstartdate , mythid )
165 heimbach 1.5 call cal_TimePassed( xx_hfluxstartdate, modelstartdate,
166     & difftime, mythid )
167     call cal_ToSeconds ( difftime, diffsecs, mythid )
168 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
169 heimbach 1.5 & xx_hfluxperiod) + 1
170 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
171 heimbach 1.5 & xx_hfluxperiod) + 2
172 heimbach 1.11 # else
173 heimbach 1.5 startrec = 1
174     endrec = 1
175 heimbach 1.11 # endif
176     diffrec = endrec - startrec + 1
177     call ctrl_init_ctrlvar (
178     & xx_hflux_file, 3, 103, diffrec, startrec, endrec,
179     & snx, sny, 1, 'c', 'xy', mythid )
180 heimbach 1.5
181     #elif (defined (ALLOW_ATEMP_CONTROL))
182     c-- Atmos. temperature
183    
184 heimbach 1.11 # ifdef ALLOW_CAL
185     call cal_FullDate( xx_atempstartdate1, xx_atempstartdate2,
186     & xx_atempstartdate , mythid )
187 heimbach 1.5 call cal_TimePassed( xx_atempstartdate, modelstartdate,
188     & difftime, mythid )
189     call cal_ToSeconds ( difftime, diffsecs, mythid )
190 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
191 heimbach 1.5 & xx_atempperiod) + 1
192 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
193 heimbach 1.5 & xx_atempperiod) + 2
194 heimbach 1.11 # else
195 heimbach 1.5 startrec = 1
196     endrec = 1
197 heimbach 1.11 # endif
198     diffrec = endrec - startrec + 1
199     call ctrl_init_ctrlvar (
200     & xx_atemp_file, 7, 107, diffrec, startrec, endrec,
201     & snx, sny, 1, 'c', 'xy', mythid )
202 heimbach 1.5
203     #elif (defined (ALLOW_HFLUX0_CONTROL))
204     c-- initial forcing only
205 heimbach 1.11 call ctrl_init_ctrlvar (
206     & xx_hflux_file, 3, 103, 1, 1, 1,
207     & snx, sny, 1, 'c', 'xy', mythid )
208 heimbach 1.1
209 heimbach 1.5 #endif /* ALLOW_HFLUX_CONTROL */
210    
211 heimbach 1.9 c-------------------------------------------------------------------------------
212 heimbach 1.6 c--
213 heimbach 1.5 #if (defined (ALLOW_SFLUX_CONTROL))
214     c-- Salt flux.
215    
216 heimbach 1.11 # ifdef ALLOW_CAL
217     call cal_FullDate( xx_sfluxstartdate1, xx_sfluxstartdate2,
218     & xx_sfluxstartdate , mythid )
219 heimbach 1.5 call cal_TimePassed( xx_sfluxstartdate, modelstartdate,
220     & difftime, mythid )
221     call cal_ToSeconds ( difftime, diffsecs, mythid )
222 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
223 heimbach 1.5 & xx_sfluxperiod) + 1
224 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
225 heimbach 1.5 & xx_sfluxperiod) + 2
226 heimbach 1.11 # else
227 heimbach 1.5 startrec = 1
228     endrec = 1
229 heimbach 1.11 # endif
230     diffrec = endrec - startrec + 1
231     call ctrl_init_ctrlvar (
232     & xx_sflux_file, 4, 104, diffrec, startrec, endrec,
233     & snx, sny, 1, 'c', 'xy', mythid )
234    
235 heimbach 1.5 #elif (defined (ALLOW_AQH_CONTROL))
236     c-- Atmos. humidity
237    
238 heimbach 1.11 # ifdef ALLOW_CAL
239     call cal_FullDate( xx_aqhstartdate1, xx_aqhstartdate2,
240     & xx_aqhstartdate , mythid )
241 heimbach 1.5 call cal_TimePassed( xx_aqhstartdate, modelstartdate,
242     & difftime, mythid )
243     call cal_ToSeconds ( difftime, diffsecs, mythid )
244 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
245 heimbach 1.5 & xx_aqhperiod) + 1
246 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
247 heimbach 1.5 & xx_aqhperiod) + 2
248 heimbach 1.11 # else
249 heimbach 1.5 startrec = 1
250     endrec = 1
251 heimbach 1.11 # endif
252     diffrec = endrec - startrec + 1
253     call ctrl_init_ctrlvar (
254     & xx_aqh_file, 8, 108, diffrec, startrec, endrec,
255     & snx, sny, 1, 'c', 'xy', mythid )
256 heimbach 1.5
257     #elif (defined (ALLOW_SFLUX0_CONTROL))
258     c-- initial forcing only
259 heimbach 1.11 call ctrl_init_ctrlvar (
260     & xx_sflux_file, 4, 104, 1, 1, 1,
261     & snx, sny, 1, 'c', 'xy', mythid )
262 heimbach 1.1
263 heimbach 1.5 #endif /* ALLOW_SFLUX_CONTROL */
264    
265 heimbach 1.9 c-------------------------------------------------------------------------------
266 heimbach 1.6 c--
267 heimbach 1.5 #if (defined (ALLOW_USTRESS_CONTROL))
268     c-- Zonal wind stress.
269    
270 heimbach 1.11 # ifdef ALLOW_CAL
271     call cal_FullDate( xx_tauustartdate1, xx_tauustartdate2,
272     & xx_tauustartdate, mythid )
273 heimbach 1.5 call cal_TimePassed( xx_tauustartdate, modelstartdate,
274     & difftime, mythid )
275     call cal_ToSeconds ( difftime, diffsecs, mythid )
276 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
277 heimbach 1.5 & xx_tauuperiod) + 1
278 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
279 heimbach 1.5 & xx_tauuperiod) + 2
280 heimbach 1.11 # else
281 heimbach 1.5 startrec = 1
282     endrec = 1
283 heimbach 1.11 # endif
284     diffrec = endrec - startrec + 1
285     call ctrl_init_ctrlvar (
286     & xx_tauu_file, 5, 105, diffrec, startrec, endrec,
287     & snx, sny, 1, 'w', 'xy', mythid )
288 heimbach 1.5
289     #elif (defined (ALLOW_UWIND_CONTROL))
290     c-- Zonal wind speed.
291    
292 heimbach 1.11 # ifdef ALLOW_CAL
293     call cal_FullDate( xx_uwindstartdate1, xx_uwindstartdate2,
294     & xx_uwindstartdate , mythid )
295 heimbach 1.5 call cal_TimePassed( xx_uwindstartdate, modelstartdate,
296     & difftime, mythid )
297     call cal_ToSeconds ( difftime, diffsecs, mythid )
298 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
299 heimbach 1.5 & xx_uwindperiod) + 1
300 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
301 heimbach 1.5 & xx_uwindperiod) + 2
302 heimbach 1.11 # else
303 heimbach 1.5 startrec = 1
304     endrec = 1
305 heimbach 1.11 # endif
306     diffrec = endrec - startrec + 1
307     call ctrl_init_ctrlvar (
308     & xx_uwind_file, 9, 109, diffrec, startrec, endrec,
309 heimbach 1.14 & snx, sny, 1, 'c', 'xy', mythid )
310 heimbach 1.5
311     #elif (defined (ALLOW_TAUU0_CONTROL))
312     c-- initial forcing only
313 heimbach 1.11 call ctrl_init_ctrlvar (
314     & xx_tauu_file, 5, 105, 1, 1, 1,
315     & snx, sny, 1, 'w', 'xy', mythid )
316 heimbach 1.1
317 heimbach 1.5 #endif /* ALLOW_USTRESS_CONTROL */
318    
319 heimbach 1.9 c-------------------------------------------------------------------------------
320 heimbach 1.6 c--
321 heimbach 1.5 #if (defined (ALLOW_VSTRESS_CONTROL))
322     c-- Meridional wind stress.
323    
324 heimbach 1.11 # ifdef ALLOW_CAL
325     call cal_FullDate( xx_tauvstartdate1, xx_tauvstartdate2,
326     & xx_tauvstartdate, mythid )
327 heimbach 1.5 call cal_TimePassed( xx_tauvstartdate, modelstartdate,
328     & difftime, mythid )
329     call cal_ToSeconds ( difftime, diffsecs, mythid )
330 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
331 heimbach 1.5 & xx_tauvperiod) + 1
332 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
333 heimbach 1.5 & xx_tauvperiod) + 2
334 heimbach 1.11 # else
335 heimbach 1.5 startrec = 1
336     endrec = 1
337 heimbach 1.11 # endif
338     diffrec = endrec - startrec + 1
339     call ctrl_init_ctrlvar (
340     & xx_tauv_file, 6, 106, diffrec, startrec, endrec,
341     & snx, sny, 1, 's', 'xy', mythid )
342 heimbach 1.5
343     #elif (defined (ALLOW_VWIND_CONTROL))
344     c-- Meridional wind speed.
345    
346 heimbach 1.11 # ifdef ALLOW_CAL
347     call cal_FullDate( xx_vwindstartdate1, xx_vwindstartdate2,
348     & xx_vwindstartdate , mythid )
349 heimbach 1.5 call cal_TimePassed( xx_vwindstartdate, modelstartdate,
350     & difftime, mythid )
351     call cal_ToSeconds ( difftime, diffsecs, mythid )
352 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
353 heimbach 1.5 & xx_vwindperiod) + 1
354 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
355 heimbach 1.5 & xx_vwindperiod) + 2
356 heimbach 1.11 # else
357 heimbach 1.5 startrec = 1
358     endrec = 1
359 heimbach 1.11 # endif
360     diffrec = endrec - startrec + 1
361     call ctrl_init_ctrlvar (
362     & xx_vwind_file, 10, 110, diffrec, startrec, endrec,
363 heimbach 1.14 & snx, sny, 1, 'c', 'xy', mythid )
364 heimbach 1.5
365     #elif (defined (ALLOW_TAUV0_CONTROL))
366     c-- initial forcing only
367 heimbach 1.11 diffrec = endrec - startrec + 1
368     call ctrl_init_ctrlvar (
369     & xx_tauv_file, 6, 106, 1, 1, 1,
370     & snx, sny, 1, 's', 'xy', mythid )
371 heimbach 1.1
372 heimbach 1.5 #endif /* ALLOW_VSTRESS_CONTROL */
373    
374 heimbach 1.11 c-- ===========================
375     c-- Open boundary contributions.
376     c-- ===========================
377    
378 heimbach 1.9 c-------------------------------------------------------------------------------
379 heimbach 1.6 c--
380 heimbach 1.5 #ifdef ALLOW_OBCSN_CONTROL
381     c-- Northern obc.
382    
383 heimbach 1.11 # ifdef ALLOW_CAL
384     call cal_FullDate( xx_obcsnstartdate1, xx_obcsnstartdate2,
385     & xx_obcsnstartdate, mythid )
386 heimbach 1.5 call cal_TimePassed( xx_obcsnstartdate, modelstartdate,
387     & difftime, mythid )
388     call cal_ToSeconds ( difftime, diffsecs, mythid )
389 heimbach 1.9 startrec = int((modelstart - diffsecs)/xx_obcsnperiod) + 1
390 heimbach 1.11 startrec = (startrec - 1)*nobcs + 1
391 heimbach 1.9 endrec = int((modelend - diffsecs)/xx_obcsnperiod) + 2
392 heimbach 1.11 endrec = (endrec - startrec + 1)*nobcs
393     # else
394 heimbach 1.5 startrec = 1
395     endrec = 1
396 heimbach 1.11 # endif
397 mlosch 1.13 diffrec = endrec
398 heimbach 1.11 call ctrl_init_ctrlvar (
399     & xx_obcsn_file, 11, 111, diffrec, startrec, endrec,
400     & snx, sny, nr, 'm', 'xz', mythid )
401 heimbach 1.5
402     #endif /* ALLOW_OBCSN_CONTROL */
403    
404 heimbach 1.11 c-------------------------------------------------------------------------------
405     c--
406 heimbach 1.5 #ifdef ALLOW_OBCSS_CONTROL
407     c-- Southern obc.
408    
409 heimbach 1.11 # ifdef ALLOW_CAL
410     call cal_FullDate( xx_obcssstartdate1, xx_obcssstartdate2,
411     & xx_obcssstartdate, mythid )
412 heimbach 1.5 call cal_TimePassed( xx_obcssstartdate, modelstartdate,
413     & difftime, mythid )
414     call cal_ToSeconds ( difftime, diffsecs, mythid )
415 heimbach 1.9 startrec = int((modelstart - diffsecs)/xx_obcssperiod) + 1
416 heimbach 1.11 startrec = (startrec - 1)*nobcs + 1
417 heimbach 1.9 endrec = int((modelend - diffsecs)/xx_obcssperiod) + 2
418 heimbach 1.11 endrec = (endrec - startrec + 1)*nobcs
419     # else
420 heimbach 1.5 startrec = 1
421     endrec = 1
422 heimbach 1.11 # endif
423 mlosch 1.13 diffrec = endrec
424 heimbach 1.11 call ctrl_init_ctrlvar (
425     & xx_obcss_file, 12, 112, diffrec, startrec, endrec,
426     & snx, sny, nr, 'm', 'xz', mythid )
427 heimbach 1.5
428     #endif /* ALLOW_OBCSS_CONTROL */
429    
430 heimbach 1.9 c-------------------------------------------------------------------------------
431 heimbach 1.6 c--
432 heimbach 1.5 #ifdef ALLOW_OBCSW_CONTROL
433     c-- Western obc.
434    
435 heimbach 1.11 # ifdef ALLOW_CAL
436     call cal_FullDate( xx_obcswstartdate1, xx_obcswstartdate2,
437     & xx_obcswstartdate, mythid )
438 heimbach 1.5 call cal_TimePassed( xx_obcswstartdate, modelstartdate,
439     & difftime, mythid )
440     call cal_ToSeconds ( difftime, diffsecs, mythid )
441 heimbach 1.9 startrec = int((modelstart - diffsecs)/xx_obcswperiod) + 1
442 heimbach 1.11 startrec = (startrec - 1)*nobcs + 1
443 heimbach 1.9 endrec = int((modelend - diffsecs)/xx_obcswperiod) + 2
444 heimbach 1.11 endrec = (endrec - startrec + 1)*nobcs
445     # else
446 heimbach 1.5 startrec = 1
447     endrec = 1
448 heimbach 1.11 # endif
449 mlosch 1.13 diffrec = endrec
450 heimbach 1.11 call ctrl_init_ctrlvar (
451     & xx_obcsw_file, 13, 113, diffrec, startrec, endrec,
452     & snx, sny, nr, 'm', 'yz', mythid )
453 heimbach 1.5
454 heimbach 1.11 #endif /* ALLOW_OBCSW_CONTROL */
455 heimbach 1.5
456 heimbach 1.9 c-------------------------------------------------------------------------------
457 heimbach 1.6 c--
458 heimbach 1.5 #ifdef ALLOW_OBCSE_CONTROL
459     c-- Eastern obc.
460    
461 heimbach 1.11 # ifdef ALLOW_CAL
462     call cal_FullDate( xx_obcsestartdate1, xx_obcsestartdate2,
463     & xx_obcsestartdate, mythid )
464 heimbach 1.5 call cal_TimePassed( xx_obcsestartdate, modelstartdate,
465     & difftime, mythid )
466     call cal_ToSeconds ( difftime, diffsecs, mythid )
467 heimbach 1.9 startrec = int((modelstart - diffsecs)/xx_obcseperiod) + 1
468 heimbach 1.11 startrec = (startrec - 1)*nobcs + 1
469 heimbach 1.9 endrec = int((modelend - diffsecs)/xx_obcseperiod) + 2
470 heimbach 1.11 endrec = (endrec - startrec + 1)*nobcs
471     # else
472 heimbach 1.5 startrec = 1
473     endrec = 1
474 heimbach 1.11 # endif
475 mlosch 1.13 diffrec = endrec
476 heimbach 1.11 call ctrl_init_ctrlvar (
477     & xx_obcse_file, 14, 114, diffrec, startrec, endrec,
478     & snx, sny, nr, 'm', 'yz', mythid )
479 heimbach 1.5
480     #endif /* ALLOW_OBCSE_CONTROL */
481 heimbach 1.2
482 heimbach 1.9 c-------------------------------------------------------------------------------
483 heimbach 1.6 c--
484 heimbach 1.3 #ifdef ALLOW_DIFFKR_CONTROL
485 heimbach 1.11 call ctrl_init_ctrlvar (
486     & xx_diffkr_file, 15, 115, 1, 1, 1,
487     & snx, sny, nr, 'c', '3d', mythid )
488 heimbach 1.3 #endif /* ALLOW_DIFFKR_CONTROL */
489    
490 heimbach 1.9 c-------------------------------------------------------------------------------
491 heimbach 1.6 c--
492 heimbach 1.3 #ifdef ALLOW_KAPGM_CONTROL
493 heimbach 1.11 call ctrl_init_ctrlvar (
494     & xx_kapgm_file, 16, 116, 1, 1, 1,
495     & snx, sny, nr, 'c', '3d', mythid )
496 heimbach 1.3 #endif /* ALLOW_KAPGM_CONTROL */
497    
498 heimbach 1.9 c-------------------------------------------------------------------------------
499 heimbach 1.6 c--
500 heimbach 1.2 #ifdef ALLOW_TR10_CONTROL
501 heimbach 1.11 call ctrl_init_ctrlvar (
502     & xx_tr1_file, 17, 117, 1, 1, 1,
503     & snx, sny, nr, 'c', '3d', mythid )
504 heimbach 1.2 #endif /* ALLOW_TR10_CONTROL */
505    
506 heimbach 1.9 c-------------------------------------------------------------------------------
507 heimbach 1.6 c--
508     #ifdef ALLOW_SST0_CONTROL
509 heimbach 1.11 call ctrl_init_ctrlvar (
510     & xx_sst_file, 18, 118, 1, 1, 1,
511     & snx, sny, 1, 'c', 'xy', mythid )
512 heimbach 1.6 #endif /* ALLOW_SST0_CONTROL */
513    
514 heimbach 1.9 c-------------------------------------------------------------------------------
515 heimbach 1.6 c--
516     #ifdef ALLOW_SSS0_CONTROL
517 heimbach 1.11 call ctrl_init_ctrlvar (
518     & xx_sss_file, 19, 119, 1, 1, 1,
519     & snx, sny, 1, 'c', 'xy', mythid )
520 heimbach 1.6 #endif /* ALLOW_SSS0_CONTROL */
521    
522 heimbach 1.9 c-------------------------------------------------------------------------------
523 heimbach 1.6 c--
524     #ifdef ALLOW_HFACC_CONTROL
525 heimbach 1.11 # ifdef ALLOW_HFACC3D_CONTROL
526     call ctrl_init_ctrlvar (
527     & xx_hfacc_file, 20, 120, 1, 1, 1,
528     & snx, sny, nr, 'c', '3d', mythid )
529     # else
530     call ctrl_init_ctrlvar (
531     & xx_hfacc_file, 20, 120, 1, 1, 1,
532     & snx, sny, 1, 'c', 'xy', mythid )
533     # endif /*ALLOW_HFACC3D_CONTROL*/
534 heimbach 1.6 #endif /* ALLOW_HFACC_CONTROL */
535    
536 heimbach 1.9 c-------------------------------------------------------------------------------
537 heimbach 1.6 c--
538 heimbach 1.5 #ifdef ALLOW_EFLUXY0_CONTROL
539 heimbach 1.11 call ctrl_init_ctrlvar (
540     & xx_efluxy_file, 21, 121, 1, 1, 1,
541     & snx, sny, nr, 's', '3d', mythid )
542 heimbach 1.5 #endif /* ALLOW_EFLUXY0_CONTROL */
543    
544 heimbach 1.9 c-------------------------------------------------------------------------------
545 heimbach 1.6 c--
546 heimbach 1.5 #ifdef ALLOW_EFLUXP0_CONTROL
547 heimbach 1.11 call ctrl_init_ctrlvar (
548     & xx_efluxp_file, 22, 122, 1, 1, 1,
549     & snx, sny, nr, 'v', '3d', mythid )
550 heimbach 1.5 #endif /* ALLOW_EFLUXP0_CONTROL */
551 heimbach 1.6
552 heimbach 1.9 c-------------------------------------------------------------------------------
553 heimbach 1.6 c--
554     #ifdef ALLOW_BOTTOMDRAG_CONTROL
555 heimbach 1.11 call ctrl_init_ctrlvar (
556     & xx_bottomdrag_file, 23, 123, 1, 1, 1,
557     & snx, sny, 1, 'c', 'xy', mythid )
558 heimbach 1.6 #endif /* ALLOW_BOTTOMDRAG_CONTROL */
559    
560 heimbach 1.9 c-------------------------------------------------------------------------------
561 heimbach 1.15 c--
562     #ifdef ALLOW_EDTAUX_CONTROL
563     call ctrl_init_ctrlvar (
564     & xx_edtaux_file, 25, 125, 1, 1, 1,
565     & snx, sny, nr, 'w', '3d', mythid )
566     #endif /* ALLOW_EDTAUX_CONTROL */
567    
568     c-------------------------------------------------------------------------------
569     c--
570     #ifdef ALLOW_EDTAUY_CONTROL
571     call ctrl_init_ctrlvar (
572     & xx_edtauy_file, 26, 126, 1, 1, 1,
573     & snx, sny, nr, 's', '3d', mythid )
574     #endif /* ALLOW_EDTAUY_CONTROL */
575    
576     c-------------------------------------------------------------------------------
577 heimbach 1.9 c-------------------------------------------------------------------------------
578     c-------------------------------------------------------------------------------
579 heimbach 1.1
580 heimbach 1.11 call ctrl_init_wet( mythid )
581    
582 heimbach 1.1 return
583     end
584    

  ViewVC Help
Powered by ViewVC 1.1.22