/[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.38 - (show annotations) (download)
Wed Apr 20 19:13:40 2011 UTC (13 years ago) by mmazloff
Branch: MAIN
CVS Tags: checkpoint62w
Changes since 1.37: +37 -1 lines
Ability to decompose obcs controls into modes using ALLOW_OBCS_CONTROL_MODES

1 C
2 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_init.F,v 1.37 2011/03/15 16:53:05 mlosch 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
90 integer itlo,ithi
91 integer jtlo,jthi
92 integer jmin,jmax
93 integer imin,imax
94
95 integer ivar
96 integer startrec
97 integer endrec
98 integer diffrec
99
100 #ifdef ALLOW_OBCS_CONTROL_MODES
101 INTEGER k,length_of_rec,dUnit
102 INTEGER MDS_RECLEN
103 EXTERNAL MDS_RECLEN
104 #endif
105
106 c == external ==
107
108 integer ilnblnk
109 external ilnblnk
110
111 c == end of interface ==
112
113 jtlo = mybylo(mythid)
114 jthi = mybyhi(mythid)
115 itlo = mybxlo(mythid)
116 ithi = mybxhi(mythid)
117 jmin = 1-oly
118 jmax = sny+oly
119 imin = 1-olx
120 imax = snx+olx
121
122 c-- Set default values.
123 do ivar = 1,maxcvars
124 ncvarindex(ivar) = -1
125 ncvarrecs(ivar) = 0
126 ncvarxmax(ivar) = 0
127 ncvarymax(ivar) = 0
128 ncvarnrmax(ivar) = 0
129 ncvargrd(ivar) = '?'
130 enddo
131
132 _BARRIER
133
134 c-- =====================
135 c-- Initial state fields.
136 c-- =====================
137
138 cph(
139 cph index 7-10 reserved for atmos. state,
140 cph index 11-14 reserved for open boundaries,
141 cph index 15-16 reserved for mixing coeff.
142 cph index 17 reserved for passive tracer TR1
143 cph index 18,19 reserved for sst, sss
144 cph index 20 for hFacC
145 cph index 21-22 for efluxy, efluxp
146 cph index 23 for bottom drag
147 cph index 24
148 cph index 25-26 for edtaux, edtauy
149 cph index 27-29 for uvel0, vvel0, etan0
150 cph index 30-31 for generic 2d, 3d field
151 cph index 32 reserved for precip (atmos. state)
152 cph index 33 reserved for swflux (atmos. state)
153 cph index 34 reserved for swdown (atmos. state)
154 cph 35 lwflux
155 cph 36 lwdown
156 cph 37 evap
157 cph 38 snowprecip
158 cph 39 apressure
159 cph 40 runoff
160 cph 41 seaice SIAREA
161 cph 42 seaice SIHEFF
162 cph 43 seaice SIHSNOW
163 cph)
164
165 c----------------------------------------------------------------------
166 c--
167 #ifdef ALLOW_THETA0_CONTROL
168 c-- Initial state temperature contribution.
169 call ctrl_init_ctrlvar (
170 & xx_theta_file, 1, 101, 1, 1, 1,
171 & snx, sny, nr, 'c', '3d', mythid )
172 #endif /* ALLOW_THETA0_CONTROL */
173
174 c----------------------------------------------------------------------
175 c--
176 #ifdef ALLOW_SALT0_CONTROL
177 c-- Initial state salinity contribution.
178 call ctrl_init_ctrlvar (
179 & xx_salt_file, 2, 102, 1, 1, 1,
180 & snx, sny, nr, 'c', '3d', mythid )
181 #endif /* ALLOW_SALT0_CONTROL */
182
183 c-- ===========================
184 c-- Surface flux contributions.
185 c-- ===========================
186
187 c----------------------------------------------------------------------
188 c--
189 #if (defined (ALLOW_HFLUX_CONTROL))
190 c-- Heat flux.
191 call ctrl_init_rec (
192 I xx_hfluxstartdate1, xx_hfluxstartdate2, xx_hfluxperiod, 1,
193 O xx_hfluxstartdate, diffrec, startrec, endrec,
194 I mythid )
195 call ctrl_init_ctrlvar (
196 & xx_hflux_file, 3, 103, diffrec, startrec, endrec,
197 & snx, sny, 1, 'c', 'xy', mythid )
198
199 #elif (defined (ALLOW_ATEMP_CONTROL))
200 c-- Atmos. temperature
201 call ctrl_init_rec (
202 I xx_atempstartdate1, xx_atempstartdate2, xx_atempperiod, 1,
203 O xx_atempstartdate, diffrec, startrec, endrec,
204 I mythid )
205 call ctrl_init_ctrlvar (
206 & xx_atemp_file, 7, 107, diffrec, startrec, endrec,
207 & snx, sny, 1, 'c', 'xy', mythid )
208
209 #elif (defined (ALLOW_HFLUX0_CONTROL))
210 c-- initial forcing only
211 call ctrl_init_ctrlvar (
212 & xx_hflux_file, 3, 103, 1, 1, 1,
213 & snx, sny, 1, 'c', 'xy', mythid )
214
215 #endif /* ALLOW_HFLUX_CONTROL */
216
217 c----------------------------------------------------------------------
218 c--
219 #if (defined (ALLOW_SFLUX_CONTROL))
220 c-- Salt flux.
221 call ctrl_init_rec (
222 I xx_sfluxstartdate1, xx_sfluxstartdate2, xx_sfluxperiod, 1,
223 O xx_sfluxstartdate, diffrec, startrec, endrec,
224 I mythid )
225 call ctrl_init_ctrlvar (
226 & xx_sflux_file, 4, 104, diffrec, startrec, endrec,
227 & snx, sny, 1, 'c', 'xy', mythid )
228
229 #elif (defined (ALLOW_AQH_CONTROL))
230 c-- Atmos. humidity
231 call ctrl_init_rec (
232 I xx_aqhstartdate1, xx_aqhstartdate2, xx_aqhperiod, 1,
233 O xx_aqhstartdate, diffrec, startrec, endrec,
234 I mythid )
235 call ctrl_init_ctrlvar (
236 & xx_aqh_file, 8, 108, diffrec, startrec, endrec,
237 & snx, sny, 1, 'c', 'xy', mythid )
238
239 #elif (defined (ALLOW_SFLUX0_CONTROL))
240 c-- initial forcing only
241 call ctrl_init_ctrlvar (
242 & xx_sflux_file, 4, 104, 1, 1, 1,
243 & snx, sny, 1, 'c', 'xy', mythid )
244
245 #endif /* ALLOW_SFLUX_CONTROL */
246
247 c----------------------------------------------------------------------
248 c--
249 #if (defined (ALLOW_USTRESS_CONTROL))
250 c-- Zonal wind stress.
251 call ctrl_init_rec (
252 I xx_tauustartdate1, xx_tauustartdate2, xx_tauuperiod, 1,
253 O xx_tauustartdate, diffrec, startrec, endrec,
254 I mythid )
255 call ctrl_init_ctrlvar (
256 & xx_tauu_file, 5, 105, diffrec, startrec, endrec,
257 #ifndef ALLOW_ROTATE_UV_CONTROLS
258 & snx, sny, 1, 'w', 'xy', mythid )
259 #else
260 & snx, sny, 1, 'c', 'xy', mythid )
261 #endif
262
263 #elif (defined (ALLOW_UWIND_CONTROL))
264 c-- Zonal wind speed.
265 call ctrl_init_rec (
266 I xx_uwindstartdate1, xx_uwindstartdate2, xx_uwindperiod, 1,
267 O xx_uwindstartdate, diffrec, startrec, endrec,
268 I mythid )
269 call ctrl_init_ctrlvar (
270 & xx_uwind_file, 9, 109, diffrec, startrec, endrec,
271 & snx, sny, 1, 'c', 'xy', mythid )
272
273 #elif (defined (ALLOW_TAUU0_CONTROL))
274 c-- initial forcing only
275 call ctrl_init_ctrlvar (
276 & xx_tauu_file, 5, 105, 1, 1, 1,
277 & snx, sny, 1, 'w', 'xy', mythid )
278
279 #endif /* ALLOW_USTRESS_CONTROL */
280
281 c----------------------------------------------------------------------
282 c--
283 #if (defined (ALLOW_VSTRESS_CONTROL))
284 c-- Meridional wind stress.
285 call ctrl_init_rec (
286 I xx_tauvstartdate1, xx_tauvstartdate2, xx_tauvperiod, 1,
287 O xx_tauvstartdate, diffrec, startrec, endrec,
288 I mythid )
289 call ctrl_init_ctrlvar (
290 & xx_tauv_file, 6, 106, diffrec, startrec, endrec,
291 #ifndef ALLOW_ROTATE_UV_CONTROLS
292 & snx, sny, 1, 's', 'xy', mythid )
293 #else
294 & snx, sny, 1, 'c', 'xy', mythid )
295 #endif
296
297 #elif (defined (ALLOW_VWIND_CONTROL))
298 c-- Meridional wind speed.
299 call ctrl_init_rec (
300 I xx_vwindstartdate1, xx_vwindstartdate2, xx_vwindperiod, 1,
301 O xx_vwindstartdate, diffrec, startrec, endrec,
302 I mythid )
303 call ctrl_init_ctrlvar (
304 & xx_vwind_file, 10, 110, diffrec, startrec, endrec,
305 & snx, sny, 1, 'c', 'xy', mythid )
306
307 #elif (defined (ALLOW_TAUV0_CONTROL))
308 c-- initial forcing only
309 call ctrl_init_ctrlvar (
310 & xx_tauv_file, 6, 106, 1, 1, 1,
311 & snx, sny, 1, 's', 'xy', mythid )
312
313 #endif /* ALLOW_VSTRESS_CONTROL */
314
315 c-- ===========================
316 c-- Open boundary contributions.
317 c-- ===========================
318
319 c----------------------------------------------------------------------
320 c--
321 #ifdef ALLOW_OBCSN_CONTROL
322 c-- Northern obc.
323 call ctrl_init_rec (
324 I xx_obcsnstartdate1, xx_obcsnstartdate2, xx_obcsnperiod, 4,
325 O xx_obcsnstartdate, diffrec, startrec, endrec,
326 I mythid )
327 call ctrl_init_ctrlvar (
328 & xx_obcsn_file, 11, 111, diffrec, startrec, endrec,
329 & snx, 1, nr, 'm', 'xz', mythid )
330 #endif /* ALLOW_OBCSN_CONTROL */
331
332 c----------------------------------------------------------------------
333 c--
334 #ifdef ALLOW_OBCSS_CONTROL
335 c-- Southern obc.
336 call ctrl_init_rec (
337 I xx_obcssstartdate1, xx_obcssstartdate2, xx_obcssperiod, 4,
338 O xx_obcssstartdate, diffrec, startrec, endrec,
339 I mythid )
340 call ctrl_init_ctrlvar (
341 & xx_obcss_file, 12, 112, diffrec, startrec, endrec,
342 & snx, 1, nr, 'm', 'xz', mythid )
343 #endif /* ALLOW_OBCSS_CONTROL */
344
345 c----------------------------------------------------------------------
346 c--
347 #ifdef ALLOW_OBCSW_CONTROL
348 c-- Western obc.
349 call ctrl_init_rec (
350 I xx_obcswstartdate1, xx_obcswstartdate2, xx_obcswperiod, 4,
351 O xx_obcswstartdate, diffrec, startrec, endrec,
352 I mythid )
353 call ctrl_init_ctrlvar (
354 & xx_obcsw_file, 13, 113, diffrec, startrec, endrec,
355 & 1, sny, nr, 'm', 'yz', mythid )
356 #endif /* ALLOW_OBCSW_CONTROL */
357
358 c----------------------------------------------------------------------
359 c--
360 #ifdef ALLOW_OBCSE_CONTROL
361 c-- Eastern obc.
362 call ctrl_init_rec (
363 I xx_obcsestartdate1, xx_obcsestartdate2, xx_obcseperiod, 4,
364 O xx_obcsestartdate, diffrec, startrec, endrec,
365 I mythid )
366 call ctrl_init_ctrlvar (
367 & xx_obcse_file, 14, 114, diffrec, startrec, endrec,
368 & 1, sny, nr, 'm', 'yz', mythid )
369 #endif /* ALLOW_OBCSE_CONTROL */
370
371 c----------------------------------------------------------------------
372 c--
373 #ifdef ALLOW_OBCS_CONTROL_MODES
374 cih Get matrices for reconstruction from barotropic-barclinic modes
375 CMM To use modes now hardcoded with ECCO_CPPOPTION. Would be good to have
376 c run-time option and also define filename=baro_invmodes.bin
377 CALL MDSFINDUNIT( dUnit, myThid )
378 length_of_rec = MDS_RECLEN( 64, NR*NR, myThid )
379 open(dUnit, file='baro_invmodes.bin', status='old',
380 & access='direct', recl=length_of_rec )
381 do j = 1,Nr
382 read(dUnit,rec=j) ((modesv(k,i,j), k=1,Nr), i=1,Nr)
383 end do
384 CLOSE( dUnit )
385 CMM double precision modesv is size [NR,NR,NR]
386 c dim one is z-space
387 c dim two is mode space
388 c dim three is the total depth for which this set of modes applies
389 c so for example modesv(:,2,nr) will be the second mode
390 c in z-space for the full model depth
391 c The modes are to be orthogonal when weighted by dz.
392 c i.e. if f_i(z) = mode i, sum_j(f_i(z_j)*f_j(z_j)*dz_j = delta_ij
393 c first mode should also be constant in depth...barotropic
394 c For a matlab code example how to construct the orthonormal modes,
395 c which are ideally the solution of planetary vertical mode equation
396 c using model mean dRho/dz, see
397 c MITgcm/verification/obcs_ctrl/input/gendata.m
398 c This code is compatible with partial cells
399 #endif
400
401 c----------------------------------------------------------------------
402 c--
403 #ifdef ALLOW_DIFFKR_CONTROL
404 call ctrl_init_ctrlvar (
405 & xx_diffkr_file, 15, 115, 1, 1, 1,
406 & snx, sny, nr, 'c', '3d', mythid )
407 #endif /* ALLOW_DIFFKR_CONTROL */
408
409 c----------------------------------------------------------------------
410 c--
411 #ifdef ALLOW_KAPGM_CONTROL
412 call ctrl_init_ctrlvar (
413 & xx_kapgm_file, 16, 116, 1, 1, 1,
414 & snx, sny, nr, 'c', '3d', mythid )
415 #endif /* ALLOW_KAPGM_CONTROL */
416
417 c----------------------------------------------------------------------
418 c--
419 #ifdef ALLOW_TR10_CONTROL
420 call ctrl_init_ctrlvar (
421 & xx_tr1_file, 17, 117, 1, 1, 1,
422 & snx, sny, nr, 'c', '3d', mythid )
423 #endif /* ALLOW_TR10_CONTROL */
424
425 c----------------------------------------------------------------------
426 c--
427 #if (defined (ALLOW_SST_CONTROL))
428 call ctrl_init_rec (
429 I xx_sststartdate1, xx_sststartdate2, xx_sstperiod, 1,
430 O xx_sststartdate, diffrec, startrec, endrec,
431 I mythid )
432 call ctrl_init_ctrlvar (
433 & xx_sst_file, 18, 118, diffrec, startrec, endrec,
434 & snx, sny, 1, 'c', 'xy', mythid )
435
436 #elif (defined (ALLOW_SST0_CONTROL))
437 call ctrl_init_ctrlvar (
438 & xx_sst_file, 18, 118, 1, 1, 1,
439 & snx, sny, 1, 'c', 'xy', mythid )
440
441 #endif /* ALLOW_SST_CONTROL */
442
443 c----------------------------------------------------------------------
444 c--
445 #if (defined (ALLOW_SSS_CONTROL))
446 call ctrl_init_rec (
447 I xx_sssstartdate1, xx_sssstartdate2, xx_sssperiod, 1,
448 O xx_sssstartdate, diffrec, startrec, endrec,
449 I mythid )
450 call ctrl_init_ctrlvar (
451 & xx_sss_file, 19, 119, diffrec, startrec, endrec,
452 & snx, sny, 1, 'c', 'xy', mythid )
453
454 #elif (defined (ALLOW_SSS0_CONTROL))
455 call ctrl_init_ctrlvar (
456 & xx_sss_file, 19, 119, 1, 1, 1,
457 & snx, sny, 1, 'c', 'xy', mythid )
458
459 #endif /* ALLOW_SSS0_CONTROL */
460
461 c----------------------------------------------------------------------
462 c--
463 #ifdef ALLOW_DEPTH_CONTROL
464 call ctrl_init_ctrlvar (
465 & xx_depth_file, 20, 120, 1, 1, 1,
466 & snx, sny, 1, 'c', 'xy', mythid )
467 #endif /* ALLOW_DEPTH_CONTROL */
468
469 c----------------------------------------------------------------------
470 c--
471 #ifdef ALLOW_EFLUXY0_CONTROL
472 call ctrl_init_ctrlvar (
473 & xx_efluxy_file, 21, 121, 1, 1, 1,
474 & snx, sny, nr, 's', '3d', mythid )
475 #endif /* ALLOW_EFLUXY0_CONTROL */
476
477 c----------------------------------------------------------------------
478 c--
479 #ifdef ALLOW_EFLUXP0_CONTROL
480 call ctrl_init_ctrlvar (
481 & xx_efluxp_file, 22, 122, 1, 1, 1,
482 & snx, sny, nr, 'v', '3d', mythid )
483 #endif /* ALLOW_EFLUXP0_CONTROL */
484
485 c----------------------------------------------------------------------
486 c--
487 #ifdef ALLOW_BOTTOMDRAG_CONTROL
488 call ctrl_init_ctrlvar (
489 & xx_bottomdrag_file, 23, 123, 1, 1, 1,
490 & snx, sny, 1, 'c', 'xy', mythid )
491 #endif /* ALLOW_BOTTOMDRAG_CONTROL */
492
493 c----------------------------------------------------------------------
494 c--
495 #ifdef ALLOW_HFLUXM_CONTROL
496 call ctrl_init_ctrlvar (
497 & xx_hfluxm_file, 24, 124, 1, 1, 1,
498 & snx, sny, 1, 'c', 'xy', mythid )
499 #endif /* ALLOW_HFLUXM_CONTROL */
500
501 c----------------------------------------------------------------------
502 c--
503 #ifdef ALLOW_EDDYPSI_CONTROL
504 call ctrl_init_ctrlvar (
505 & xx_edtaux_file, 25, 125, 1, 1, 1,
506 & snx, sny, nr, 'w', '3d', mythid )
507
508 call ctrl_init_ctrlvar (
509 & xx_edtauy_file, 26, 126, 1, 1, 1,
510 & snx, sny, nr, 's', '3d', mythid )
511 #endif /* ALLOW_EDDYPSI_CONTROL */
512
513 c----------------------------------------------------------------------
514 c--
515 #ifdef ALLOW_UVEL0_CONTROL
516 call ctrl_init_ctrlvar (
517 & xx_uvel_file, 27, 127, 1, 1, 1,
518 & snx, sny, nr, 'w', '3d', mythid )
519 #endif /* ALLOW_UVEL0_CONTROL */
520
521 c----------------------------------------------------------------------
522 c--
523 #ifdef ALLOW_VVEL0_CONTROL
524 call ctrl_init_ctrlvar (
525 & xx_vvel_file, 28, 128, 1, 1, 1,
526 & snx, sny, nr, 's', '3d', mythid )
527 #endif /* ALLOW_VVEL0_CONTROL */
528
529 c----------------------------------------------------------------------
530 c--
531 #ifdef ALLOW_ETAN0_CONTROL
532 call ctrl_init_ctrlvar (
533 & xx_etan_file, 29, 129, 1, 1, 1,
534 & snx, sny, 1, 'c', 'xy', mythid )
535 #endif /* ALLOW_VVEL0_CONTROL */
536
537 c----------------------------------------------------------------------
538 c--
539 #ifdef ALLOW_GEN2D_CONTROL
540 call ctrl_init_ctrlvar (
541 & xx_gen2d_file, 30, 130, 1, 1, 1,
542 & snx, sny, 1, 'c', 'xy', mythid )
543 #endif /* ALLOW_GEN2D_CONTROL */
544
545 c----------------------------------------------------------------------
546 c--
547 #ifdef ALLOW_GEN3D_CONTROL
548 call ctrl_init_ctrlvar (
549 & xx_gen3d_file, 31, 131, 1, 1, 1,
550 & snx, sny, nr, 'c', '3d', mythid )
551 #endif /* ALLOW_GEN3D_CONTROL */
552
553 c----------------------------------------------------------------------
554 c--
555 #ifdef ALLOW_PRECIP_CONTROL
556 c-- Atmos. precipitation
557 call ctrl_init_rec (
558 I xx_precipstartdate1, xx_precipstartdate2, xx_precipperiod,1,
559 O xx_precipstartdate, diffrec, startrec, endrec,
560 I mythid )
561 call ctrl_init_ctrlvar (
562 & xx_precip_file, 32, 132, diffrec, startrec, endrec,
563 & snx, sny, 1, 'c', 'xy', mythid )
564
565 #endif /* ALLOW_PRECIP_CONTROL */
566
567 c----------------------------------------------------------------------
568 c--
569 #ifdef ALLOW_SWFLUX_CONTROL
570 c-- Atmos. swflux
571 call ctrl_init_rec (
572 I xx_swfluxstartdate1, xx_swfluxstartdate2, xx_swfluxperiod, 1,
573 O xx_swfluxstartdate, diffrec, startrec, endrec,
574 I mythid )
575 call ctrl_init_ctrlvar (
576 & xx_swflux_file, 33, 133, diffrec, startrec, endrec,
577 & snx, sny, 1, 'c', 'xy', mythid )
578
579 #endif /* ALLOW_SWFLUX_CONTROL */
580
581 c----------------------------------------------------------------------
582 c--
583 #ifdef ALLOW_SWDOWN_CONTROL
584 c-- Atmos. swdown
585 call ctrl_init_rec (
586 I xx_swdownstartdate1, xx_swdownstartdate2, xx_swdownperiod, 1,
587 O xx_swdownstartdate, diffrec, startrec, endrec,
588 I mythid )
589 call ctrl_init_ctrlvar (
590 & xx_swdown_file, 34, 134, diffrec, startrec, endrec,
591 & snx, sny, 1, 'c', 'xy', mythid )
592
593 #endif /* ALLOW_SWDOWN_CONTROL */
594
595 c----------------------------------------------------------------------
596 c--
597 #ifdef ALLOW_LWFLUX_CONTROL
598 c-- Atmos. lwflux
599 call ctrl_init_rec (
600 I xx_lwfluxstartdate1, xx_lwfluxstartdate2, xx_lwfluxperiod, 1,
601 O xx_lwfluxstartdate, diffrec, startrec, endrec,
602 I mythid )
603 call ctrl_init_ctrlvar (
604 & xx_lwflux_file, 35, 135, diffrec, startrec, endrec,
605 & snx, sny, 1, 'c', 'xy', mythid )
606
607 #endif /* ALLOW_LWFLUX_CONTROL */
608
609 c----------------------------------------------------------------------
610 c--
611 #ifdef ALLOW_LWDOWN_CONTROL
612 c-- Atmos. lwdown
613 call ctrl_init_rec (
614 I xx_lwdownstartdate1, xx_lwdownstartdate2, xx_lwdownperiod, 1,
615 O xx_lwdownstartdate, diffrec, startrec, endrec,
616 I mythid )
617 call ctrl_init_ctrlvar (
618 & xx_lwdown_file, 36, 136, diffrec, startrec, endrec,
619 & snx, sny, 1, 'c', 'xy', mythid )
620
621 #endif /* ALLOW_LWDOWN_CONTROL */
622
623 c----------------------------------------------------------------------
624 c--
625 #ifdef ALLOW_EVAP_CONTROL
626 c-- Atmos. evap
627 call ctrl_init_rec (
628 I xx_evapstartdate1, xx_evapstartdate2, xx_evapperiod, 1,
629 O xx_evapstartdate, diffrec, startrec, endrec,
630 I mythid )
631 call ctrl_init_ctrlvar (
632 & xx_evap_file, 37, 137, diffrec, startrec, endrec,
633 & snx, sny, 1, 'c', 'xy', mythid )
634
635 #endif /* ALLOW_EVAP_CONTROL */
636
637 c----------------------------------------------------------------------
638 c--
639 #ifdef ALLOW_SNOWPRECIP_CONTROL
640 c-- Atmos. snowprecip
641 call ctrl_init_rec (
642 I xx_snowprecipstartdate1, xx_snowprecipstartdate2,
643 I xx_snowprecipperiod, 1,
644 O xx_snowprecipstartdate, diffrec, startrec, endrec,
645 I mythid )
646 call ctrl_init_ctrlvar (
647 & xx_snowprecip_file, 38, 138, diffrec, startrec, endrec,
648 & snx, sny, 1, 'c', 'xy', mythid )
649
650 #endif /* ALLOW_SNOWPRECIP_CONTROL */
651
652 c----------------------------------------------------------------------
653 c--
654 #ifdef ALLOW_APRESSURE_CONTROL
655 c-- Atmos. apressure
656 call ctrl_init_rec (
657 I xx_apressurestartdate1, xx_apressurestartdate2,
658 I xx_apressureperiod, 1,
659 O xx_apressurestartdate, diffrec, startrec, endrec,
660 I mythid )
661 call ctrl_init_ctrlvar (
662 & xx_apressure_file, 39, 139, diffrec, startrec, endrec,
663 & snx, sny, 1, 'c', 'xy', mythid )
664
665 #endif /* ALLOW_APRESSURE_CONTROL */
666
667 c----------------------------------------------------------------------
668 c--
669 #ifdef ALLOW_RUNOFF_CONTROL
670 c-- Atmos. runoff
671 call ctrl_init_rec (
672 I xx_runoffstartdate1, xx_runoffstartdate2, xx_runoffperiod, 1,
673 O xx_runoffstartdate, diffrec, startrec, endrec,
674 I mythid )
675 call ctrl_init_ctrlvar (
676 & xx_runoff_file, 40, 140, diffrec, startrec, endrec,
677 & snx, sny, 1, 'c', 'xy', mythid )
678 #endif /* ALLOW_RUNOFF_CONTROL */
679
680 c----------------------------------------------------------------------
681 c--
682 #ifdef ALLOW_SIAREA_CONTROL
683 C-- so far there are no xx_siareastartdate1, etc., so we need to fudge it.
684 CML call ctrl_init_rec (
685 CML I xx_siareastartdate1, xx_siareastartdate2, xx_siareaperiod, 1,
686 CML O xx_siareastartdate, diffrec, startrec, endrec,
687 CML I mythid )
688 startrec = 1
689 endrec = 1
690 diffrec = endrec - startrec + 1
691 call ctrl_init_ctrlvar (
692 & xx_siarea_file, 41, 141, diffrec, startrec, endrec,
693 & snx, sny, 1, 'c', 'xy', mythid )
694 #endif /* ALLOW_siarea_CONTROL */
695
696 c----------------------------------------------------------------------
697 c--
698 #ifdef ALLOW_SIHEFF_CONTROL
699 C-- so far there are no xx_siheffstartdate1, etc., so we need to fudge it.
700 CML call ctrl_init_rec (
701 CML I xx_siheffstartdate1, xx_siheffstartdate2, xx_siheffperiod, 1,
702 CML O xx_siheffstartdate, diffrec, startrec, endrec,
703 CML I mythid )
704 startrec = 1
705 endrec = 1
706 diffrec = endrec - startrec + 1
707 call ctrl_init_ctrlvar (
708 & xx_siheff_file, 42, 142, diffrec, startrec, endrec,
709 & snx, sny, 1, 'c', 'xy', mythid )
710 #endif /* ALLOW_siheff_CONTROL */
711
712 c----------------------------------------------------------------------
713 c--
714 #ifdef ALLOW_SIHSNOW_CONTROL
715 C-- so far there are no xx_sihsnowstartdate1, etc., so we need to fudge it.
716 CML call ctrl_init_rec (
717 CML I xx_sihsnowstartdate1, xx_sihsnowstartdate2, xx_sihsnowperiod, 1,
718 CML O xx_sihsnowstartdate, diffrec, startrec, endrec,
719 CML I mythid )
720 startrec = 1
721 endrec = 1
722 diffrec = endrec - startrec + 1
723 call ctrl_init_ctrlvar (
724 & xx_sihsnow_file, 43, 143, diffrec, startrec, endrec,
725 & snx, sny, 1, 'c', 'xy', mythid )
726 #endif /* ALLOW_sihsnow_CONTROL */
727
728
729 c----------------------------------------------------------------------
730 c--
731 #ifdef ALLOW_KAPREDI_CONTROL
732 call ctrl_init_ctrlvar (
733 & xx_kapredi_file, 44, 144, 1, 1, 1,
734 & snx, sny, nr, 'c', '3d', mythid )
735 #endif /* ALLOW_KAPREDI_CONTROL */
736
737 c----------------------------------------------------------------------
738 c----------------------------------------------------------------------
739
740 call ctrl_init_wet( mythid )
741
742 c----------------------------------------------------------------------
743 c----------------------------------------------------------------------
744
745 #ifdef ALLOW_DIC_CONTROL
746 do i = 1, dic_n_control
747 xx_dic(i) = 0. _d 0
748 enddo
749 #endif
750
751 c----------------------------------------------------------------------
752 c----------------------------------------------------------------------
753
754 do bj = jtlo,jthi
755 do bi = itlo,ithi
756 do j = jmin,jmax
757 do i = imin,imax
758 wareaunit (i,j,bi,bj) = 1.0
759 #ifndef ALLOW_ECCO
760 whflux (i,j,bi,bj) = maskC(i,j,1,bi,bj)
761 wsflux (i,j,bi,bj) = maskC(i,j,1,bi,bj)
762 wtauu (i,j,bi,bj) = maskW(i,j,1,bi,bj)
763 wtauv (i,j,bi,bj) = maskS(i,j,1,bi,bj)
764 watemp (i,j,bi,bj) = maskC(i,j,1,bi,bj)
765 waqh (i,j,bi,bj) = maskC(i,j,1,bi,bj)
766 wprecip (i,j,bi,bj) = maskC(i,j,1,bi,bj)
767 wswflux (i,j,bi,bj) = maskC(i,j,1,bi,bj)
768 wswdown (i,j,bi,bj) = maskC(i,j,1,bi,bj)
769 wuwind (i,j,bi,bj) = maskC(i,j,1,bi,bj)
770 wvwind (i,j,bi,bj) = maskC(i,j,1,bi,bj)
771 wlwflux (i,j,bi,bj) = maskC(i,j,1,bi,bj)
772 wlwdown (i,j,bi,bj) = maskC(i,j,1,bi,bj)
773 wevap (i,j,bi,bj) = maskC(i,j,1,bi,bj)
774 wsnowprecip(i,j,bi,bj) = maskC(i,j,1,bi,bj)
775 wapressure(i,j,bi,bj) = maskC(i,j,1,bi,bj)
776 wrunoff (i,j,bi,bj) = maskC(i,j,1,bi,bj)
777 wsst (i,j,bi,bj) = maskC(i,j,1,bi,bj)
778 wsss (i,j,bi,bj) = maskC(i,j,1,bi,bj)
779 #endif
780 enddo
781 enddo
782 enddo
783 enddo
784
785 return
786 end
787
788 subroutine ctrl_init_rec(
789 I fldstartdate1, fldstartdate2, fldperiod, nfac,
790 O fldstartdate, diffrec, startrec, endrec,
791 I mythid )
792
793 c ==================================================================
794 c SUBROUTINE ctrl_init_rec
795 c ==================================================================
796 c
797 c helper routine to compute the first and last record of a
798 c time dependent control variable
799 c
800 c Martin.Losch@awi.de, 2011-Mar-15
801 c
802 c ==================================================================
803 c SUBROUTINE ctrl_init_rec
804 c ==================================================================
805
806 implicit none
807
808 c == global variables ==
809 #include "SIZE.h"
810 #include "EEPARAMS.h"
811 #include "PARAMS.h"
812 #ifdef ALLOW_CAL
813 # include "cal.h"
814 #endif
815
816 c == input variables ==
817 c fldstartdate1/2 : start time (date/time) of fld
818 c fldperod : sampling interval of fld
819 c nfac : factor for the case that fld is an obcs variable
820 c in this case nfac = 4, otherwise nfac = 1
821 c mythid : thread ID of this instance
822 integer fldstartdate1
823 integer fldstartdate2
824 _RL fldperiod
825 integer nfac
826 integer mythid
827
828 c == output variables ==
829 c fldstartdate : full date from fldstartdate1 and 2
830 c startrec : first record of ctrl variable
831 c startrec : last record of ctrl variable
832 c diffrec : difference between first and last record of ctrl variable
833 integer fldstartdate(4)
834 integer startrec
835 integer endrec
836 integer diffrec
837
838 c == local variables ==
839 integer i
840 #ifdef ALLOW_CAL
841 integer difftime(4)
842 _RL diffsecs
843 #endif /* ALLOW_CAL */
844
845 c initialise some output
846 do i = 1,4
847 fldstartdate(i) = 0
848 end do
849 startrec = 0
850 endrec = 0
851 diffrec = 0
852 # ifdef ALLOW_CAL
853 call cal_FullDate( fldstartdate1, fldstartdate2,
854 & fldstartdate , mythid )
855 call cal_TimePassed( fldstartdate, modelstartdate,
856 & difftime, mythid )
857 call cal_ToSeconds ( difftime, diffsecs, mythid )
858 if ( fldperiod .EQ. -12. ) then
859 startrec = 1
860 endrec = 12*nfac
861 elseif ( fldperiod .EQ. 0. ) then
862 startrec = 1
863 endrec = 1*nfac
864 else
865 startrec = int((modelstart + startTime - diffsecs)/
866 & fldperiod) + 1
867 endrec = int((modelend + startTime - diffsecs + modelstep/2)/
868 & fldperiod) + 2
869 if ( nfac .ne. 1 ) then
870 c This is the case of obcs. I am not sure that this is correct, but
871 c it seems to work in most configurations.
872 endrec = (endrec - startrec + 1)*nfac
873 startrec = (startrec - 1)*nfac + 1
874 endif
875 endif
876 # else /* ndef ALLOW_CAL */
877 if ( fldperiod .EQ. 0. ) then
878 startrec = 1
879 endrec = 1*nfac
880 else
881 startrec = 1
882 endrec = (int((endTime - startTime)/fldperiod) + 1)*nfac
883 endif
884 #endif /* ALLOW_CAL */
885 diffrec = endrec - startrec + 1
886
887 return
888 end

  ViewVC Help
Powered by ViewVC 1.1.22