/[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.44 - (show annotations) (download)
Wed Apr 11 12:11:30 2012 UTC (12 years, 1 month ago) by mlosch
Branch: MAIN
Changes since 1.43: +62 -40 lines
 - add a formal parameter to s/r ctrl_init_rec to print
   xx_$(ctrl_valiable)_file for easier debugging
 - move call of cal_* routines into the if-block where they are really
   needed so that xx_$(ctrl_variable)_period can be 0 or -12 without having
   to specify a startdate

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

  ViewVC Help
Powered by ViewVC 1.1.22