/[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.30 - (hide annotations) (download)
Fri May 30 02:48:28 2008 UTC (16 years ago) by gforget
Branch: MAIN
CVS Tags: checkpoint60, checkpoint61, checkpoint61f, checkpoint61n, checkpoint61e, checkpoint61g, checkpoint61d, checkpoint61b, checkpoint61c, checkpoint61a, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61j, checkpoint61k, checkpoint61h, checkpoint61i
Changes since 1.29: +3 -7 lines
o bridging the gap between eddy stress and GM.
  -> eddyTau is replaced with eddyPsi (eddyTau = f x rho0 x eddyPsi)
      along with a change in CPP option (now ALLOW_EDDYPSI).
  -> when using GM w/ GM_AdvForm:
      The total eddy streamfunction (Psi = eddyPsi + K x Slope)
      is applied either in the tracer Eq. or in momentum Eq.
      depending on data.gmredi (intro. GM_InMomAsStress).
  -> ALLOW_EDDYPSI_CONTROL for estimation purpose.
  The key modifications are in model/src/taueddy_external_forcing.F
  pkg/gmredi/gmredi_calc_*F pkg/gmredi/gmredi_*transport.F

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

  ViewVC Help
Powered by ViewVC 1.1.22