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