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