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

1 edhill 1.10 C
2 heimbach 1.17 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_init.F,v 1.16 2005/04/07 23:38:43 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 heimbach 1.17 cph index 17 reserved for precip. (atmos. state)
129 heimbach 1.5 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 heimbach 1.17 cph index 23 for bottom drag
133     cph index 24
134     cph index 25-26 for edtaux, edtauy
135     cph index 27-29 for uvel0, vvel0, etan0
136     cph index 30-31 for relax. SST, SSS
137     cph index 32 reserved for passive tracer TR1
138    
139 heimbach 1.5 cph)
140    
141 heimbach 1.16 c----------------------------------------------------------------------
142 heimbach 1.6 c--
143 heimbach 1.1 #ifdef ALLOW_THETA0_CONTROL
144 heimbach 1.5 c-- Initial state temperature contribution.
145 heimbach 1.11 call ctrl_init_ctrlvar (
146     & xx_theta_file, 1, 101, 1, 1, 1,
147     & snx, sny, nr, 'c', '3d', mythid )
148 heimbach 1.1 #endif /* ALLOW_THETA0_CONTROL */
149    
150 heimbach 1.16 c----------------------------------------------------------------------
151 heimbach 1.6 c--
152 heimbach 1.1 #ifdef ALLOW_SALT0_CONTROL
153 heimbach 1.5 c-- Initial state salinity contribution.
154 heimbach 1.11 call ctrl_init_ctrlvar (
155     & xx_salt_file, 2, 102, 1, 1, 1,
156     & snx, sny, nr, 'c', '3d', mythid )
157 heimbach 1.1 #endif /* ALLOW_SALT0_CONTROL */
158    
159 heimbach 1.5 c-- ===========================
160     c-- Surface flux contributions.
161     c-- ===========================
162    
163 heimbach 1.16 c----------------------------------------------------------------------
164 heimbach 1.6 c--
165 heimbach 1.5 #if (defined (ALLOW_HFLUX_CONTROL))
166     c-- Heat flux.
167    
168 heimbach 1.11 # ifdef ALLOW_CAL
169     call cal_FullDate( xx_hfluxstartdate1, xx_hfluxstartdate2,
170     & xx_hfluxstartdate , mythid )
171 heimbach 1.5 call cal_TimePassed( xx_hfluxstartdate, modelstartdate,
172     & difftime, mythid )
173     call cal_ToSeconds ( difftime, diffsecs, mythid )
174 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
175 heimbach 1.5 & xx_hfluxperiod) + 1
176 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
177 heimbach 1.5 & xx_hfluxperiod) + 2
178 heimbach 1.11 # else
179 heimbach 1.5 startrec = 1
180     endrec = 1
181 heimbach 1.11 # endif
182     diffrec = endrec - startrec + 1
183     call ctrl_init_ctrlvar (
184     & xx_hflux_file, 3, 103, diffrec, startrec, endrec,
185     & snx, sny, 1, 'c', 'xy', mythid )
186 heimbach 1.5
187     #elif (defined (ALLOW_ATEMP_CONTROL))
188     c-- Atmos. temperature
189    
190 heimbach 1.11 # ifdef ALLOW_CAL
191     call cal_FullDate( xx_atempstartdate1, xx_atempstartdate2,
192     & xx_atempstartdate , mythid )
193 heimbach 1.5 call cal_TimePassed( xx_atempstartdate, modelstartdate,
194     & difftime, mythid )
195     call cal_ToSeconds ( difftime, diffsecs, mythid )
196 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
197 heimbach 1.5 & xx_atempperiod) + 1
198 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
199 heimbach 1.5 & xx_atempperiod) + 2
200 heimbach 1.11 # else
201 heimbach 1.5 startrec = 1
202     endrec = 1
203 heimbach 1.11 # endif
204     diffrec = endrec - startrec + 1
205     call ctrl_init_ctrlvar (
206     & xx_atemp_file, 7, 107, diffrec, startrec, endrec,
207     & snx, sny, 1, 'c', 'xy', mythid )
208 heimbach 1.5
209     #elif (defined (ALLOW_HFLUX0_CONTROL))
210     c-- initial forcing only
211 heimbach 1.11 call ctrl_init_ctrlvar (
212     & xx_hflux_file, 3, 103, 1, 1, 1,
213     & snx, sny, 1, 'c', 'xy', mythid )
214 heimbach 1.1
215 heimbach 1.5 #endif /* ALLOW_HFLUX_CONTROL */
216    
217 heimbach 1.16 c----------------------------------------------------------------------
218 heimbach 1.6 c--
219 heimbach 1.5 #if (defined (ALLOW_SFLUX_CONTROL))
220     c-- Salt flux.
221    
222 heimbach 1.11 # ifdef ALLOW_CAL
223     call cal_FullDate( xx_sfluxstartdate1, xx_sfluxstartdate2,
224     & xx_sfluxstartdate , mythid )
225 heimbach 1.5 call cal_TimePassed( xx_sfluxstartdate, modelstartdate,
226     & difftime, mythid )
227     call cal_ToSeconds ( difftime, diffsecs, mythid )
228 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
229 heimbach 1.5 & xx_sfluxperiod) + 1
230 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
231 heimbach 1.5 & xx_sfluxperiod) + 2
232 heimbach 1.11 # else
233 heimbach 1.5 startrec = 1
234     endrec = 1
235 heimbach 1.11 # endif
236     diffrec = endrec - startrec + 1
237     call ctrl_init_ctrlvar (
238     & xx_sflux_file, 4, 104, diffrec, startrec, endrec,
239     & snx, sny, 1, 'c', 'xy', mythid )
240    
241 heimbach 1.5 #elif (defined (ALLOW_AQH_CONTROL))
242     c-- Atmos. humidity
243    
244 heimbach 1.11 # ifdef ALLOW_CAL
245     call cal_FullDate( xx_aqhstartdate1, xx_aqhstartdate2,
246     & xx_aqhstartdate , mythid )
247 heimbach 1.5 call cal_TimePassed( xx_aqhstartdate, modelstartdate,
248     & difftime, mythid )
249     call cal_ToSeconds ( difftime, diffsecs, mythid )
250 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
251 heimbach 1.5 & xx_aqhperiod) + 1
252 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
253 heimbach 1.5 & xx_aqhperiod) + 2
254 heimbach 1.11 # else
255 heimbach 1.5 startrec = 1
256     endrec = 1
257 heimbach 1.11 # endif
258     diffrec = endrec - startrec + 1
259     call ctrl_init_ctrlvar (
260     & xx_aqh_file, 8, 108, diffrec, startrec, endrec,
261     & snx, sny, 1, 'c', 'xy', mythid )
262 heimbach 1.5
263     #elif (defined (ALLOW_SFLUX0_CONTROL))
264     c-- initial forcing only
265 heimbach 1.11 call ctrl_init_ctrlvar (
266     & xx_sflux_file, 4, 104, 1, 1, 1,
267     & snx, sny, 1, 'c', 'xy', mythid )
268 heimbach 1.1
269 heimbach 1.5 #endif /* ALLOW_SFLUX_CONTROL */
270    
271 heimbach 1.16 c----------------------------------------------------------------------
272 heimbach 1.6 c--
273 heimbach 1.5 #if (defined (ALLOW_USTRESS_CONTROL))
274     c-- Zonal wind stress.
275    
276 heimbach 1.11 # ifdef ALLOW_CAL
277     call cal_FullDate( xx_tauustartdate1, xx_tauustartdate2,
278     & xx_tauustartdate, mythid )
279 heimbach 1.5 call cal_TimePassed( xx_tauustartdate, modelstartdate,
280     & difftime, mythid )
281     call cal_ToSeconds ( difftime, diffsecs, mythid )
282 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
283 heimbach 1.5 & xx_tauuperiod) + 1
284 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
285 heimbach 1.5 & xx_tauuperiod) + 2
286 heimbach 1.11 # else
287 heimbach 1.5 startrec = 1
288     endrec = 1
289 heimbach 1.11 # endif
290     diffrec = endrec - startrec + 1
291     call ctrl_init_ctrlvar (
292     & xx_tauu_file, 5, 105, diffrec, startrec, endrec,
293     & snx, sny, 1, 'w', 'xy', mythid )
294 heimbach 1.5
295     #elif (defined (ALLOW_UWIND_CONTROL))
296     c-- Zonal wind speed.
297    
298 heimbach 1.11 # ifdef ALLOW_CAL
299     call cal_FullDate( xx_uwindstartdate1, xx_uwindstartdate2,
300     & xx_uwindstartdate , mythid )
301 heimbach 1.5 call cal_TimePassed( xx_uwindstartdate, modelstartdate,
302     & difftime, mythid )
303     call cal_ToSeconds ( difftime, diffsecs, mythid )
304 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
305 heimbach 1.5 & xx_uwindperiod) + 1
306 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
307 heimbach 1.5 & xx_uwindperiod) + 2
308 heimbach 1.11 # else
309 heimbach 1.5 startrec = 1
310     endrec = 1
311 heimbach 1.11 # endif
312     diffrec = endrec - startrec + 1
313     call ctrl_init_ctrlvar (
314     & xx_uwind_file, 9, 109, diffrec, startrec, endrec,
315 heimbach 1.14 & snx, sny, 1, 'c', 'xy', mythid )
316 heimbach 1.5
317     #elif (defined (ALLOW_TAUU0_CONTROL))
318     c-- initial forcing only
319 heimbach 1.11 call ctrl_init_ctrlvar (
320     & xx_tauu_file, 5, 105, 1, 1, 1,
321     & snx, sny, 1, 'w', 'xy', mythid )
322 heimbach 1.1
323 heimbach 1.5 #endif /* ALLOW_USTRESS_CONTROL */
324    
325 heimbach 1.16 c----------------------------------------------------------------------
326 heimbach 1.6 c--
327 heimbach 1.5 #if (defined (ALLOW_VSTRESS_CONTROL))
328     c-- Meridional wind stress.
329    
330 heimbach 1.11 # ifdef ALLOW_CAL
331     call cal_FullDate( xx_tauvstartdate1, xx_tauvstartdate2,
332     & xx_tauvstartdate, mythid )
333 heimbach 1.5 call cal_TimePassed( xx_tauvstartdate, modelstartdate,
334     & difftime, mythid )
335     call cal_ToSeconds ( difftime, diffsecs, mythid )
336 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
337 heimbach 1.5 & xx_tauvperiod) + 1
338 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
339 heimbach 1.5 & xx_tauvperiod) + 2
340 heimbach 1.11 # else
341 heimbach 1.5 startrec = 1
342     endrec = 1
343 heimbach 1.11 # endif
344     diffrec = endrec - startrec + 1
345     call ctrl_init_ctrlvar (
346     & xx_tauv_file, 6, 106, diffrec, startrec, endrec,
347     & snx, sny, 1, 's', 'xy', mythid )
348 heimbach 1.5
349     #elif (defined (ALLOW_VWIND_CONTROL))
350     c-- Meridional wind speed.
351    
352 heimbach 1.11 # ifdef ALLOW_CAL
353     call cal_FullDate( xx_vwindstartdate1, xx_vwindstartdate2,
354     & xx_vwindstartdate , mythid )
355 heimbach 1.5 call cal_TimePassed( xx_vwindstartdate, modelstartdate,
356     & difftime, mythid )
357     call cal_ToSeconds ( difftime, diffsecs, mythid )
358 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
359 heimbach 1.5 & xx_vwindperiod) + 1
360 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
361 heimbach 1.5 & xx_vwindperiod) + 2
362 heimbach 1.11 # else
363 heimbach 1.5 startrec = 1
364     endrec = 1
365 heimbach 1.11 # endif
366     diffrec = endrec - startrec + 1
367     call ctrl_init_ctrlvar (
368     & xx_vwind_file, 10, 110, diffrec, startrec, endrec,
369 heimbach 1.14 & snx, sny, 1, 'c', 'xy', mythid )
370 heimbach 1.5
371     #elif (defined (ALLOW_TAUV0_CONTROL))
372     c-- initial forcing only
373 heimbach 1.11 diffrec = endrec - startrec + 1
374     call ctrl_init_ctrlvar (
375     & xx_tauv_file, 6, 106, 1, 1, 1,
376     & snx, sny, 1, 's', 'xy', mythid )
377 heimbach 1.1
378 heimbach 1.5 #endif /* ALLOW_VSTRESS_CONTROL */
379    
380 heimbach 1.11 c-- ===========================
381     c-- Open boundary contributions.
382     c-- ===========================
383    
384 heimbach 1.16 c----------------------------------------------------------------------
385 heimbach 1.6 c--
386 heimbach 1.5 #ifdef ALLOW_OBCSN_CONTROL
387     c-- Northern obc.
388    
389 heimbach 1.11 # ifdef ALLOW_CAL
390     call cal_FullDate( xx_obcsnstartdate1, xx_obcsnstartdate2,
391     & xx_obcsnstartdate, mythid )
392 heimbach 1.5 call cal_TimePassed( xx_obcsnstartdate, modelstartdate,
393     & difftime, mythid )
394     call cal_ToSeconds ( difftime, diffsecs, mythid )
395 heimbach 1.9 startrec = int((modelstart - diffsecs)/xx_obcsnperiod) + 1
396 heimbach 1.11 startrec = (startrec - 1)*nobcs + 1
397 heimbach 1.9 endrec = int((modelend - diffsecs)/xx_obcsnperiod) + 2
398 heimbach 1.11 endrec = (endrec - startrec + 1)*nobcs
399     # else
400 heimbach 1.5 startrec = 1
401     endrec = 1
402 heimbach 1.11 # endif
403 mlosch 1.13 diffrec = endrec
404 heimbach 1.11 call ctrl_init_ctrlvar (
405     & xx_obcsn_file, 11, 111, diffrec, startrec, endrec,
406     & snx, sny, nr, 'm', 'xz', mythid )
407 heimbach 1.5
408     #endif /* ALLOW_OBCSN_CONTROL */
409    
410 heimbach 1.16 c----------------------------------------------------------------------
411 heimbach 1.11 c--
412 heimbach 1.5 #ifdef ALLOW_OBCSS_CONTROL
413     c-- Southern obc.
414    
415 heimbach 1.11 # ifdef ALLOW_CAL
416     call cal_FullDate( xx_obcssstartdate1, xx_obcssstartdate2,
417     & xx_obcssstartdate, mythid )
418 heimbach 1.5 call cal_TimePassed( xx_obcssstartdate, modelstartdate,
419     & difftime, mythid )
420     call cal_ToSeconds ( difftime, diffsecs, mythid )
421 heimbach 1.9 startrec = int((modelstart - diffsecs)/xx_obcssperiod) + 1
422 heimbach 1.11 startrec = (startrec - 1)*nobcs + 1
423 heimbach 1.9 endrec = int((modelend - diffsecs)/xx_obcssperiod) + 2
424 heimbach 1.11 endrec = (endrec - startrec + 1)*nobcs
425     # else
426 heimbach 1.5 startrec = 1
427     endrec = 1
428 heimbach 1.11 # endif
429 mlosch 1.13 diffrec = endrec
430 heimbach 1.11 call ctrl_init_ctrlvar (
431     & xx_obcss_file, 12, 112, diffrec, startrec, endrec,
432     & snx, sny, nr, 'm', 'xz', mythid )
433 heimbach 1.5
434     #endif /* ALLOW_OBCSS_CONTROL */
435    
436 heimbach 1.16 c----------------------------------------------------------------------
437 heimbach 1.6 c--
438 heimbach 1.5 #ifdef ALLOW_OBCSW_CONTROL
439     c-- Western obc.
440    
441 heimbach 1.11 # ifdef ALLOW_CAL
442     call cal_FullDate( xx_obcswstartdate1, xx_obcswstartdate2,
443     & xx_obcswstartdate, mythid )
444 heimbach 1.5 call cal_TimePassed( xx_obcswstartdate, modelstartdate,
445     & difftime, mythid )
446     call cal_ToSeconds ( difftime, diffsecs, mythid )
447 heimbach 1.9 startrec = int((modelstart - diffsecs)/xx_obcswperiod) + 1
448 heimbach 1.11 startrec = (startrec - 1)*nobcs + 1
449 heimbach 1.9 endrec = int((modelend - diffsecs)/xx_obcswperiod) + 2
450 heimbach 1.11 endrec = (endrec - startrec + 1)*nobcs
451     # else
452 heimbach 1.5 startrec = 1
453     endrec = 1
454 heimbach 1.11 # endif
455 mlosch 1.13 diffrec = endrec
456 heimbach 1.11 call ctrl_init_ctrlvar (
457     & xx_obcsw_file, 13, 113, diffrec, startrec, endrec,
458     & snx, sny, nr, 'm', 'yz', mythid )
459 heimbach 1.5
460 heimbach 1.11 #endif /* ALLOW_OBCSW_CONTROL */
461 heimbach 1.5
462 heimbach 1.16 c----------------------------------------------------------------------
463 heimbach 1.6 c--
464 heimbach 1.5 #ifdef ALLOW_OBCSE_CONTROL
465     c-- Eastern obc.
466    
467 heimbach 1.11 # ifdef ALLOW_CAL
468     call cal_FullDate( xx_obcsestartdate1, xx_obcsestartdate2,
469     & xx_obcsestartdate, mythid )
470 heimbach 1.5 call cal_TimePassed( xx_obcsestartdate, modelstartdate,
471     & difftime, mythid )
472     call cal_ToSeconds ( difftime, diffsecs, mythid )
473 heimbach 1.9 startrec = int((modelstart - diffsecs)/xx_obcseperiod) + 1
474 heimbach 1.11 startrec = (startrec - 1)*nobcs + 1
475 heimbach 1.9 endrec = int((modelend - diffsecs)/xx_obcseperiod) + 2
476 heimbach 1.11 endrec = (endrec - startrec + 1)*nobcs
477     # else
478 heimbach 1.5 startrec = 1
479     endrec = 1
480 heimbach 1.11 # endif
481 mlosch 1.13 diffrec = endrec
482 heimbach 1.11 call ctrl_init_ctrlvar (
483     & xx_obcse_file, 14, 114, diffrec, startrec, endrec,
484     & snx, sny, nr, 'm', 'yz', mythid )
485 heimbach 1.5
486     #endif /* ALLOW_OBCSE_CONTROL */
487 heimbach 1.2
488 heimbach 1.16 c----------------------------------------------------------------------
489 heimbach 1.6 c--
490 heimbach 1.3 #ifdef ALLOW_DIFFKR_CONTROL
491 heimbach 1.11 call ctrl_init_ctrlvar (
492     & xx_diffkr_file, 15, 115, 1, 1, 1,
493     & snx, sny, nr, 'c', '3d', mythid )
494 heimbach 1.3 #endif /* ALLOW_DIFFKR_CONTROL */
495    
496 heimbach 1.16 c----------------------------------------------------------------------
497 heimbach 1.6 c--
498 heimbach 1.3 #ifdef ALLOW_KAPGM_CONTROL
499 heimbach 1.11 call ctrl_init_ctrlvar (
500     & xx_kapgm_file, 16, 116, 1, 1, 1,
501     & snx, sny, nr, 'c', '3d', mythid )
502 heimbach 1.3 #endif /* ALLOW_KAPGM_CONTROL */
503    
504 heimbach 1.16 c----------------------------------------------------------------------
505 heimbach 1.6 c--
506 heimbach 1.17 #ifdef ALLOW_PRECIP_CONTROL
507     c-- Atmos. precipitation
508    
509     # ifdef ALLOW_CAL
510     call cal_FullDate( xx_precipstartdate1, xx_precipstartdate2,
511     & xx_precipstartdate , mythid )
512     call cal_TimePassed( xx_precipstartdate, modelstartdate,
513     & difftime, mythid )
514     call cal_ToSeconds ( difftime, diffsecs, mythid )
515     startrec = int((modelstart + startTime - diffsecs)/
516     & xx_precipperiod) + 1
517     endrec = int((modelend + startTime - diffsecs + modelstep/2)/
518     & xx_precipperiod) + 2
519     # else
520     startrec = 1
521     endrec = 1
522     # endif
523     diffrec = endrec - startrec + 1
524 heimbach 1.11 call ctrl_init_ctrlvar (
525 heimbach 1.17 & xx_precip_file, 17, 117, diffrec, startrec, endrec,
526     & snx, sny, 1, 'c', 'xy', mythid )
527    
528     #endif /* ALLOW_PRECIP_CONTROL */
529 heimbach 1.2
530 heimbach 1.16 c----------------------------------------------------------------------
531 heimbach 1.6 c--
532     #ifdef ALLOW_SST0_CONTROL
533 heimbach 1.11 call ctrl_init_ctrlvar (
534     & xx_sst_file, 18, 118, 1, 1, 1,
535     & snx, sny, 1, 'c', 'xy', mythid )
536 heimbach 1.6 #endif /* ALLOW_SST0_CONTROL */
537    
538 heimbach 1.16 c----------------------------------------------------------------------
539 heimbach 1.6 c--
540     #ifdef ALLOW_SSS0_CONTROL
541 heimbach 1.11 call ctrl_init_ctrlvar (
542     & xx_sss_file, 19, 119, 1, 1, 1,
543     & snx, sny, 1, 'c', 'xy', mythid )
544 heimbach 1.6 #endif /* ALLOW_SSS0_CONTROL */
545    
546 heimbach 1.16 c----------------------------------------------------------------------
547 heimbach 1.6 c--
548     #ifdef ALLOW_HFACC_CONTROL
549 heimbach 1.11 # ifdef ALLOW_HFACC3D_CONTROL
550     call ctrl_init_ctrlvar (
551     & xx_hfacc_file, 20, 120, 1, 1, 1,
552     & snx, sny, nr, 'c', '3d', mythid )
553     # else
554     call ctrl_init_ctrlvar (
555     & xx_hfacc_file, 20, 120, 1, 1, 1,
556     & snx, sny, 1, 'c', 'xy', mythid )
557     # endif /*ALLOW_HFACC3D_CONTROL*/
558 heimbach 1.6 #endif /* ALLOW_HFACC_CONTROL */
559    
560 heimbach 1.16 c----------------------------------------------------------------------
561 heimbach 1.6 c--
562 heimbach 1.5 #ifdef ALLOW_EFLUXY0_CONTROL
563 heimbach 1.11 call ctrl_init_ctrlvar (
564     & xx_efluxy_file, 21, 121, 1, 1, 1,
565     & snx, sny, nr, 's', '3d', mythid )
566 heimbach 1.5 #endif /* ALLOW_EFLUXY0_CONTROL */
567    
568 heimbach 1.16 c----------------------------------------------------------------------
569 heimbach 1.6 c--
570 heimbach 1.5 #ifdef ALLOW_EFLUXP0_CONTROL
571 heimbach 1.11 call ctrl_init_ctrlvar (
572     & xx_efluxp_file, 22, 122, 1, 1, 1,
573     & snx, sny, nr, 'v', '3d', mythid )
574 heimbach 1.5 #endif /* ALLOW_EFLUXP0_CONTROL */
575 heimbach 1.6
576 heimbach 1.16 c----------------------------------------------------------------------
577 heimbach 1.6 c--
578     #ifdef ALLOW_BOTTOMDRAG_CONTROL
579 heimbach 1.11 call ctrl_init_ctrlvar (
580     & xx_bottomdrag_file, 23, 123, 1, 1, 1,
581     & snx, sny, 1, 'c', 'xy', mythid )
582 heimbach 1.6 #endif /* ALLOW_BOTTOMDRAG_CONTROL */
583    
584 heimbach 1.16 c----------------------------------------------------------------------
585 heimbach 1.15 c--
586     #ifdef ALLOW_EDTAUX_CONTROL
587     call ctrl_init_ctrlvar (
588     & xx_edtaux_file, 25, 125, 1, 1, 1,
589     & snx, sny, nr, 'w', '3d', mythid )
590     #endif /* ALLOW_EDTAUX_CONTROL */
591    
592 heimbach 1.16 c----------------------------------------------------------------------
593 heimbach 1.15 c--
594     #ifdef ALLOW_EDTAUY_CONTROL
595     call ctrl_init_ctrlvar (
596     & xx_edtauy_file, 26, 126, 1, 1, 1,
597     & snx, sny, nr, 's', '3d', mythid )
598     #endif /* ALLOW_EDTAUY_CONTROL */
599    
600 heimbach 1.16 c----------------------------------------------------------------------
601     c--
602     #ifdef ALLOW_UVEL0_CONTROL
603     call ctrl_init_ctrlvar (
604     & xx_uvel_file, 27, 127, 1, 1, 1,
605     & snx, sny, nr, 'w', '3d', mythid )
606     #endif /* ALLOW_UVEL0_CONTROL */
607    
608     c----------------------------------------------------------------------
609     c--
610     #ifdef ALLOW_VVEL0_CONTROL
611     call ctrl_init_ctrlvar (
612     & xx_vvel_file, 28, 128, 1, 1, 1,
613     & snx, sny, nr, 's', '3d', mythid )
614     #endif /* ALLOW_VVEL0_CONTROL */
615    
616     c----------------------------------------------------------------------
617     c--
618     #ifdef ALLOW_ETAN0_CONTROL
619     call ctrl_init_ctrlvar (
620     & xx_etan_file, 29, 129, 1, 1, 1,
621     & snx, sny, 1, 'c', 'xy', mythid )
622     #endif /* ALLOW_VVEL0_CONTROL */
623    
624     c----------------------------------------------------------------------
625     c--
626     #ifdef ALLOW_RELAXSST_CONTROL
627     call ctrl_init_ctrlvar (
628     & xx_relaxsst_file, 30, 130, 1, 1, 1,
629     & snx, sny, 1, 'c', 'xy', mythid )
630     #endif /* ALLOW_RELAXSST_CONTROL */
631    
632     c----------------------------------------------------------------------
633     c--
634     #ifdef ALLOW_RELAXSSS_CONTROL
635     call ctrl_init_ctrlvar (
636     & xx_relaxsss_file, 31, 131, 1, 1, 1,
637     & snx, sny, 1, 'c', 'xy', mythid )
638     #endif /* ALLOW_RELAXSSS_CONTROL */
639    
640     c----------------------------------------------------------------------
641 heimbach 1.17 c--
642     #ifdef ALLOW_TR10_CONTROL
643     call ctrl_init_ctrlvar (
644     & xx_tr1_file, 32, 132, 1, 1, 1,
645     & snx, sny, nr, 'c', '3d', mythid )
646     #endif /* ALLOW_TR10_CONTROL */
647    
648     c----------------------------------------------------------------------
649 heimbach 1.16 c----------------------------------------------------------------------
650     c----------------------------------------------------------------------
651 heimbach 1.1
652 heimbach 1.11 call ctrl_init_wet( mythid )
653    
654 heimbach 1.1 return
655     end
656    

  ViewVC Help
Powered by ViewVC 1.1.22