/[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.46 - (show annotations) (download)
Fri Jul 27 18:22:42 2012 UTC (11 years, 10 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint63q
Changes since 1.45: +24 -1 lines
Enable generic arrays of control variables

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

  ViewVC Help
Powered by ViewVC 1.1.22