/[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.23 - (hide annotations) (download)
Wed Jun 7 01:55:13 2006 UTC (17 years, 11 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint58m_post, checkpoint58o_post, checkpoint58p_post, checkpoint58q_post, checkpoint58n_post, checkpoint58k_post, checkpoint58l_post, checkpoint58h_post, checkpoint58j_post, checkpoint58i_post
Changes since 1.22: +4 -10 lines
Modifications for bottom topography control
o replace hFacC by _hFacC at various places
o replace ALLOW_HFACC_CONTROL by ALLOW_DEPTH_CONTROL
o add non-self-adjoint cg2d_nsa
o update autodiff support routines
o re-initialise hfac after ctrl_depth_ini
o works for 5x5 box, doesnt work for global_ocean.90x40x15

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

  ViewVC Help
Powered by ViewVC 1.1.22