/[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.54 - (show annotations) (download)
Tue Feb 5 22:45:56 2013 UTC (11 years, 3 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65d, checkpoint65, checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f
Changes since 1.53: +3 -2 lines
- ctrl_init.F : fix previous modification.
- ctrl_readparms.F : read/init mult_gentim2d.

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

  ViewVC Help
Powered by ViewVC 1.1.22