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