/[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.35 - (show annotations) (download)
Wed Jan 19 08:29:35 2011 UTC (13 years, 4 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint62u, checkpoint62t, checkpoint62s, checkpoint62r
Changes since 1.34: +69 -29 lines
- add code that allows the use of ALLOW_OBCS?_CONTROL without exf/cal
- fix initialisation (computation of diffrec) for obcs-ctrl
- handle case of obcs?period == 0

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

  ViewVC Help
Powered by ViewVC 1.1.22