/[MITgcm]/MITgcm/pkg/ctrl/ctrl_init.F
ViewVC logotype

Contents of /MITgcm/pkg/ctrl/ctrl_init.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.31 - (show annotations) (download)
Fri May 29 06:12:05 2009 UTC (15 years ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint61q, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p
Changes since 1.30: +1 -2 lines
Bug fix (reported by jmc): remove un-assigned vars.

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

  ViewVC Help
Powered by ViewVC 1.1.22