/[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.21 - (hide annotations) (download)
Thu Jan 5 17:05:42 2006 UTC (18 years, 4 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint58b_post, checkpoint58, checkpoint58a_post, checkpoint58c_post
Changes since 1.20: +46 -4 lines
Bring ctrl up-to-data, in particular w.r.t. wthetaLev, wsaltLev

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

  ViewVC Help
Powered by ViewVC 1.1.22