/[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.50 - (show annotations) (download)
Tue Aug 28 19:18:26 2012 UTC (11 years, 9 months ago) by gforget
Branch: MAIN
Changes since 1.49: +21 -19 lines
- pkg/exf : added run time switch useAtmWind to replace ALLOW_ATM_WIND
  cpp switch. ALLOW_ATM_WIND now just sets the useAtmWind default (see
  exf_readparms.F) and force defines ALLOW_BULKFORMULAE (EXF_OPTIONS.h).
- pkg/exf, autodiff, ctrl, ecco and seaice : remove ALLOW_ATM_WIND
  brackets, or replace them with useAtmWind ones.
- pkg/ctrl, ecco : allow to compile both ALLOW_U/VSTRESS_CONTROL and
  ALLOW_U/VWIND_CONTROL. Depending on useAtmWind, one is inactive,
  and the other is active (see exf_getffields.F/exf_getsurfacefluxes.F).

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

  ViewVC Help
Powered by ViewVC 1.1.22