/[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.47 - (show annotations) (download)
Tue Jul 31 16:05:57 2012 UTC (11 years, 10 months ago) by heimbach
Branch: MAIN
Changes since 1.46: +3 -1 lines
Attempt at adding CTRL_SIZE.h

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

  ViewVC Help
Powered by ViewVC 1.1.22