/[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.20 - (show annotations) (download)
Wed Aug 31 00:03:30 2005 UTC (18 years, 9 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57s_post, checkpoint57y_post, checkpoint57r_post, checkpoint57z_post, checkpoint57t_post, checkpoint57v_post, checkpoint57y_pre, checkpint57u_post, checkpoint57w_post, checkpoint57x_post
Changes since 1.19: +50 -4 lines
Adding time-dependent SST, SSS control.

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

  ViewVC Help
Powered by ViewVC 1.1.22