/[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.24 - (hide annotations) (download)
Fri Oct 27 05:16:54 2006 UTC (17 years, 6 months ago) by heimbach
Branch: MAIN
Changes since 1.23: +180 -1 lines
Adding new control variables:
lwflux, lwdown, evap, snowprecip, apressure, runoff.

1 edhill 1.10 C
2 heimbach 1.24 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_init.F,v 1.23 2006/06/07 01:55:13 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 heimbach 1.21 integer bi,bj
86     integer i,j,k
87     integer itlo,ithi
88     integer jtlo,jthi
89     integer jmin,jmax
90     integer imin,imax
91    
92 heimbach 1.1 integer ntmp
93 heimbach 1.11 integer ivar
94 heimbach 1.5 integer iobcs
95 heimbach 1.1 integer il
96     integer errio
97     integer startrec
98     integer endrec
99 heimbach 1.11 integer diffrec
100 heimbach 1.5 integer difftime(4)
101     _RL diffsecs
102 heimbach 1.1
103     character*(max_len_prec) record
104     character*(max_len_mbuf) msgbuf
105 heimbach 1.11 character*2 whichxyz
106 heimbach 1.5
107 heimbach 1.1 c == external ==
108    
109     integer ilnblnk
110     external ilnblnk
111    
112     c == end of interface ==
113    
114 heimbach 1.21 jtlo = mybylo(mythid)
115     jthi = mybyhi(mythid)
116     itlo = mybxlo(mythid)
117     ithi = mybxhi(mythid)
118     jmin = 1-oly
119     jmax = sny+oly
120     imin = 1-olx
121     imax = snx+olx
122    
123 heimbach 1.1 c-- Set default values.
124 heimbach 1.11 do ivar = 1,maxcvars
125     ncvarindex(ivar) = -1
126     ncvarrecs(ivar) = 0
127     ncvarxmax(ivar) = 0
128     ncvarymax(ivar) = 0
129     ncvarnrmax(ivar) = 0
130     ncvargrd(ivar) = '?'
131     enddo
132 heimbach 1.1
133     _BARRIER
134    
135     c-- =====================
136     c-- Initial state fields.
137     c-- =====================
138    
139 heimbach 1.5 cph(
140     cph index 7-10 reserved for atmos. state,
141     cph index 11-14 reserved for open boundaries,
142     cph index 15-16 reserved for mixing coeff.
143 heimbach 1.18 cph index 17 reserved for passive tracer TR1
144 heimbach 1.5 cph index 18,19 reserved for sst, sss
145     cph index 20 for hFacC
146 heimbach 1.6 cph index 21-22 for efluxy, efluxp
147 heimbach 1.17 cph index 23 for bottom drag
148     cph index 24
149     cph index 25-26 for edtaux, edtauy
150     cph index 27-29 for uvel0, vvel0, etan0
151     cph index 30-31 for relax. SST, SSS
152 heimbach 1.18 cph index 32 reserved for precip (atmos. state)
153     cph index 33 reserved for swflux (atmos. state)
154 heimbach 1.19 cph index 34 reserved for swdown (atmos. state)
155 heimbach 1.24 cph 35 lwflux
156     cph 36 lwdown
157     cph 37 evap
158     cph 38 snowprecip
159     cph 39 apressure
160     cph 40 runoff
161 heimbach 1.5 cph)
162    
163 heimbach 1.16 c----------------------------------------------------------------------
164 heimbach 1.6 c--
165 heimbach 1.1 #ifdef ALLOW_THETA0_CONTROL
166 heimbach 1.5 c-- Initial state temperature contribution.
167 heimbach 1.11 call ctrl_init_ctrlvar (
168     & xx_theta_file, 1, 101, 1, 1, 1,
169     & snx, sny, nr, 'c', '3d', mythid )
170 heimbach 1.1 #endif /* ALLOW_THETA0_CONTROL */
171    
172 heimbach 1.16 c----------------------------------------------------------------------
173 heimbach 1.6 c--
174 heimbach 1.1 #ifdef ALLOW_SALT0_CONTROL
175 heimbach 1.5 c-- Initial state salinity contribution.
176 heimbach 1.11 call ctrl_init_ctrlvar (
177     & xx_salt_file, 2, 102, 1, 1, 1,
178     & snx, sny, nr, 'c', '3d', mythid )
179 heimbach 1.1 #endif /* ALLOW_SALT0_CONTROL */
180    
181 heimbach 1.5 c-- ===========================
182     c-- Surface flux contributions.
183     c-- ===========================
184    
185 heimbach 1.16 c----------------------------------------------------------------------
186 heimbach 1.6 c--
187 heimbach 1.5 #if (defined (ALLOW_HFLUX_CONTROL))
188     c-- Heat flux.
189    
190 heimbach 1.11 # ifdef ALLOW_CAL
191     call cal_FullDate( xx_hfluxstartdate1, xx_hfluxstartdate2,
192     & xx_hfluxstartdate , mythid )
193 heimbach 1.5 call cal_TimePassed( xx_hfluxstartdate, modelstartdate,
194     & difftime, mythid )
195     call cal_ToSeconds ( difftime, diffsecs, mythid )
196 gforget 1.22 if ( xx_hfluxperiod .EQ. 0 ) then
197     startrec=1
198     endrec=12
199     else
200 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
201 heimbach 1.5 & xx_hfluxperiod) + 1
202 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
203 heimbach 1.5 & xx_hfluxperiod) + 2
204 gforget 1.22 endif
205 heimbach 1.11 # else
206 heimbach 1.5 startrec = 1
207     endrec = 1
208 heimbach 1.11 # endif
209     diffrec = endrec - startrec + 1
210     call ctrl_init_ctrlvar (
211     & xx_hflux_file, 3, 103, diffrec, startrec, endrec,
212     & snx, sny, 1, 'c', 'xy', mythid )
213 heimbach 1.5
214     #elif (defined (ALLOW_ATEMP_CONTROL))
215     c-- Atmos. temperature
216    
217 heimbach 1.11 # ifdef ALLOW_CAL
218     call cal_FullDate( xx_atempstartdate1, xx_atempstartdate2,
219     & xx_atempstartdate , mythid )
220 heimbach 1.5 call cal_TimePassed( xx_atempstartdate, modelstartdate,
221     & difftime, mythid )
222     call cal_ToSeconds ( difftime, diffsecs, mythid )
223 gforget 1.22 if ( xx_atempperiod .EQ. 0 ) then
224     startrec=1
225     endrec=12
226     else
227 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
228 heimbach 1.5 & xx_atempperiod) + 1
229 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
230 heimbach 1.5 & xx_atempperiod) + 2
231 gforget 1.22 endif
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_atemp_file, 7, 107, diffrec, startrec, endrec,
239     & snx, sny, 1, 'c', 'xy', mythid )
240 heimbach 1.5
241     #elif (defined (ALLOW_HFLUX0_CONTROL))
242     c-- initial forcing only
243 heimbach 1.11 call ctrl_init_ctrlvar (
244     & xx_hflux_file, 3, 103, 1, 1, 1,
245     & snx, sny, 1, 'c', 'xy', mythid )
246 heimbach 1.1
247 heimbach 1.5 #endif /* ALLOW_HFLUX_CONTROL */
248    
249 heimbach 1.16 c----------------------------------------------------------------------
250 heimbach 1.6 c--
251 heimbach 1.5 #if (defined (ALLOW_SFLUX_CONTROL))
252     c-- Salt flux.
253    
254 heimbach 1.11 # ifdef ALLOW_CAL
255     call cal_FullDate( xx_sfluxstartdate1, xx_sfluxstartdate2,
256     & xx_sfluxstartdate , mythid )
257 heimbach 1.5 call cal_TimePassed( xx_sfluxstartdate, modelstartdate,
258     & difftime, mythid )
259     call cal_ToSeconds ( difftime, diffsecs, mythid )
260 gforget 1.22 if ( xx_sfluxperiod .EQ. 0 ) then
261     startrec=1
262     endrec=12
263     else
264 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
265 heimbach 1.5 & xx_sfluxperiod) + 1
266 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
267 heimbach 1.5 & xx_sfluxperiod) + 2
268 gforget 1.22 endif
269 heimbach 1.11 # else
270 heimbach 1.5 startrec = 1
271     endrec = 1
272 heimbach 1.11 # endif
273     diffrec = endrec - startrec + 1
274     call ctrl_init_ctrlvar (
275     & xx_sflux_file, 4, 104, diffrec, startrec, endrec,
276     & snx, sny, 1, 'c', 'xy', mythid )
277    
278 heimbach 1.5 #elif (defined (ALLOW_AQH_CONTROL))
279     c-- Atmos. humidity
280    
281 heimbach 1.11 # ifdef ALLOW_CAL
282     call cal_FullDate( xx_aqhstartdate1, xx_aqhstartdate2,
283     & xx_aqhstartdate , mythid )
284 heimbach 1.5 call cal_TimePassed( xx_aqhstartdate, modelstartdate,
285     & difftime, mythid )
286     call cal_ToSeconds ( difftime, diffsecs, mythid )
287 gforget 1.22 if ( xx_aqhperiod .EQ. 0 ) then
288     startrec=1
289     endrec=12
290     else
291 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
292 heimbach 1.5 & xx_aqhperiod) + 1
293 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
294 heimbach 1.5 & xx_aqhperiod) + 2
295 gforget 1.22 endif
296 heimbach 1.11 # else
297 heimbach 1.5 startrec = 1
298     endrec = 1
299 heimbach 1.11 # endif
300     diffrec = endrec - startrec + 1
301     call ctrl_init_ctrlvar (
302     & xx_aqh_file, 8, 108, diffrec, startrec, endrec,
303     & snx, sny, 1, 'c', 'xy', mythid )
304 heimbach 1.5
305     #elif (defined (ALLOW_SFLUX0_CONTROL))
306     c-- initial forcing only
307 heimbach 1.11 call ctrl_init_ctrlvar (
308     & xx_sflux_file, 4, 104, 1, 1, 1,
309     & snx, sny, 1, 'c', 'xy', mythid )
310 heimbach 1.1
311 heimbach 1.5 #endif /* ALLOW_SFLUX_CONTROL */
312    
313 heimbach 1.16 c----------------------------------------------------------------------
314 heimbach 1.6 c--
315 heimbach 1.5 #if (defined (ALLOW_USTRESS_CONTROL))
316     c-- Zonal wind stress.
317    
318 heimbach 1.11 # ifdef ALLOW_CAL
319     call cal_FullDate( xx_tauustartdate1, xx_tauustartdate2,
320     & xx_tauustartdate, mythid )
321 heimbach 1.5 call cal_TimePassed( xx_tauustartdate, modelstartdate,
322     & difftime, mythid )
323     call cal_ToSeconds ( difftime, diffsecs, mythid )
324 gforget 1.22 if ( xx_tauuperiod .EQ. 0 ) then
325     startrec=1
326     endrec=12
327     else
328 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
329 heimbach 1.5 & xx_tauuperiod) + 1
330 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
331 heimbach 1.5 & xx_tauuperiod) + 2
332 gforget 1.22 endif
333 heimbach 1.11 # else
334 heimbach 1.5 startrec = 1
335     endrec = 1
336 heimbach 1.11 # endif
337     diffrec = endrec - startrec + 1
338     call ctrl_init_ctrlvar (
339     & xx_tauu_file, 5, 105, diffrec, startrec, endrec,
340     & snx, sny, 1, 'w', 'xy', mythid )
341 heimbach 1.5
342     #elif (defined (ALLOW_UWIND_CONTROL))
343     c-- Zonal wind speed.
344    
345 heimbach 1.11 # ifdef ALLOW_CAL
346     call cal_FullDate( xx_uwindstartdate1, xx_uwindstartdate2,
347     & xx_uwindstartdate , mythid )
348 heimbach 1.5 call cal_TimePassed( xx_uwindstartdate, modelstartdate,
349     & difftime, mythid )
350     call cal_ToSeconds ( difftime, diffsecs, mythid )
351 gforget 1.22 if ( xx_uwindperiod .EQ. 0 ) then
352     startrec=1
353     endrec=12
354     else
355 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
356 heimbach 1.5 & xx_uwindperiod) + 1
357 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
358 heimbach 1.5 & xx_uwindperiod) + 2
359 gforget 1.22 endif
360 heimbach 1.11 # else
361 heimbach 1.5 startrec = 1
362     endrec = 1
363 heimbach 1.11 # endif
364     diffrec = endrec - startrec + 1
365     call ctrl_init_ctrlvar (
366     & xx_uwind_file, 9, 109, diffrec, startrec, endrec,
367 heimbach 1.14 & snx, sny, 1, 'c', 'xy', mythid )
368 heimbach 1.5
369     #elif (defined (ALLOW_TAUU0_CONTROL))
370     c-- initial forcing only
371 heimbach 1.11 call ctrl_init_ctrlvar (
372     & xx_tauu_file, 5, 105, 1, 1, 1,
373     & snx, sny, 1, 'w', 'xy', mythid )
374 heimbach 1.1
375 heimbach 1.5 #endif /* ALLOW_USTRESS_CONTROL */
376    
377 heimbach 1.16 c----------------------------------------------------------------------
378 heimbach 1.6 c--
379 heimbach 1.5 #if (defined (ALLOW_VSTRESS_CONTROL))
380     c-- Meridional wind stress.
381    
382 heimbach 1.11 # ifdef ALLOW_CAL
383     call cal_FullDate( xx_tauvstartdate1, xx_tauvstartdate2,
384     & xx_tauvstartdate, mythid )
385 heimbach 1.5 call cal_TimePassed( xx_tauvstartdate, modelstartdate,
386     & difftime, mythid )
387     call cal_ToSeconds ( difftime, diffsecs, mythid )
388 gforget 1.22 if ( xx_tauvperiod .EQ. 0 ) then
389     startrec=1
390     endrec=12
391     else
392 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
393 heimbach 1.5 & xx_tauvperiod) + 1
394 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
395 heimbach 1.5 & xx_tauvperiod) + 2
396 gforget 1.22 endif
397 heimbach 1.11 # else
398 heimbach 1.5 startrec = 1
399     endrec = 1
400 heimbach 1.11 # endif
401     diffrec = endrec - startrec + 1
402     call ctrl_init_ctrlvar (
403     & xx_tauv_file, 6, 106, diffrec, startrec, endrec,
404     & snx, sny, 1, 's', 'xy', mythid )
405 heimbach 1.5
406     #elif (defined (ALLOW_VWIND_CONTROL))
407     c-- Meridional wind speed.
408    
409 heimbach 1.11 # ifdef ALLOW_CAL
410     call cal_FullDate( xx_vwindstartdate1, xx_vwindstartdate2,
411     & xx_vwindstartdate , mythid )
412 heimbach 1.5 call cal_TimePassed( xx_vwindstartdate, modelstartdate,
413     & difftime, mythid )
414     call cal_ToSeconds ( difftime, diffsecs, mythid )
415 gforget 1.22 if ( xx_vwindperiod .EQ. 0 ) then
416     startrec=1
417     endrec=12
418     else
419 heimbach 1.9 startrec = int((modelstart + startTime - diffsecs)/
420 heimbach 1.5 & xx_vwindperiod) + 1
421 heimbach 1.9 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
422 heimbach 1.5 & xx_vwindperiod) + 2
423 gforget 1.22 endif
424 heimbach 1.11 # else
425 heimbach 1.5 startrec = 1
426     endrec = 1
427 heimbach 1.11 # endif
428     diffrec = endrec - startrec + 1
429     call ctrl_init_ctrlvar (
430     & xx_vwind_file, 10, 110, diffrec, startrec, endrec,
431 heimbach 1.14 & snx, sny, 1, 'c', 'xy', mythid )
432 heimbach 1.5
433     #elif (defined (ALLOW_TAUV0_CONTROL))
434     c-- initial forcing only
435 heimbach 1.11 diffrec = endrec - startrec + 1
436     call ctrl_init_ctrlvar (
437     & xx_tauv_file, 6, 106, 1, 1, 1,
438     & snx, sny, 1, 's', 'xy', mythid )
439 heimbach 1.1
440 heimbach 1.5 #endif /* ALLOW_VSTRESS_CONTROL */
441    
442 heimbach 1.11 c-- ===========================
443     c-- Open boundary contributions.
444     c-- ===========================
445    
446 heimbach 1.16 c----------------------------------------------------------------------
447 heimbach 1.6 c--
448 heimbach 1.5 #ifdef ALLOW_OBCSN_CONTROL
449     c-- Northern obc.
450    
451 heimbach 1.11 # ifdef ALLOW_CAL
452     call cal_FullDate( xx_obcsnstartdate1, xx_obcsnstartdate2,
453     & xx_obcsnstartdate, mythid )
454 heimbach 1.5 call cal_TimePassed( xx_obcsnstartdate, modelstartdate,
455     & difftime, mythid )
456     call cal_ToSeconds ( difftime, diffsecs, mythid )
457 heimbach 1.9 startrec = int((modelstart - diffsecs)/xx_obcsnperiod) + 1
458 heimbach 1.11 startrec = (startrec - 1)*nobcs + 1
459 heimbach 1.9 endrec = int((modelend - diffsecs)/xx_obcsnperiod) + 2
460 heimbach 1.11 endrec = (endrec - startrec + 1)*nobcs
461     # else
462 heimbach 1.5 startrec = 1
463     endrec = 1
464 heimbach 1.11 # endif
465 mlosch 1.13 diffrec = endrec
466 heimbach 1.11 call ctrl_init_ctrlvar (
467     & xx_obcsn_file, 11, 111, diffrec, startrec, endrec,
468     & snx, sny, nr, 'm', 'xz', mythid )
469 heimbach 1.5
470     #endif /* ALLOW_OBCSN_CONTROL */
471    
472 heimbach 1.16 c----------------------------------------------------------------------
473 heimbach 1.11 c--
474 heimbach 1.5 #ifdef ALLOW_OBCSS_CONTROL
475     c-- Southern obc.
476    
477 heimbach 1.11 # ifdef ALLOW_CAL
478     call cal_FullDate( xx_obcssstartdate1, xx_obcssstartdate2,
479     & xx_obcssstartdate, mythid )
480 heimbach 1.5 call cal_TimePassed( xx_obcssstartdate, modelstartdate,
481     & difftime, mythid )
482     call cal_ToSeconds ( difftime, diffsecs, mythid )
483 heimbach 1.9 startrec = int((modelstart - diffsecs)/xx_obcssperiod) + 1
484 heimbach 1.11 startrec = (startrec - 1)*nobcs + 1
485 heimbach 1.9 endrec = int((modelend - diffsecs)/xx_obcssperiod) + 2
486 heimbach 1.11 endrec = (endrec - startrec + 1)*nobcs
487     # else
488 heimbach 1.5 startrec = 1
489     endrec = 1
490 heimbach 1.11 # endif
491 mlosch 1.13 diffrec = endrec
492 heimbach 1.11 call ctrl_init_ctrlvar (
493     & xx_obcss_file, 12, 112, diffrec, startrec, endrec,
494     & snx, sny, nr, 'm', 'xz', mythid )
495 heimbach 1.5
496     #endif /* ALLOW_OBCSS_CONTROL */
497    
498 heimbach 1.16 c----------------------------------------------------------------------
499 heimbach 1.6 c--
500 heimbach 1.5 #ifdef ALLOW_OBCSW_CONTROL
501     c-- Western obc.
502    
503 heimbach 1.11 # ifdef ALLOW_CAL
504     call cal_FullDate( xx_obcswstartdate1, xx_obcswstartdate2,
505     & xx_obcswstartdate, mythid )
506 heimbach 1.5 call cal_TimePassed( xx_obcswstartdate, modelstartdate,
507     & difftime, mythid )
508     call cal_ToSeconds ( difftime, diffsecs, mythid )
509 heimbach 1.9 startrec = int((modelstart - diffsecs)/xx_obcswperiod) + 1
510 heimbach 1.11 startrec = (startrec - 1)*nobcs + 1
511 heimbach 1.9 endrec = int((modelend - diffsecs)/xx_obcswperiod) + 2
512 heimbach 1.11 endrec = (endrec - startrec + 1)*nobcs
513     # else
514 heimbach 1.5 startrec = 1
515     endrec = 1
516 heimbach 1.11 # endif
517 mlosch 1.13 diffrec = endrec
518 heimbach 1.11 call ctrl_init_ctrlvar (
519     & xx_obcsw_file, 13, 113, diffrec, startrec, endrec,
520     & snx, sny, nr, 'm', 'yz', mythid )
521 heimbach 1.5
522 heimbach 1.11 #endif /* ALLOW_OBCSW_CONTROL */
523 heimbach 1.5
524 heimbach 1.16 c----------------------------------------------------------------------
525 heimbach 1.6 c--
526 heimbach 1.5 #ifdef ALLOW_OBCSE_CONTROL
527     c-- Eastern obc.
528    
529 heimbach 1.11 # ifdef ALLOW_CAL
530     call cal_FullDate( xx_obcsestartdate1, xx_obcsestartdate2,
531     & xx_obcsestartdate, mythid )
532 heimbach 1.5 call cal_TimePassed( xx_obcsestartdate, modelstartdate,
533     & difftime, mythid )
534     call cal_ToSeconds ( difftime, diffsecs, mythid )
535 heimbach 1.9 startrec = int((modelstart - diffsecs)/xx_obcseperiod) + 1
536 heimbach 1.11 startrec = (startrec - 1)*nobcs + 1
537 heimbach 1.9 endrec = int((modelend - diffsecs)/xx_obcseperiod) + 2
538 heimbach 1.11 endrec = (endrec - startrec + 1)*nobcs
539     # else
540 heimbach 1.5 startrec = 1
541     endrec = 1
542 heimbach 1.11 # endif
543 mlosch 1.13 diffrec = endrec
544 heimbach 1.11 call ctrl_init_ctrlvar (
545     & xx_obcse_file, 14, 114, diffrec, startrec, endrec,
546     & snx, sny, nr, 'm', 'yz', mythid )
547 heimbach 1.5
548     #endif /* ALLOW_OBCSE_CONTROL */
549 heimbach 1.2
550 heimbach 1.16 c----------------------------------------------------------------------
551 heimbach 1.6 c--
552 heimbach 1.3 #ifdef ALLOW_DIFFKR_CONTROL
553 heimbach 1.11 call ctrl_init_ctrlvar (
554     & xx_diffkr_file, 15, 115, 1, 1, 1,
555     & snx, sny, nr, 'c', '3d', mythid )
556 heimbach 1.3 #endif /* ALLOW_DIFFKR_CONTROL */
557    
558 heimbach 1.16 c----------------------------------------------------------------------
559 heimbach 1.6 c--
560 heimbach 1.3 #ifdef ALLOW_KAPGM_CONTROL
561 heimbach 1.11 call ctrl_init_ctrlvar (
562     & xx_kapgm_file, 16, 116, 1, 1, 1,
563     & snx, sny, nr, 'c', '3d', mythid )
564 heimbach 1.3 #endif /* ALLOW_KAPGM_CONTROL */
565    
566 heimbach 1.16 c----------------------------------------------------------------------
567 heimbach 1.6 c--
568 heimbach 1.18 #ifdef ALLOW_TR10_CONTROL
569 heimbach 1.11 call ctrl_init_ctrlvar (
570 heimbach 1.18 & xx_tr1_file, 17, 117, 1, 1, 1,
571     & snx, sny, nr, 'c', '3d', mythid )
572     #endif /* ALLOW_TR10_CONTROL */
573 heimbach 1.2
574 heimbach 1.16 c----------------------------------------------------------------------
575 heimbach 1.6 c--
576 heimbach 1.20 #if (defined (ALLOW_SST_CONTROL))
577    
578     # ifdef ALLOW_CAL
579     call cal_FullDate( xx_sststartdate1, xx_sststartdate2,
580     & xx_sststartdate , mythid )
581     call cal_TimePassed( xx_sststartdate, modelstartdate,
582     & difftime, mythid )
583     call cal_ToSeconds ( difftime, diffsecs, mythid )
584 gforget 1.22 if ( xx_sstperiod .EQ. 0 ) then
585     startrec=1
586     endrec=12
587     else
588 heimbach 1.20 startrec = int((modelstart + startTime - diffsecs)/
589     & xx_sstperiod) + 1
590     endrec = int((modelend + startTime - diffsecs + modelstep/2)/
591     & xx_sstperiod) + 2
592 gforget 1.22 endif
593 heimbach 1.20 # else
594     startrec = 1
595     endrec = 1
596     # endif
597     diffrec = endrec - startrec + 1
598     call ctrl_init_ctrlvar (
599     & xx_sst_file, 18, 118, diffrec, startrec, endrec,
600     & snx, sny, 1, 'c', 'xy', mythid )
601    
602     #elif (defined (ALLOW_SST0_CONTROL))
603    
604 heimbach 1.11 call ctrl_init_ctrlvar (
605     & xx_sst_file, 18, 118, 1, 1, 1,
606     & snx, sny, 1, 'c', 'xy', mythid )
607 heimbach 1.20
608     #endif /* ALLOW_SST_CONTROL */
609 heimbach 1.6
610 heimbach 1.16 c----------------------------------------------------------------------
611 heimbach 1.6 c--
612 heimbach 1.20 #if (defined (ALLOW_SSS_CONTROL))
613    
614     # ifdef ALLOW_CAL
615     call cal_FullDate( xx_sssstartdate1, xx_sssstartdate2,
616     & xx_sssstartdate , mythid )
617     call cal_TimePassed( xx_sssstartdate, modelstartdate,
618     & difftime, mythid )
619     call cal_ToSeconds ( difftime, diffsecs, mythid )
620 gforget 1.22 if ( xx_sssperiod .EQ. 0 ) then
621     startrec=1
622     endrec=12
623     else
624 heimbach 1.20 startrec = int((modelstart + startTime - diffsecs)/
625     & xx_sssperiod) + 1
626     endrec = int((modelend + startTime - diffsecs + modelstep/2)/
627     & xx_sssperiod) + 2
628 gforget 1.22 endif
629 heimbach 1.20 # else
630     startrec = 1
631     endrec = 1
632     # endif
633     diffrec = endrec - startrec + 1
634     call ctrl_init_ctrlvar (
635     & xx_sss_file, 19, 119, diffrec, startrec, endrec,
636     & snx, sny, 1, 'c', 'xy', mythid )
637    
638     #elif (defined (ALLOW_SSS0_CONTROL))
639    
640 heimbach 1.11 call ctrl_init_ctrlvar (
641     & xx_sss_file, 19, 119, 1, 1, 1,
642     & snx, sny, 1, 'c', 'xy', mythid )
643 heimbach 1.20
644 heimbach 1.6 #endif /* ALLOW_SSS0_CONTROL */
645    
646 heimbach 1.16 c----------------------------------------------------------------------
647 heimbach 1.6 c--
648 heimbach 1.23 #ifdef ALLOW_DEPTH_CONTROL
649 heimbach 1.11 call ctrl_init_ctrlvar (
650 heimbach 1.23 & xx_depth_file, 20, 120, 1, 1, 1,
651 heimbach 1.11 & snx, sny, 1, 'c', 'xy', mythid )
652 heimbach 1.23 #endif /* ALLOW_DEPTH_CONTROL */
653 heimbach 1.6
654 heimbach 1.16 c----------------------------------------------------------------------
655 heimbach 1.6 c--
656 heimbach 1.5 #ifdef ALLOW_EFLUXY0_CONTROL
657 heimbach 1.11 call ctrl_init_ctrlvar (
658     & xx_efluxy_file, 21, 121, 1, 1, 1,
659     & snx, sny, nr, 's', '3d', mythid )
660 heimbach 1.5 #endif /* ALLOW_EFLUXY0_CONTROL */
661    
662 heimbach 1.16 c----------------------------------------------------------------------
663 heimbach 1.6 c--
664 heimbach 1.5 #ifdef ALLOW_EFLUXP0_CONTROL
665 heimbach 1.11 call ctrl_init_ctrlvar (
666     & xx_efluxp_file, 22, 122, 1, 1, 1,
667     & snx, sny, nr, 'v', '3d', mythid )
668 heimbach 1.5 #endif /* ALLOW_EFLUXP0_CONTROL */
669 heimbach 1.6
670 heimbach 1.16 c----------------------------------------------------------------------
671 heimbach 1.6 c--
672     #ifdef ALLOW_BOTTOMDRAG_CONTROL
673 heimbach 1.11 call ctrl_init_ctrlvar (
674     & xx_bottomdrag_file, 23, 123, 1, 1, 1,
675     & snx, sny, 1, 'c', 'xy', mythid )
676 heimbach 1.6 #endif /* ALLOW_BOTTOMDRAG_CONTROL */
677    
678 heimbach 1.16 c----------------------------------------------------------------------
679 heimbach 1.15 c--
680     #ifdef ALLOW_EDTAUX_CONTROL
681     call ctrl_init_ctrlvar (
682     & xx_edtaux_file, 25, 125, 1, 1, 1,
683     & snx, sny, nr, 'w', '3d', mythid )
684     #endif /* ALLOW_EDTAUX_CONTROL */
685    
686 heimbach 1.16 c----------------------------------------------------------------------
687 heimbach 1.15 c--
688     #ifdef ALLOW_EDTAUY_CONTROL
689     call ctrl_init_ctrlvar (
690     & xx_edtauy_file, 26, 126, 1, 1, 1,
691     & snx, sny, nr, 's', '3d', mythid )
692     #endif /* ALLOW_EDTAUY_CONTROL */
693    
694 heimbach 1.16 c----------------------------------------------------------------------
695     c--
696     #ifdef ALLOW_UVEL0_CONTROL
697     call ctrl_init_ctrlvar (
698     & xx_uvel_file, 27, 127, 1, 1, 1,
699     & snx, sny, nr, 'w', '3d', mythid )
700     #endif /* ALLOW_UVEL0_CONTROL */
701    
702     c----------------------------------------------------------------------
703     c--
704     #ifdef ALLOW_VVEL0_CONTROL
705     call ctrl_init_ctrlvar (
706     & xx_vvel_file, 28, 128, 1, 1, 1,
707     & snx, sny, nr, 's', '3d', mythid )
708     #endif /* ALLOW_VVEL0_CONTROL */
709    
710     c----------------------------------------------------------------------
711     c--
712     #ifdef ALLOW_ETAN0_CONTROL
713     call ctrl_init_ctrlvar (
714     & xx_etan_file, 29, 129, 1, 1, 1,
715     & snx, sny, 1, 'c', 'xy', mythid )
716     #endif /* ALLOW_VVEL0_CONTROL */
717    
718     c----------------------------------------------------------------------
719     c--
720     #ifdef ALLOW_RELAXSST_CONTROL
721     call ctrl_init_ctrlvar (
722     & xx_relaxsst_file, 30, 130, 1, 1, 1,
723     & snx, sny, 1, 'c', 'xy', mythid )
724     #endif /* ALLOW_RELAXSST_CONTROL */
725    
726     c----------------------------------------------------------------------
727     c--
728     #ifdef ALLOW_RELAXSSS_CONTROL
729     call ctrl_init_ctrlvar (
730     & xx_relaxsss_file, 31, 131, 1, 1, 1,
731     & snx, sny, 1, 'c', 'xy', mythid )
732     #endif /* ALLOW_RELAXSSS_CONTROL */
733    
734     c----------------------------------------------------------------------
735 heimbach 1.17 c--
736 heimbach 1.18 #ifdef ALLOW_PRECIP_CONTROL
737     c-- Atmos. precipitation
738    
739     # ifdef ALLOW_CAL
740     call cal_FullDate( xx_precipstartdate1, xx_precipstartdate2,
741     & xx_precipstartdate , mythid )
742     call cal_TimePassed( xx_precipstartdate, modelstartdate,
743     & difftime, mythid )
744     call cal_ToSeconds ( difftime, diffsecs, mythid )
745 gforget 1.22 if ( xx_precipperiod .EQ. 0 ) then
746     startrec=1
747     endrec=12
748     else
749 heimbach 1.18 startrec = int((modelstart + startTime - diffsecs)/
750     & xx_precipperiod) + 1
751     endrec = int((modelend + startTime - diffsecs + modelstep/2)/
752     & xx_precipperiod) + 2
753 gforget 1.22 endif
754 heimbach 1.18 # else
755     startrec = 1
756     endrec = 1
757     # endif
758     diffrec = endrec - startrec + 1
759     call ctrl_init_ctrlvar (
760     & xx_precip_file, 32, 132, diffrec, startrec, endrec,
761     & snx, sny, 1, 'c', 'xy', mythid )
762    
763     #endif /* ALLOW_PRECIP_CONTROL */
764    
765     c----------------------------------------------------------------------
766     c--
767     #ifdef ALLOW_SWFLUX_CONTROL
768     c-- Atmos. swflux
769    
770     # ifdef ALLOW_CAL
771     call cal_FullDate( xx_swfluxstartdate1, xx_swfluxstartdate2,
772     & xx_swfluxstartdate , mythid )
773     call cal_TimePassed( xx_swfluxstartdate, modelstartdate,
774     & difftime, mythid )
775     call cal_ToSeconds ( difftime, diffsecs, mythid )
776 gforget 1.22 if ( xx_swfluxperiod .EQ. 0 ) then
777     startrec=1
778     endrec=12
779     else
780 heimbach 1.18 startrec = int((modelstart + startTime - diffsecs)/
781     & xx_swfluxperiod) + 1
782     endrec = int((modelend + startTime - diffsecs + modelstep/2)/
783     & xx_swfluxperiod) + 2
784 gforget 1.22 endif
785 heimbach 1.18 # else
786     startrec = 1
787     endrec = 1
788     # endif
789     diffrec = endrec - startrec + 1
790 heimbach 1.17 call ctrl_init_ctrlvar (
791 heimbach 1.18 & xx_swflux_file, 33, 133, diffrec, startrec, endrec,
792     & snx, sny, 1, 'c', 'xy', mythid )
793    
794     #endif /* ALLOW_SWFLUX_CONTROL */
795 heimbach 1.17
796     c----------------------------------------------------------------------
797 heimbach 1.19 c--
798     #ifdef ALLOW_SWDOWN_CONTROL
799     c-- Atmos. swdown
800    
801     # ifdef ALLOW_CAL
802     call cal_FullDate( xx_swdownstartdate1, xx_swdownstartdate2,
803     & xx_swdownstartdate , mythid )
804     call cal_TimePassed( xx_swdownstartdate, modelstartdate,
805     & difftime, mythid )
806     call cal_ToSeconds ( difftime, diffsecs, mythid )
807 gforget 1.22 if ( xx_swdownperiod .EQ. 0 ) then
808     startrec=1
809     endrec=12
810     else
811 heimbach 1.19 startrec = int((modelstart + startTime - diffsecs)/
812     & xx_swdownperiod) + 1
813     endrec = int((modelend + startTime - diffsecs + modelstep/2)/
814     & xx_swdownperiod) + 2
815 gforget 1.22 endif
816 heimbach 1.19 # else
817     startrec = 1
818     endrec = 1
819     # endif
820     diffrec = endrec - startrec + 1
821     call ctrl_init_ctrlvar (
822     & xx_swdown_file, 34, 134, diffrec, startrec, endrec,
823     & snx, sny, 1, 'c', 'xy', mythid )
824    
825     #endif /* ALLOW_SWDOWN_CONTROL */
826    
827     c----------------------------------------------------------------------
828 heimbach 1.24 c--
829     #ifdef ALLOW_LWFLUX_CONTROL
830     c-- Atmos. lwflux
831    
832     # ifdef ALLOW_CAL
833     call cal_FullDate( xx_lwfluxstartdate1, xx_lwfluxstartdate2,
834     & xx_lwfluxstartdate , mythid )
835     call cal_TimePassed( xx_lwfluxstartdate, modelstartdate,
836     & difftime, mythid )
837     call cal_ToSeconds ( difftime, diffsecs, mythid )
838     if ( xx_lwfluxperiod .EQ. 0 ) then
839     startrec=1
840     endrec=12
841     else
842     startrec = int((modelstart + startTime - diffsecs)/
843     & xx_lwfluxperiod) + 1
844     endrec = int((modelend + startTime - diffsecs + modelstep/2)/
845     & xx_lwfluxperiod) + 2
846     endif
847     # else
848     startrec = 1
849     endrec = 1
850     # endif
851     diffrec = endrec - startrec + 1
852     call ctrl_init_ctrlvar (
853     & xx_lwflux_file, 35, 135, diffrec, startrec, endrec,
854     & snx, sny, 1, 'c', 'xy', mythid )
855    
856     #endif /* ALLOW_LWFLUX_CONTROL */
857    
858     c----------------------------------------------------------------------
859     c--
860     #ifdef ALLOW_LWDOWN_CONTROL
861     c-- Atmos. lwdown
862    
863     # ifdef ALLOW_CAL
864     call cal_FullDate( xx_lwdownstartdate1, xx_lwdownstartdate2,
865     & xx_lwdownstartdate , mythid )
866     call cal_TimePassed( xx_lwdownstartdate, modelstartdate,
867     & difftime, mythid )
868     call cal_ToSeconds ( difftime, diffsecs, mythid )
869     if ( xx_lwdownperiod .EQ. 0 ) then
870     startrec=1
871     endrec=12
872     else
873     startrec = int((modelstart + startTime - diffsecs)/
874     & xx_lwdownperiod) + 1
875     endrec = int((modelend + startTime - diffsecs + modelstep/2)/
876     & xx_lwdownperiod) + 2
877     endif
878     # else
879     startrec = 1
880     endrec = 1
881     # endif
882     diffrec = endrec - startrec + 1
883     call ctrl_init_ctrlvar (
884     & xx_lwdown_file, 36, 136, diffrec, startrec, endrec,
885     & snx, sny, 1, 'c', 'xy', mythid )
886    
887     #endif /* ALLOW_LWDOWN_CONTROL */
888    
889     c----------------------------------------------------------------------
890     c--
891     #ifdef ALLOW_EVAP_CONTROL
892     c-- Atmos. evap
893    
894     # ifdef ALLOW_CAL
895     call cal_FullDate( xx_evapstartdate1, xx_evapstartdate2,
896     & xx_evapstartdate , mythid )
897     call cal_TimePassed( xx_evapstartdate, modelstartdate,
898     & difftime, mythid )
899     call cal_ToSeconds ( difftime, diffsecs, mythid )
900     if ( xx_evapperiod .EQ. 0 ) then
901     startrec=1
902     endrec=12
903     else
904     startrec = int((modelstart + startTime - diffsecs)/
905     & xx_evapperiod) + 1
906     endrec = int((modelend + startTime - diffsecs + modelstep/2)/
907     & xx_evapperiod) + 2
908     endif
909     # else
910     startrec = 1
911     endrec = 1
912     # endif
913     diffrec = endrec - startrec + 1
914     call ctrl_init_ctrlvar (
915     & xx_evap_file, 37, 137, diffrec, startrec, endrec,
916     & snx, sny, 1, 'c', 'xy', mythid )
917    
918     #endif /* ALLOW_EVAP_CONTROL */
919    
920     c----------------------------------------------------------------------
921     c--
922     #ifdef ALLOW_SNOWPRECIP_CONTROL
923     c-- Atmos. snowprecip
924    
925     # ifdef ALLOW_CAL
926     call cal_FullDate( xx_snowprecipstartdate1,
927     & xx_snowprecipstartdate2, xx_snowprecipstartdate , mythid )
928     call cal_TimePassed( xx_snowprecipstartdate, modelstartdate,
929     & difftime, mythid )
930     call cal_ToSeconds ( difftime, diffsecs, mythid )
931     if ( xx_snowprecipperiod .EQ. 0 ) then
932     startrec=1
933     endrec=12
934     else
935     startrec = int((modelstart + startTime - diffsecs)/
936     & xx_snowprecipperiod) + 1
937     endrec = int((modelend + startTime - diffsecs + modelstep/2)/
938     & xx_snowprecipperiod) + 2
939     endif
940     # else
941     startrec = 1
942     endrec = 1
943     # endif
944     diffrec = endrec - startrec + 1
945     call ctrl_init_ctrlvar (
946     & xx_snowprecip_file, 38, 138, diffrec, startrec, endrec,
947     & snx, sny, 1, 'c', 'xy', mythid )
948    
949     #endif /* ALLOW_SNOWPRECIP_CONTROL */
950    
951     c----------------------------------------------------------------------
952     c--
953     #ifdef ALLOW_APRESSURE_CONTROL
954     c-- Atmos. apressure
955    
956     # ifdef ALLOW_CAL
957     call cal_FullDate( xx_apressurestartdate1,
958     & xx_apressurestartdate2, xx_apressurestartdate , mythid )
959     call cal_TimePassed( xx_apressurestartdate, modelstartdate,
960     & difftime, mythid )
961     call cal_ToSeconds ( difftime, diffsecs, mythid )
962     if ( xx_apressureperiod .EQ. 0 ) then
963     startrec=1
964     endrec=12
965     else
966     startrec = int((modelstart + startTime - diffsecs)/
967     & xx_apressureperiod) + 1
968     endrec = int((modelend + startTime - diffsecs + modelstep/2)/
969     & xx_apressureperiod) + 2
970     endif
971     # else
972     startrec = 1
973     endrec = 1
974     # endif
975     diffrec = endrec - startrec + 1
976     call ctrl_init_ctrlvar (
977     & xx_apressure_file, 39, 139, diffrec, startrec, endrec,
978     & snx, sny, 1, 'c', 'xy', mythid )
979    
980     #endif /* ALLOW_APRESSURE_CONTROL */
981    
982     c----------------------------------------------------------------------
983     c--
984     #ifdef ALLOW_RUNOFF_CONTROL
985     c-- Atmos. runoff
986     startrec = 1
987     endrec = 1
988     diffrec = endrec - startrec + 1
989     call ctrl_init_ctrlvar (
990     & xx_runoff_file, 40, 140, diffrec, startrec, endrec,
991     & snx, sny, 1, 'c', 'xy', mythid )
992     #endif /* ALLOW_RUNOFF_CONTROL */
993    
994     c----------------------------------------------------------------------
995 heimbach 1.16 c----------------------------------------------------------------------
996 heimbach 1.21
997     call ctrl_init_wet( mythid )
998    
999     c----------------------------------------------------------------------
1000 heimbach 1.16 c----------------------------------------------------------------------
1001 heimbach 1.1
1002 heimbach 1.21 do bj = jtlo,jthi
1003     do bi = itlo,ithi
1004     do j = jmin,jmax
1005     do i = imin,imax
1006     wareaunit (i,j,bi,bj) = 1.0
1007     #ifndef ALLOW_ECCO
1008     whflux (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1009     wsflux (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1010     wtauu (i,j,bi,bj) = maskW(i,j,1,bi,bj)
1011     wtauv (i,j,bi,bj) = maskS(i,j,1,bi,bj)
1012     watemp (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1013     waqh (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1014     wprecip (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1015     wswflux (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1016     wswdown (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1017     wuwind (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1018     wvwind (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1019 heimbach 1.24 wlwflux (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1020     wlwdown (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1021     wevap (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1022     wsnowprecip(i,j,bi,bj) = maskC(i,j,1,bi,bj)
1023     wapressure(i,j,bi,bj) = maskC(i,j,1,bi,bj)
1024     wrunoff (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1025 heimbach 1.21 wsst (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1026     wsss (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1027     #endif
1028     enddo
1029     enddo
1030     enddo
1031     enddo
1032    
1033 heimbach 1.1 return
1034     end
1035    

  ViewVC Help
Powered by ViewVC 1.1.22