/[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.60 - (show annotations) (download)
Sun Nov 9 13:10:03 2014 UTC (9 years, 6 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint65g
Changes since 1.59: +16 -1 lines
- carry out smooth_correl2D at initialization
  stage rather than during time stepping.

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

  ViewVC Help
Powered by ViewVC 1.1.22