/[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.34 - (show annotations) (download)
Tue Jul 13 00:02:06 2010 UTC (13 years, 10 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint62q, checkpoint62p, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l
Changes since 1.33: +9 -1 lines
ALLOW_ROTATE_UV_CONTROLS: when defined, we
rotate wind/stress controls adjustments
from Eastward/Northward to model grid directions.

1 C
2 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_init.F,v 1.33 2009/10/15 05:22:29 heimbach 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 startrec = int((modelstart - diffsecs)/xx_obcsnperiod) + 1
471 startrec = (startrec - 1)*nobcs + 1
472 endrec = int((modelend - diffsecs)/xx_obcsnperiod) + 2
473 endrec = (endrec - startrec + 1)*nobcs
474 # else
475 startrec = 1
476 endrec = 1
477 # endif
478 diffrec = endrec
479 call ctrl_init_ctrlvar (
480 & xx_obcsn_file, 11, 111, diffrec, startrec, endrec,
481 & snx, 1, nr, 'm', 'xz', mythid )
482
483 #endif /* ALLOW_OBCSN_CONTROL */
484
485 c----------------------------------------------------------------------
486 c--
487 #ifdef ALLOW_OBCSS_CONTROL
488 c-- Southern obc.
489
490 # ifdef ALLOW_CAL
491 call cal_FullDate( xx_obcssstartdate1, xx_obcssstartdate2,
492 & xx_obcssstartdate, mythid )
493 call cal_TimePassed( xx_obcssstartdate, modelstartdate,
494 & difftime, mythid )
495 call cal_ToSeconds ( difftime, diffsecs, mythid )
496 startrec = int((modelstart - diffsecs)/xx_obcssperiod) + 1
497 startrec = (startrec - 1)*nobcs + 1
498 endrec = int((modelend - diffsecs)/xx_obcssperiod) + 2
499 endrec = (endrec - startrec + 1)*nobcs
500 # else
501 startrec = 1
502 endrec = 1
503 # endif
504 diffrec = endrec
505 call ctrl_init_ctrlvar (
506 & xx_obcss_file, 12, 112, diffrec, startrec, endrec,
507 & snx, 1, nr, 'm', 'xz', mythid )
508
509 #endif /* ALLOW_OBCSS_CONTROL */
510
511 c----------------------------------------------------------------------
512 c--
513 #ifdef ALLOW_OBCSW_CONTROL
514 c-- Western obc.
515
516 # ifdef ALLOW_CAL
517 call cal_FullDate( xx_obcswstartdate1, xx_obcswstartdate2,
518 & xx_obcswstartdate, mythid )
519 call cal_TimePassed( xx_obcswstartdate, modelstartdate,
520 & difftime, mythid )
521 call cal_ToSeconds ( difftime, diffsecs, mythid )
522 startrec = int((modelstart - diffsecs)/xx_obcswperiod) + 1
523 startrec = (startrec - 1)*nobcs + 1
524 endrec = int((modelend - diffsecs)/xx_obcswperiod) + 2
525 endrec = (endrec - startrec + 1)*nobcs
526 # else
527 startrec = 1
528 endrec = 1
529 # endif
530 diffrec = endrec
531 call ctrl_init_ctrlvar (
532 & xx_obcsw_file, 13, 113, diffrec, startrec, endrec,
533 & 1, sny, nr, 'm', 'yz', mythid )
534
535 #endif /* ALLOW_OBCSW_CONTROL */
536
537 c----------------------------------------------------------------------
538 c--
539 #ifdef ALLOW_OBCSE_CONTROL
540 c-- Eastern obc.
541
542 # ifdef ALLOW_CAL
543 call cal_FullDate( xx_obcsestartdate1, xx_obcsestartdate2,
544 & xx_obcsestartdate, mythid )
545 call cal_TimePassed( xx_obcsestartdate, modelstartdate,
546 & difftime, mythid )
547 call cal_ToSeconds ( difftime, diffsecs, mythid )
548 startrec = int((modelstart - diffsecs)/xx_obcseperiod) + 1
549 startrec = (startrec - 1)*nobcs + 1
550 endrec = int((modelend - diffsecs)/xx_obcseperiod) + 2
551 endrec = (endrec - startrec + 1)*nobcs
552 # else
553 startrec = 1
554 endrec = 1
555 # endif
556 diffrec = endrec
557 call ctrl_init_ctrlvar (
558 & xx_obcse_file, 14, 114, diffrec, startrec, endrec,
559 & 1, sny, nr, 'm', 'yz', mythid )
560
561 #endif /* ALLOW_OBCSE_CONTROL */
562
563 c----------------------------------------------------------------------
564 c--
565 #ifdef ALLOW_DIFFKR_CONTROL
566 call ctrl_init_ctrlvar (
567 & xx_diffkr_file, 15, 115, 1, 1, 1,
568 & snx, sny, nr, 'c', '3d', mythid )
569 #endif /* ALLOW_DIFFKR_CONTROL */
570
571 c----------------------------------------------------------------------
572 c--
573 #ifdef ALLOW_KAPGM_CONTROL
574 call ctrl_init_ctrlvar (
575 & xx_kapgm_file, 16, 116, 1, 1, 1,
576 & snx, sny, nr, 'c', '3d', mythid )
577 #endif /* ALLOW_KAPGM_CONTROL */
578
579 c----------------------------------------------------------------------
580 c--
581 #ifdef ALLOW_TR10_CONTROL
582 call ctrl_init_ctrlvar (
583 & xx_tr1_file, 17, 117, 1, 1, 1,
584 & snx, sny, nr, 'c', '3d', mythid )
585 #endif /* ALLOW_TR10_CONTROL */
586
587 c----------------------------------------------------------------------
588 c--
589 #if (defined (ALLOW_SST_CONTROL))
590
591 # ifdef ALLOW_CAL
592 call cal_FullDate( xx_sststartdate1, xx_sststartdate2,
593 & xx_sststartdate , mythid )
594 call cal_TimePassed( xx_sststartdate, modelstartdate,
595 & difftime, mythid )
596 call cal_ToSeconds ( difftime, diffsecs, mythid )
597 if ( xx_sstperiod .EQ. 0 ) then
598 startrec=1
599 endrec=12
600 else
601 startrec = int((modelstart + startTime - diffsecs)/
602 & xx_sstperiod) + 1
603 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
604 & xx_sstperiod) + 2
605 endif
606 # else
607 startrec = 1
608 endrec = 1
609 # endif
610 diffrec = endrec - startrec + 1
611 call ctrl_init_ctrlvar (
612 & xx_sst_file, 18, 118, diffrec, startrec, endrec,
613 & snx, sny, 1, 'c', 'xy', mythid )
614
615 #elif (defined (ALLOW_SST0_CONTROL))
616
617 call ctrl_init_ctrlvar (
618 & xx_sst_file, 18, 118, 1, 1, 1,
619 & snx, sny, 1, 'c', 'xy', mythid )
620
621 #endif /* ALLOW_SST_CONTROL */
622
623 c----------------------------------------------------------------------
624 c--
625 #if (defined (ALLOW_SSS_CONTROL))
626
627 # ifdef ALLOW_CAL
628 call cal_FullDate( xx_sssstartdate1, xx_sssstartdate2,
629 & xx_sssstartdate , mythid )
630 call cal_TimePassed( xx_sssstartdate, modelstartdate,
631 & difftime, mythid )
632 call cal_ToSeconds ( difftime, diffsecs, mythid )
633 if ( xx_sssperiod .EQ. 0 ) then
634 startrec=1
635 endrec=12
636 else
637 startrec = int((modelstart + startTime - diffsecs)/
638 & xx_sssperiod) + 1
639 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
640 & xx_sssperiod) + 2
641 endif
642 # else
643 startrec = 1
644 endrec = 1
645 # endif
646 diffrec = endrec - startrec + 1
647 call ctrl_init_ctrlvar (
648 & xx_sss_file, 19, 119, diffrec, startrec, endrec,
649 & snx, sny, 1, 'c', 'xy', mythid )
650
651 #elif (defined (ALLOW_SSS0_CONTROL))
652
653 call ctrl_init_ctrlvar (
654 & xx_sss_file, 19, 119, 1, 1, 1,
655 & snx, sny, 1, 'c', 'xy', mythid )
656
657 #endif /* ALLOW_SSS0_CONTROL */
658
659 c----------------------------------------------------------------------
660 c--
661 #ifdef ALLOW_DEPTH_CONTROL
662 call ctrl_init_ctrlvar (
663 & xx_depth_file, 20, 120, 1, 1, 1,
664 & snx, sny, 1, 'c', 'xy', mythid )
665 #endif /* ALLOW_DEPTH_CONTROL */
666
667 c----------------------------------------------------------------------
668 c--
669 #ifdef ALLOW_EFLUXY0_CONTROL
670 call ctrl_init_ctrlvar (
671 & xx_efluxy_file, 21, 121, 1, 1, 1,
672 & snx, sny, nr, 's', '3d', mythid )
673 #endif /* ALLOW_EFLUXY0_CONTROL */
674
675 c----------------------------------------------------------------------
676 c--
677 #ifdef ALLOW_EFLUXP0_CONTROL
678 call ctrl_init_ctrlvar (
679 & xx_efluxp_file, 22, 122, 1, 1, 1,
680 & snx, sny, nr, 'v', '3d', mythid )
681 #endif /* ALLOW_EFLUXP0_CONTROL */
682
683 c----------------------------------------------------------------------
684 c--
685 #ifdef ALLOW_BOTTOMDRAG_CONTROL
686 call ctrl_init_ctrlvar (
687 & xx_bottomdrag_file, 23, 123, 1, 1, 1,
688 & snx, sny, 1, 'c', 'xy', mythid )
689 #endif /* ALLOW_BOTTOMDRAG_CONTROL */
690
691 c----------------------------------------------------------------------
692 c--
693 #ifdef ALLOW_HFLUXM_CONTROL
694 call ctrl_init_ctrlvar (
695 & xx_hfluxm_file, 24, 124, 1, 1, 1,
696 & snx, sny, 1, 'c', 'xy', mythid )
697 #endif /* ALLOW_HFLUXM_CONTROL */
698
699 c----------------------------------------------------------------------
700 c--
701 #ifdef ALLOW_EDDYPSI_CONTROL
702 call ctrl_init_ctrlvar (
703 & xx_edtaux_file, 25, 125, 1, 1, 1,
704 & snx, sny, nr, 'w', '3d', mythid )
705
706 call ctrl_init_ctrlvar (
707 & xx_edtauy_file, 26, 126, 1, 1, 1,
708 & snx, sny, nr, 's', '3d', mythid )
709 #endif /* ALLOW_EDDYPSI_CONTROL */
710
711 c----------------------------------------------------------------------
712 c--
713 #ifdef ALLOW_UVEL0_CONTROL
714 call ctrl_init_ctrlvar (
715 & xx_uvel_file, 27, 127, 1, 1, 1,
716 & snx, sny, nr, 'w', '3d', mythid )
717 #endif /* ALLOW_UVEL0_CONTROL */
718
719 c----------------------------------------------------------------------
720 c--
721 #ifdef ALLOW_VVEL0_CONTROL
722 call ctrl_init_ctrlvar (
723 & xx_vvel_file, 28, 128, 1, 1, 1,
724 & snx, sny, nr, 's', '3d', mythid )
725 #endif /* ALLOW_VVEL0_CONTROL */
726
727 c----------------------------------------------------------------------
728 c--
729 #ifdef ALLOW_ETAN0_CONTROL
730 call ctrl_init_ctrlvar (
731 & xx_etan_file, 29, 129, 1, 1, 1,
732 & snx, sny, 1, 'c', 'xy', mythid )
733 #endif /* ALLOW_VVEL0_CONTROL */
734
735 c----------------------------------------------------------------------
736 c--
737 #ifdef ALLOW_GEN2D_CONTROL
738 call ctrl_init_ctrlvar (
739 & xx_gen2d_file, 30, 130, 1, 1, 1,
740 & snx, sny, 1, 'c', 'xy', mythid )
741 #endif /* ALLOW_GEN2D_CONTROL */
742
743 c----------------------------------------------------------------------
744 c--
745 #ifdef ALLOW_GEN3D_CONTROL
746 call ctrl_init_ctrlvar (
747 & xx_gen3d_file, 31, 131, 1, 1, 1,
748 & snx, sny, nr, 'c', '3d', mythid )
749 #endif /* ALLOW_GEN3D_CONTROL */
750
751 c----------------------------------------------------------------------
752 c--
753 #ifdef ALLOW_PRECIP_CONTROL
754 c-- Atmos. precipitation
755
756 # ifdef ALLOW_CAL
757 call cal_FullDate( xx_precipstartdate1, xx_precipstartdate2,
758 & xx_precipstartdate , mythid )
759 call cal_TimePassed( xx_precipstartdate, modelstartdate,
760 & difftime, mythid )
761 call cal_ToSeconds ( difftime, diffsecs, mythid )
762 if ( xx_precipperiod .EQ. 0 ) then
763 startrec=1
764 endrec=12
765 else
766 startrec = int((modelstart + startTime - diffsecs)/
767 & xx_precipperiod) + 1
768 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
769 & xx_precipperiod) + 2
770 endif
771 # else
772 startrec = 1
773 endrec = 1
774 # endif
775 diffrec = endrec - startrec + 1
776 call ctrl_init_ctrlvar (
777 & xx_precip_file, 32, 132, diffrec, startrec, endrec,
778 & snx, sny, 1, 'c', 'xy', mythid )
779
780 #endif /* ALLOW_PRECIP_CONTROL */
781
782 c----------------------------------------------------------------------
783 c--
784 #ifdef ALLOW_SWFLUX_CONTROL
785 c-- Atmos. swflux
786
787 # ifdef ALLOW_CAL
788 call cal_FullDate( xx_swfluxstartdate1, xx_swfluxstartdate2,
789 & xx_swfluxstartdate , mythid )
790 call cal_TimePassed( xx_swfluxstartdate, modelstartdate,
791 & difftime, mythid )
792 call cal_ToSeconds ( difftime, diffsecs, mythid )
793 if ( xx_swfluxperiod .EQ. 0 ) then
794 startrec=1
795 endrec=12
796 else
797 startrec = int((modelstart + startTime - diffsecs)/
798 & xx_swfluxperiod) + 1
799 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
800 & xx_swfluxperiod) + 2
801 endif
802 # else
803 startrec = 1
804 endrec = 1
805 # endif
806 diffrec = endrec - startrec + 1
807 call ctrl_init_ctrlvar (
808 & xx_swflux_file, 33, 133, diffrec, startrec, endrec,
809 & snx, sny, 1, 'c', 'xy', mythid )
810
811 #endif /* ALLOW_SWFLUX_CONTROL */
812
813 c----------------------------------------------------------------------
814 c--
815 #ifdef ALLOW_SWDOWN_CONTROL
816 c-- Atmos. swdown
817
818 # ifdef ALLOW_CAL
819 call cal_FullDate( xx_swdownstartdate1, xx_swdownstartdate2,
820 & xx_swdownstartdate , mythid )
821 call cal_TimePassed( xx_swdownstartdate, modelstartdate,
822 & difftime, mythid )
823 call cal_ToSeconds ( difftime, diffsecs, mythid )
824 if ( xx_swdownperiod .EQ. 0 ) then
825 startrec=1
826 endrec=12
827 else
828 startrec = int((modelstart + startTime - diffsecs)/
829 & xx_swdownperiod) + 1
830 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
831 & xx_swdownperiod) + 2
832 endif
833 # else
834 startrec = 1
835 endrec = 1
836 # endif
837 diffrec = endrec - startrec + 1
838 call ctrl_init_ctrlvar (
839 & xx_swdown_file, 34, 134, diffrec, startrec, endrec,
840 & snx, sny, 1, 'c', 'xy', mythid )
841
842 #endif /* ALLOW_SWDOWN_CONTROL */
843
844 c----------------------------------------------------------------------
845 c--
846 #ifdef ALLOW_LWFLUX_CONTROL
847 c-- Atmos. lwflux
848
849 # ifdef ALLOW_CAL
850 call cal_FullDate( xx_lwfluxstartdate1, xx_lwfluxstartdate2,
851 & xx_lwfluxstartdate , mythid )
852 call cal_TimePassed( xx_lwfluxstartdate, modelstartdate,
853 & difftime, mythid )
854 call cal_ToSeconds ( difftime, diffsecs, mythid )
855 if ( xx_lwfluxperiod .EQ. 0 ) then
856 startrec=1
857 endrec=12
858 else
859 startrec = int((modelstart + startTime - diffsecs)/
860 & xx_lwfluxperiod) + 1
861 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
862 & xx_lwfluxperiod) + 2
863 endif
864 # else
865 startrec = 1
866 endrec = 1
867 # endif
868 diffrec = endrec - startrec + 1
869 call ctrl_init_ctrlvar (
870 & xx_lwflux_file, 35, 135, diffrec, startrec, endrec,
871 & snx, sny, 1, 'c', 'xy', mythid )
872
873 #endif /* ALLOW_LWFLUX_CONTROL */
874
875 c----------------------------------------------------------------------
876 c--
877 #ifdef ALLOW_LWDOWN_CONTROL
878 c-- Atmos. lwdown
879
880 # ifdef ALLOW_CAL
881 call cal_FullDate( xx_lwdownstartdate1, xx_lwdownstartdate2,
882 & xx_lwdownstartdate , mythid )
883 call cal_TimePassed( xx_lwdownstartdate, modelstartdate,
884 & difftime, mythid )
885 call cal_ToSeconds ( difftime, diffsecs, mythid )
886 if ( xx_lwdownperiod .EQ. 0 ) then
887 startrec=1
888 endrec=12
889 else
890 startrec = int((modelstart + startTime - diffsecs)/
891 & xx_lwdownperiod) + 1
892 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
893 & xx_lwdownperiod) + 2
894 endif
895 # else
896 startrec = 1
897 endrec = 1
898 # endif
899 diffrec = endrec - startrec + 1
900 call ctrl_init_ctrlvar (
901 & xx_lwdown_file, 36, 136, diffrec, startrec, endrec,
902 & snx, sny, 1, 'c', 'xy', mythid )
903
904 #endif /* ALLOW_LWDOWN_CONTROL */
905
906 c----------------------------------------------------------------------
907 c--
908 #ifdef ALLOW_EVAP_CONTROL
909 c-- Atmos. evap
910
911 # ifdef ALLOW_CAL
912 call cal_FullDate( xx_evapstartdate1, xx_evapstartdate2,
913 & xx_evapstartdate , mythid )
914 call cal_TimePassed( xx_evapstartdate, modelstartdate,
915 & difftime, mythid )
916 call cal_ToSeconds ( difftime, diffsecs, mythid )
917 if ( xx_evapperiod .EQ. 0 ) then
918 startrec=1
919 endrec=12
920 else
921 startrec = int((modelstart + startTime - diffsecs)/
922 & xx_evapperiod) + 1
923 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
924 & xx_evapperiod) + 2
925 endif
926 # else
927 startrec = 1
928 endrec = 1
929 # endif
930 diffrec = endrec - startrec + 1
931 call ctrl_init_ctrlvar (
932 & xx_evap_file, 37, 137, diffrec, startrec, endrec,
933 & snx, sny, 1, 'c', 'xy', mythid )
934
935 #endif /* ALLOW_EVAP_CONTROL */
936
937 c----------------------------------------------------------------------
938 c--
939 #ifdef ALLOW_SNOWPRECIP_CONTROL
940 c-- Atmos. snowprecip
941
942 # ifdef ALLOW_CAL
943 call cal_FullDate( xx_snowprecipstartdate1,
944 & xx_snowprecipstartdate2, xx_snowprecipstartdate , mythid )
945 call cal_TimePassed( xx_snowprecipstartdate, modelstartdate,
946 & difftime, mythid )
947 call cal_ToSeconds ( difftime, diffsecs, mythid )
948 if ( xx_snowprecipperiod .EQ. 0 ) then
949 startrec=1
950 endrec=12
951 else
952 startrec = int((modelstart + startTime - diffsecs)/
953 & xx_snowprecipperiod) + 1
954 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
955 & xx_snowprecipperiod) + 2
956 endif
957 # else
958 startrec = 1
959 endrec = 1
960 # endif
961 diffrec = endrec - startrec + 1
962 call ctrl_init_ctrlvar (
963 & xx_snowprecip_file, 38, 138, diffrec, startrec, endrec,
964 & snx, sny, 1, 'c', 'xy', mythid )
965
966 #endif /* ALLOW_SNOWPRECIP_CONTROL */
967
968 c----------------------------------------------------------------------
969 c--
970 #ifdef ALLOW_APRESSURE_CONTROL
971 c-- Atmos. apressure
972
973 # ifdef ALLOW_CAL
974 call cal_FullDate( xx_apressurestartdate1,
975 & xx_apressurestartdate2, xx_apressurestartdate , mythid )
976 call cal_TimePassed( xx_apressurestartdate, modelstartdate,
977 & difftime, mythid )
978 call cal_ToSeconds ( difftime, diffsecs, mythid )
979 if ( xx_apressureperiod .EQ. 0 ) then
980 startrec=1
981 endrec=12
982 else
983 startrec = int((modelstart + startTime - diffsecs)/
984 & xx_apressureperiod) + 1
985 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
986 & xx_apressureperiod) + 2
987 endif
988 # else
989 startrec = 1
990 endrec = 1
991 # endif
992 diffrec = endrec - startrec + 1
993 call ctrl_init_ctrlvar (
994 & xx_apressure_file, 39, 139, diffrec, startrec, endrec,
995 & snx, sny, 1, 'c', 'xy', mythid )
996
997 #endif /* ALLOW_APRESSURE_CONTROL */
998
999 c----------------------------------------------------------------------
1000 c--
1001 #ifdef ALLOW_RUNOFF_CONTROL
1002 c-- Atmos. runoff
1003 startrec = 1
1004 endrec = 1
1005 diffrec = endrec - startrec + 1
1006 call ctrl_init_ctrlvar (
1007 & xx_runoff_file, 40, 140, diffrec, startrec, endrec,
1008 & snx, sny, 1, 'c', 'xy', mythid )
1009 #endif /* ALLOW_RUNOFF_CONTROL */
1010
1011 c----------------------------------------------------------------------
1012 c--
1013 #ifdef ALLOW_SIAREA_CONTROL
1014 startrec = 1
1015 endrec = 1
1016 diffrec = endrec - startrec + 1
1017 call ctrl_init_ctrlvar (
1018 & xx_siarea_file, 41, 141, diffrec, startrec, endrec,
1019 & snx, sny, 1, 'c', 'xy', mythid )
1020 #endif /* ALLOW_siarea_CONTROL */
1021
1022 c----------------------------------------------------------------------
1023 c--
1024 #ifdef ALLOW_SIHEFF_CONTROL
1025 startrec = 1
1026 endrec = 1
1027 diffrec = endrec - startrec + 1
1028 call ctrl_init_ctrlvar (
1029 & xx_siheff_file, 42, 142, diffrec, startrec, endrec,
1030 & snx, sny, 1, 'c', 'xy', mythid )
1031 #endif /* ALLOW_siheff_CONTROL */
1032
1033 c----------------------------------------------------------------------
1034 c--
1035 #ifdef ALLOW_SIHSNOW_CONTROL
1036 startrec = 1
1037 endrec = 1
1038 diffrec = endrec - startrec + 1
1039 call ctrl_init_ctrlvar (
1040 & xx_sihsnow_file, 43, 143, diffrec, startrec, endrec,
1041 & snx, sny, 1, 'c', 'xy', mythid )
1042 #endif /* ALLOW_sihsnow_CONTROL */
1043
1044
1045 c----------------------------------------------------------------------
1046 c--
1047 #ifdef ALLOW_KAPREDI_CONTROL
1048 call ctrl_init_ctrlvar (
1049 & xx_kapredi_file, 44, 144, 1, 1, 1,
1050 & snx, sny, nr, 'c', '3d', mythid )
1051 #endif /* ALLOW_KAPREDI_CONTROL */
1052
1053 c----------------------------------------------------------------------
1054 c----------------------------------------------------------------------
1055
1056 call ctrl_init_wet( mythid )
1057
1058 c----------------------------------------------------------------------
1059 c----------------------------------------------------------------------
1060
1061 #ifdef ALLOW_DIC_CONTROL
1062 do i = 1, dic_n_control
1063 xx_dic(i) = 0. _d 0
1064 enddo
1065 #endif
1066
1067 c----------------------------------------------------------------------
1068 c----------------------------------------------------------------------
1069
1070 do bj = jtlo,jthi
1071 do bi = itlo,ithi
1072 do j = jmin,jmax
1073 do i = imin,imax
1074 wareaunit (i,j,bi,bj) = 1.0
1075 #ifndef ALLOW_ECCO
1076 whflux (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1077 wsflux (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1078 wtauu (i,j,bi,bj) = maskW(i,j,1,bi,bj)
1079 wtauv (i,j,bi,bj) = maskS(i,j,1,bi,bj)
1080 watemp (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1081 waqh (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1082 wprecip (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1083 wswflux (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1084 wswdown (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1085 wuwind (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1086 wvwind (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1087 wlwflux (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1088 wlwdown (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1089 wevap (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1090 wsnowprecip(i,j,bi,bj) = maskC(i,j,1,bi,bj)
1091 wapressure(i,j,bi,bj) = maskC(i,j,1,bi,bj)
1092 wrunoff (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1093 wsst (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1094 wsss (i,j,bi,bj) = maskC(i,j,1,bi,bj)
1095 #endif
1096 enddo
1097 enddo
1098 enddo
1099 enddo
1100
1101 return
1102 end
1103

  ViewVC Help
Powered by ViewVC 1.1.22