/[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.33 - (hide annotations) (download)
Thu Oct 15 05:22:29 2009 UTC (14 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint62c, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62h, checkpoint62, checkpoint62b, checkpoint61z, checkpoint61x, checkpoint61y
Changes since 1.32: +13 -1 lines
Initialise xx_dic (still missing adxx_dic)

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

  ViewVC Help
Powered by ViewVC 1.1.22