/[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.20 - (hide annotations) (download)
Wed Aug 31 00:03:30 2005 UTC (18 years, 9 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57s_post, checkpoint57y_post, checkpoint57r_post, checkpoint57z_post, checkpoint57t_post, checkpoint57v_post, checkpoint57y_pre, checkpint57u_post, checkpoint57w_post, checkpoint57x_post
Changes since 1.19: +50 -4 lines
Adding time-dependent SST, SSS control.

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

  ViewVC Help
Powered by ViewVC 1.1.22