/[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.40 - (show annotations) (download)
Tue May 24 20:56:54 2011 UTC (13 years ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62z, checkpoint62y, checkpoint63
Changes since 1.39: +96 -100 lines
remove #include "OBCS.h" (don't seem to be needed)

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

  ViewVC Help
Powered by ViewVC 1.1.22