/[MITgcm]/MITgcm_contrib/dcarroll/highres_darwin/code/exf_getffields.F
ViewVC logotype

Contents of /MITgcm_contrib/dcarroll/highres_darwin/code/exf_getffields.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.1 - (show annotations) (download)
Sun Sep 22 21:23:46 2019 UTC (5 years, 10 months ago) by dcarroll
Branch: MAIN
CVS Tags: HEAD
Initial check in of high resolution Darwin simulation code

1 C $Header: /u/gcmpack/MITgcm_contrib/ecco_darwin/v4_llc270/code_darwin/exf_getffields.F,v 1.1 2017/12/01 19:02:11 dimitri Exp $
2 C $Name: $
3
4 #include "EXF_OPTIONS.h"
5 #ifdef ALLOW_CTRL
6 # include "CTRL_OPTIONS.h"
7 #endif
8 #ifdef ALLOW_ECCO
9 # include "ECCO_OPTIONS.h"
10 #endif
11
12 SUBROUTINE EXF_GETFFIELDS( myTime, myIter, myThid )
13
14 C ==================================================================
15 C SUBROUTINE exf_getffields
16 C ==================================================================
17 C
18 C o Read-in atmospheric state and/or surface fluxes from files.
19 C
20 C heimbach@mit.edu, 23-May-2003 totally re-structured
21 C 5-Aug-2003: added USE_EXF_INTERPOLATION for arbitrary input grid
22 C
23 C ==================================================================
24 C SUBROUTINE exf_getffields
25 C ==================================================================
26
27 IMPLICIT NONE
28
29 C == global variables ==
30
31 #include "EEPARAMS.h"
32 #include "SIZE.h"
33 #include "PARAMS.h"
34 #include "DYNVARS.h"
35 #include "GRID.h"
36
37 #include "EXF_PARAM.h"
38 #include "EXF_FIELDS.h"
39 #include "EXF_CONSTANTS.h"
40
41 #ifdef ALLOW_CTRL
42 # include "CTRL_SIZE.h"
43 # include "ctrl.h"
44 # include "ctrl_dummy.h"
45 # ifdef ALLOW_GENTIM2D_CONTROL
46 # include "CTRL_GENARR.h"
47 # endif
48 #endif
49 #if (defined (ALLOW_ECCO) && defined (ECCO_CTRL_DEPRECATED))
50 # include "ecco_cost.h"
51 #endif
52
53 C == routine arguments ==
54 _RL myTime
55 INTEGER myIter
56 INTEGER myThid
57
58 C == local variables ==
59 INTEGER i, j, bi, bj
60 #ifdef ALLOW_ROTATE_UV_CONTROLS
61 _RL tmpUE(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
62 _RL tmpVN(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
63 _RL tmpUX(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
64 _RL tmpVY(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
65 #endif
66 #if (defined (ALLOW_CTRL) && \
67 defined (ALLOW_GENTIM2D_CONTROL))
68 INTEGER iarr
69 #endif
70
71 C == end of interface ==
72
73 c Atmospheric carbon dioxide concentration
74 CALL EXF_SET_FLD(
75 I 'apco2', apco2file, apco2mask,
76 I apco2StartTime, apco2period, apco2RepCycle,
77 I exf_inscal_apco2,
78 I apco2_exfremo_intercept, apco2_exfremo_slope,
79 U apco2, apco20, apco21,
80 #ifdef USE_EXF_INTERPOLATION
81 I apco2_lon0, apco2_lon_inc, apco2_lat0, apco2_lat_inc,
82 I apco2_nlon, apco2_nlat, xC, yC, apco2_interpMethod,
83 #endif
84 I myTime, myIter, myThid )
85
86 C-- read forcing fields from files and temporal interpolation
87
88 C- Zonal and meridional wind stress.
89 IF ( .NOT.useAtmWind ) THEN
90 CALL EXF_SET_UV(
91 I 'ustress', ustressfile, ustressmask,
92 I ustressStartTime, ustressperiod, ustressRepCycle,
93 I exf_inscal_ustress,
94 I ustress_exfremo_intercept, ustress_exfremo_slope,
95 U ustress, ustress0, ustress1,
96 I 'vstress', vstressfile, vstressmask,
97 I vstressStartTime, vstressperiod, vstressRepCycle,
98 I exf_inscal_vstress,
99 I vstress_exfremo_intercept, vstress_exfremo_slope,
100 U vstress, vstress0, vstress1,
101 #ifdef USE_EXF_INTERPOLATION
102 I ustress_lon0, ustress_lon_inc, ustress_lat0, ustress_lat_inc,
103 I ustress_nlon, ustress_nlat, ustress_interpMethod,
104 I vstress_lon0, vstress_lon_inc, vstress_lat0, vstress_lat_inc,
105 I vstress_nlon, vstress_nlat, vstress_interpMethod,
106 I uvInterp_stress,
107 #endif /* USE_EXF_INTERPOLATION */
108 I myTime, myIter, myThid )
109 ELSE
110 DO bj = myByLo(myThid),myByHi(myThid)
111 DO bi = myBxLo(myThid),mybxhi(myThid)
112 DO j=1-OLy,sNy+OLy
113 DO i=1-OLx,sNx+OLx
114 ustress(i,j,bi,bj) = 0. _d 0
115 vstress(i,j,bi,bj) = 0. _d 0
116 ENDDO
117 ENDDO
118 ENDDO
119 ENDDO
120 ENDIF
121
122 C- Wind speed
123 CALL EXF_SET_FLD(
124 I 'wspeed', wspeedfile, wspeedmask,
125 I wspeedStartTime, wspeedperiod, wspeedRepCycle,
126 I exf_inscal_wspeed,
127 I wspeed_exfremo_intercept, wspeed_exfremo_slope,
128 U wspeed, wspeed0, wspeed1,
129 #ifdef USE_EXF_INTERPOLATION
130 I wspeed_lon0, wspeed_lon_inc,
131 I wspeed_lat0, wspeed_lat_inc,
132 I wspeed_nlon, wspeed_nlat, xC, yC, wspeed_interpMethod,
133 #endif
134 I myTime, myIter, myThid )
135
136 C- Zonal and meridional wind.
137 IF ( useAtmWind ) THEN
138 CALL EXF_SET_UV(
139 I 'uwind', uwindfile, uwindmask,
140 I uwindStartTime, uwindperiod, uwindRepCycle,
141 I exf_inscal_uwind,
142 I uwind_exfremo_intercept, uwind_exfremo_slope,
143 U uwind, uwind0, uwind1,
144 I 'vwind', vwindfile, vwindmask,
145 I vwindStartTime, vwindperiod, vwindRepCycle,
146 I exf_inscal_vwind,
147 I vwind_exfremo_intercept, vwind_exfremo_slope,
148 U vwind, vwind0, vwind1,
149 #ifdef USE_EXF_INTERPOLATION
150 I uwind_lon0, uwind_lon_inc, uwind_lat0, uwind_lat_inc,
151 I uwind_nlon, uwind_nlat, uwind_interpMethod,
152 I vwind_lon0, vwind_lon_inc, vwind_lat0, vwind_lat_inc,
153 I vwind_nlon, vwind_nlat, vwind_interpMethod, uvInterp_wind,
154 #endif /* USE_EXF_INTERPOLATION */
155 I myTime, myIter, myThid )
156
157 IF (useRelativeWind) THEN
158 C Subtract UVEL and VVEL from UWIND and VWIND.
159 DO bj = myByLo(myThid),myByHi(myThid)
160 DO bi = myBxLo(myThid),mybxhi(myThid)
161 DO j = 1,sNy
162 DO i = 1,sNx
163 uwind(i,j,bi,bj) = uwind(i,j,bi,bj) - 0.5 _d 0
164 & * (uVel(i,j,1,bi,bj)+uVel(i+1,j,1,bi,bj))
165 vwind(i,j,bi,bj) = vwind(i,j,bi,bj) - 0.5 _d 0
166 & * (vVel(i,j,1,bi,bj)+vVel(i,j+1,1,bi,bj))
167 ENDDO
168 ENDDO
169 ENDDO
170 ENDDO
171 ENDIF
172
173 ELSE
174 DO bj = myByLo(myThid),myByHi(myThid)
175 DO bi = myBxLo(myThid),mybxhi(myThid)
176 DO j=1-OLy,sNy+OLy
177 DO i=1-OLx,sNx+OLx
178 uwind(i,j,bi,bj) = 0. _d 0
179 vwind(i,j,bi,bj) = 0. _d 0
180 ENDDO
181 ENDDO
182 ENDDO
183 ENDDO
184 ENDIF
185
186 C- Atmospheric heat flux.
187 CALL EXF_SET_FLD(
188 I 'hflux', hfluxfile, hfluxmask,
189 I hfluxStartTime, hfluxperiod, hfluxRepCycle,
190 I exf_inscal_hflux,
191 I hflux_exfremo_intercept, hflux_exfremo_slope,
192 U hflux, hflux0, hflux1,
193 #ifdef USE_EXF_INTERPOLATION
194 I hflux_lon0, hflux_lon_inc, hflux_lat0, hflux_lat_inc,
195 I hflux_nlon, hflux_nlat, xC, yC, hflux_interpMethod,
196 #endif
197 I myTime, myIter, myThid )
198
199 C- Freshwater flux.
200 CALL EXF_SET_FLD(
201 I 'sflux', sfluxfile, sfluxmask,
202 I sfluxStartTime, sfluxperiod, sfluxRepCycle,
203 I exf_inscal_sflux,
204 I sflux_exfremo_intercept, sflux_exfremo_slope,
205 U sflux, sflux0, sflux1,
206 #ifdef USE_EXF_INTERPOLATION
207 I sflux_lon0, sflux_lon_inc, sflux_lat0, sflux_lat_inc,
208 I sflux_nlon, sflux_nlat, xC, yC, sflux_interpMethod,
209 #endif
210 I myTime, myIter, myThid )
211
212 #ifdef ALLOW_ATM_TEMP
213
214 C- Atmospheric temperature.
215 CALL EXF_SET_FLD(
216 I 'atemp', atempfile, atempmask,
217 I atempStartTime, atempperiod, atempRepCycle,
218 I exf_inscal_atemp,
219 I atemp_exfremo_intercept, atemp_exfremo_slope,
220 U atemp, atemp0, atemp1,
221 #ifdef USE_EXF_INTERPOLATION
222 I atemp_lon0, atemp_lon_inc, atemp_lat0, atemp_lat_inc,
223 I atemp_nlon, atemp_nlat, xC, yC, atemp_interpMethod,
224 #endif
225 I myTime, myIter, myThid )
226 DO bj = myByLo(myThid),myByHi(myThid)
227 DO bi = myBxLo(myThid),mybxhi(myThid)
228 DO j = 1,sNy
229 DO i = 1,sNx
230 atemp(i,j,bi,bj) = atemp(i,j,bi,bj) + exf_offset_atemp
231 ENDDO
232 ENDDO
233 ENDDO
234 ENDDO
235
236 C- Atmospheric humidity.
237 CALL EXF_SET_FLD(
238 I 'aqh', aqhfile, aqhmask,
239 I aqhStartTime, aqhperiod, aqhRepCycle,
240 I exf_inscal_aqh,
241 I aqh_exfremo_intercept, aqh_exfremo_slope,
242 U aqh, aqh0, aqh1,
243 #ifdef USE_EXF_INTERPOLATION
244 I aqh_lon0, aqh_lon_inc, aqh_lat0, aqh_lat_inc,
245 I aqh_nlon, aqh_nlat, xC, yC, aqh_interpMethod,
246 #endif
247 I myTime, myIter, myThid )
248
249 # ifdef ALLOW_READ_TURBFLUXES
250
251 C- Sensible Heat flux
252 CALL EXF_SET_FLD(
253 I 'hs', hs_file, hs_mask,
254 I hs_StartTime, hs_period, hs_RepCycle,
255 I exf_inscal_hs,
256 I hs_exfremo_intercept, hs_exfremo_slope,
257 U hs, hs0, hs1,
258 # ifdef USE_EXF_INTERPOLATION
259 I hs_lon0, hs_lon_inc, hs_lat0, hs_lat_inc,
260 I hs_nlon, hs_nlat, xC, yC, hs_interpMethod,
261 # endif
262 I myTime, myIter, myThid )
263
264 C- Latent Heat flux
265 CALL EXF_SET_FLD(
266 I 'hl', hl_file, hl_mask,
267 I hl_StartTime, hl_period, hl_RepCycle,
268 I exf_inscal_hl,
269 I hl_exfremo_intercept, hl_exfremo_slope,
270 U hl, hl0, hl1,
271 # ifdef USE_EXF_INTERPOLATION
272 I hl_lon0, hl_lon_inc, hl_lat0, hl_lat_inc,
273 I hl_nlon, hl_nlat, xC, yC, hl_interpMethod,
274 # endif
275 I myTime, myIter, myThid )
276
277 # endif /* ALLOW_READ_TURBFLUXES */
278
279 C- Net long wave radiative flux.
280 CALL EXF_SET_FLD(
281 I 'lwflux', lwfluxfile, lwfluxmask,
282 I lwfluxStartTime, lwfluxperiod, lwfluxRepCycle,
283 I exf_inscal_lwflux,
284 I lwflux_exfremo_intercept, lwflux_exfremo_slope,
285 U lwflux, lwflux0, lwflux1,
286 #ifdef USE_EXF_INTERPOLATION
287 I lwflux_lon0, lwflux_lon_inc, lwflux_lat0, lwflux_lat_inc,
288 I lwflux_nlon, lwflux_nlat, xC, yC, lwflux_interpMethod,
289 #endif
290 I myTime, myIter, myThid )
291
292 #ifdef EXF_READ_EVAP
293 C- Evaporation
294 CALL EXF_SET_FLD(
295 I 'evap', evapfile, evapmask,
296 I evapStartTime, evapperiod, evapRepCycle,
297 I exf_inscal_evap,
298 I evap_exfremo_intercept, evap_exfremo_slope,
299 U evap, evap0, evap1,
300 #ifdef USE_EXF_INTERPOLATION
301 I evap_lon0, evap_lon_inc, evap_lat0, evap_lat_inc,
302 I evap_nlon, evap_nlat, xC, yC, evap_interpMethod,
303 #endif
304 I myTime, myIter, myThid )
305 #endif /* EXF_READ_EVAP */
306
307 C- Precipitation.
308 CALL EXF_SET_FLD(
309 I 'precip', precipfile, precipmask,
310 I precipStartTime, precipperiod, precipRepCycle,
311 I exf_inscal_precip,
312 I precip_exfremo_intercept, precip_exfremo_slope,
313 U precip, precip0, precip1,
314 #ifdef USE_EXF_INTERPOLATION
315 I precip_lon0, precip_lon_inc, precip_lat0, precip_lat_inc,
316 I precip_nlon, precip_nlat, xC, yC, precip_interpMethod,
317 #endif
318 I myTime, myIter, myThid )
319
320 C- Snow.
321 CALL EXF_SET_FLD(
322 I 'snowprecip', snowprecipfile, snowprecipmask,
323 I snowprecipStartTime, snowprecipperiod, snowprecipRepCycle,
324 I exf_inscal_snowprecip,
325 I snowprecip_exfremo_intercept, snowprecip_exfremo_slope,
326 U snowprecip, snowprecip0, snowprecip1,
327 #ifdef USE_EXF_INTERPOLATION
328 I snowprecip_lon0, snowprecip_lon_inc,
329 I snowprecip_lat0, snowprecip_lat_inc,
330 I snowprecip_nlon, snowprecip_nlat, xC, yC,
331 I snowprecip_interpMethod,
332 #endif
333 I myTime, myIter, myThid )
334 C Take care of case where total precip is not defined
335 IF ( snowPrecipFile .NE. ' ' ) THEN
336 DO bj = myByLo(myThid),myByHi(myThid)
337 DO bi = myBxLo(myThid),mybxhi(myThid)
338 DO j = 1,sNy
339 DO i = 1,sNx
340 precip(i,j,bi,bj) =
341 & MAX( precip(i,j,bi,bj), snowPrecip(i,j,bi,bj) )
342 ENDDO
343 ENDDO
344 ENDDO
345 ENDDO
346 ENDIF
347
348 #endif /* ALLOW_ATM_TEMP */
349
350 #if defined(ALLOW_ATM_TEMP) || defined(SHORTWAVE_HEATING)
351 C- Net short wave radiative flux.
352 CALL EXF_SET_FLD(
353 I 'swflux', swfluxfile, swfluxmask,
354 I swfluxStartTime, swfluxperiod, swfluxRepCycle,
355 I exf_inscal_swflux,
356 I swflux_exfremo_intercept, swflux_exfremo_slope,
357 U swflux, swflux0, swflux1,
358 #ifdef USE_EXF_INTERPOLATION
359 I swflux_lon0, swflux_lon_inc, swflux_lat0, swflux_lat_inc,
360 I swflux_nlon, swflux_nlat, xC, yC, swflux_interpMethod,
361 #endif
362 I myTime, myIter, myThid )
363 #endif
364
365 #ifdef ALLOW_DOWNWARD_RADIATION
366
367 C- Downward shortwave radiation.
368 CALL EXF_SET_FLD(
369 I 'swdown', swdownfile, swdownmask,
370 I swdownStartTime, swdownperiod, swdownRepCycle,
371 I exf_inscal_swdown,
372 I swdown_exfremo_intercept, swdown_exfremo_slope,
373 U swdown, swdown0, swdown1,
374 #ifdef USE_EXF_INTERPOLATION
375 I swdown_lon0, swdown_lon_inc, swdown_lat0, swdown_lat_inc,
376 I swdown_nlon, swdown_nlat, xC, yC, swdown_interpMethod,
377 #endif
378 I myTime, myIter, myThid )
379
380 C- Downward longwave radiation.
381 CALL EXF_SET_FLD(
382 I 'lwdown', lwdownfile, lwdownmask,
383 I lwdownStartTime, lwdownperiod, lwdownRepCycle,
384 I exf_inscal_lwdown,
385 I lwdown_exfremo_intercept, lwdown_exfremo_slope,
386 U lwdown, lwdown0, lwdown1,
387 #ifdef USE_EXF_INTERPOLATION
388 I lwdown_lon0, lwdown_lon_inc, lwdown_lat0, lwdown_lat_inc,
389 I lwdown_nlon, lwdown_nlat, xC, yC, lwdown_interpMethod,
390 #endif
391 I myTime, myIter, myThid )
392
393 #endif /* ALLOW_DOWNWARD_RADIATION */
394
395 #ifdef ATMOSPHERIC_LOADING
396 C- Atmos. pressure forcing
397 CALL EXF_SET_FLD(
398 I 'apressure', apressurefile, apressuremask,
399 I apressureStartTime, apressureperiod, apressureRepCycle,
400 I exf_inscal_apressure,
401 I apressure_exfremo_intercept, apressure_exfremo_slope,
402 U apressure, apressure0, apressure1,
403 #ifdef USE_EXF_INTERPOLATION
404 I apressure_lon0, apressure_lon_inc,
405 I apressure_lat0, apressure_lat_inc,
406 I apressure_nlon,apressure_nlat,xC,yC, apressure_interpMethod,
407 #endif
408 I myTime, myIter, myThid )
409 #endif
410
411 #ifdef EXF_ALLOW_TIDES
412 C- Tidal geopotential
413 CALL EXF_SET_FLD(
414 I 'tidePot', tidePotFile, tidePotMask,
415 I tidePotStartTime, tidePotPeriod, tidePotRepCycle,
416 I exf_inscal_tidePot,
417 I tidePot_exfremo_intercept, tidePot_exfremo_slope,
418 U tidePot, tidePot0, tidePot1,
419 #ifdef USE_EXF_INTERPOLATION
420 I tidePot_lon0, tidePot_lon_inc,
421 I tidePot_lat0, tidePot_lat_inc,
422 I tidePot_nlon, tidePot_nlat, xC, yC, tidePot_interpMethod,
423 #endif
424 I myTime, myIter, myThid )
425 #endif /* EXF_ALLOW_TIDES */
426
427 #ifdef EXF_SEAICE_FRACTION
428 C- fractional ice-covered area mask
429 CALL EXF_SET_FLD(
430 I 'areamask', areamaskfile, areamaskmask,
431 I areamaskStartTime, areamaskperiod, areamaskRepCycle,
432 I exf_inscal_areamask,
433 I areamask_exfremo_intercept, areamask_exfremo_slope,
434 U areamask, areamask0, areamask1,
435 #ifdef USE_EXF_INTERPOLATION
436 I areamask_lon0, areamask_lon_inc,
437 I areamask_lat0, areamask_lat_inc,
438 I areamask_nlon, areamask_nlat, xC, yC, areamask_interpMethod,
439 #endif
440 I myTime, myIter, myThid )
441 #endif
442
443 #ifdef ALLOW_RUNOFF
444 C- Runoff
445 CALL EXF_SET_FLD(
446 I 'runoff', runofffile, runoffmask,
447 I runoffStartTime, runoffperiod, runoffRepCycle,
448 I exf_inscal_runoff,
449 I runoff_exfremo_intercept, runoff_exfremo_slope,
450 U runoff, runoff0, runoff1,
451 #ifdef USE_EXF_INTERPOLATION
452 I runoff_lon0, runoff_lon_inc, runoff_lat0, runoff_lat_inc,
453 I runoff_nlon, runoff_nlat, xC, yC, runoff_interpMethod,
454 #endif
455 I myTime, myIter, myThid )
456 #endif /* ALLOW_RUNOFF */
457
458 #ifdef ALLOW_RUNOFTEMP
459 C- Runoff temperature
460 CALL EXF_SET_FLD(
461 I 'runoftemp', runoftempfile, runoffmask,
462 I runoffStartTime, runoffperiod, runoffRepCycle,
463 I exf_inscal_runoftemp,
464 I runoftemp_exfremo_intercept, runoftemp_exfremo_slope,
465 U runoftemp, runoftemp0, runoftemp1,
466 #ifdef USE_EXF_INTERPOLATION
467 I runoff_lon0, runoff_lon_inc, runoff_lat0, runoff_lat_inc,
468 I runoff_nlon, runoff_nlat, xC, yC, runoff_interpMethod,
469 #endif
470 I myTime, myIter, myThid )
471 #endif /* ALLOW_RUNOFTEMP */
472
473 #ifdef ALLOW_SALTFLX
474 C- Salt flux
475 CALL EXF_SET_FLD(
476 I 'saltflx', saltflxfile, saltflxmask,
477 I saltflxStartTime, saltflxperiod, saltflxRepCycle,
478 I exf_inscal_saltflx,
479 I saltflx_exfremo_intercept, saltflx_exfremo_slope,
480 U saltflx, saltflx0, saltflx1,
481 #ifdef USE_EXF_INTERPOLATION
482 I saltflx_lon0, saltflx_lon_inc,
483 I saltflx_lat0, saltflx_lat_inc,
484 I saltflx_nlon, saltflx_nlat, xC, yC, saltflx_interpMethod,
485 #endif
486 I myTime, myIter, myThid )
487 #endif
488
489 #ifdef ALLOW_ROTATE_UV_CONTROLS
490 IF ( useCTRL ) THEN
491 DO bj = myByLo(myThid),myByHi(myThid)
492 DO bi = myBxLo(myThid),mybxhi(myThid)
493 DO j = 1-OLy,sNy+OLy
494 DO i = 1-OLx,sNx+OLx
495 tmpUE(i,j,bi,bj) = 0. _d 0
496 tmpVN(i,j,bi,bj) = 0. _d 0
497 tmpUX(i,j,bi,bj) = 0. _d 0
498 tmpVY(i,j,bi,bj) = 0. _d 0
499 ENDDO
500 ENDDO
501 ENDDO
502 ENDDO
503 ENDIF
504 #endif
505
506 # if (!defined (ALLOW_ECCO) || defined (ECCO_CTRL_DEPRECATED))
507
508 C-- Control variables for atmos. state
509 #ifdef ALLOW_CTRL
510 IF (.NOT.ctrlUseGen) THEN
511
512 #ifdef ALLOW_ATEMP_CONTROL
513 CALL CTRL_GET_GEN (
514 & xx_atemp_file, xx_atempstartdate, xx_atempperiod,
515 & maskc, atemp, xx_atemp0, xx_atemp1, xx_atemp_dummy,
516 & xx_atemp_remo_intercept, xx_atemp_remo_slope,
517 & watemp, myTime, myIter, myThid )
518 #endif
519
520 #ifdef ALLOW_AQH_CONTROL
521 CALL CTRL_GET_GEN (
522 & xx_aqh_file, xx_aqhstartdate, xx_aqhperiod,
523 & maskc, aqh, xx_aqh0, xx_aqh1, xx_aqh_dummy,
524 & xx_aqh_remo_intercept, xx_aqh_remo_slope,
525 & waqh, myTime, myIter, myThid )
526 #endif
527
528 #ifdef ALLOW_PRECIP_CONTROL
529 CALL CTRL_GET_GEN (
530 & xx_precip_file, xx_precipstartdate, xx_precipperiod,
531 & maskc, precip, xx_precip0, xx_precip1, xx_precip_dummy,
532 & xx_precip_remo_intercept, xx_precip_remo_slope,
533 & wprecip, myTime, myIter, myThid )
534 #endif
535
536 #ifdef ALLOW_SWDOWN_CONTROL
537 CALL CTRL_GET_GEN (
538 & xx_swdown_file, xx_swdownstartdate, xx_swdownperiod,
539 & maskc, swdown, xx_swdown0, xx_swdown1, xx_swdown_dummy,
540 & xx_swdown_remo_intercept, xx_swdown_remo_slope,
541 & wswdown, myTime, myIter, myThid )
542 #endif
543
544 #ifdef ALLOW_LWDOWN_CONTROL
545 CALL CTRL_GET_GEN (
546 & xx_lwdown_file, xx_lwdownstartdate, xx_lwdownperiod,
547 & maskc, lwdown, xx_lwdown0, xx_lwdown1, xx_lwdown_dummy,
548 & xx_lwdown_remo_intercept, xx_lwdown_remo_slope,
549 & wlwdown, myTime, myIter, myThid )
550 #endif
551
552 ENDIF !if (.NOT.ctrlUseGen) then
553
554 #ifdef ALLOW_SWFLUX_CONTROL
555 CALL CTRL_GET_GEN (
556 & xx_swflux_file, xx_swfluxstartdate, xx_swfluxperiod,
557 & maskc, swflux, xx_swflux0, xx_swflux1, xx_swflux_dummy,
558 & xx_swflux_remo_intercept, xx_swflux_remo_slope,
559 & wswflux, myTime, myIter, myThid )
560 #endif
561
562 #ifdef ALLOW_LWFLUX_CONTROL
563 CALL CTRL_GET_GEN (
564 & xx_lwflux_file, xx_lwfluxstartdate, xx_lwfluxperiod,
565 & maskc, lwflux, xx_lwflux0, xx_lwflux1, xx_lwflux_dummy,
566 & xx_lwflux_remo_intercept, xx_lwflux_remo_slope,
567 & wswflux, myTime, myIter, myThid )
568 #endif
569
570 #ifdef ALLOW_EVAP_CONTROL
571 CALL CTRL_GET_GEN (
572 & xx_evap_file, xx_evapstartdate, xx_evapperiod,
573 & maskc, evap, xx_evap0, xx_evap1, xx_evap_dummy,
574 & xx_evap_remo_intercept, xx_evap_remo_slope,
575 & wevap, myTime, myIter, myThid )
576 #endif
577
578 #ifdef ALLOW_SNOWPRECIP_CONTROL
579 CALL CTRL_GET_GEN (
580 & xx_snowprecip_file, xx_snowprecipstartdate,
581 & xx_snowprecipperiod,
582 & maskc, snowprecip, xx_snowprecip0, xx_snowprecip1,
583 & xx_snowprecip_dummy,
584 & xx_snowprecip_remo_intercept, xx_snowprecip_remo_slope,
585 & wsnowprecip, myTime, myIter, myThid )
586 #endif
587
588 #ifdef ALLOW_APRESSURE_CONTROL
589 CALL CTRL_GET_GEN (
590 & xx_apressure_file, xx_apressurestartdate,
591 & xx_apressureperiod,
592 & maskc, apressure, xx_apressure0, xx_apressure1,
593 & xx_apressure_dummy,
594 & xx_apressure_remo_intercept, xx_apressure_remo_slope,
595 & wapressure, myTime, myIter, myThid )
596 #endif
597
598 IF ( useAtmWind ) THEN
599 #ifndef ALLOW_ROTATE_UV_CONTROLS
600
601 #ifdef ALLOW_UWIND_CONTROL
602 CALL CTRL_GET_GEN (
603 & xx_uwind_file, xx_uwindstartdate, xx_uwindperiod,
604 & maskc, uwind, xx_uwind0, xx_uwind1, xx_uwind_dummy,
605 & xx_uwind_remo_intercept, xx_uwind_remo_slope,
606 & wuwind, myTime, myIter, myThid )
607 #endif /* ALLOW_UWIND_CONTROL */
608
609 #ifdef ALLOW_VWIND_CONTROL
610 CALL CTRL_GET_GEN (
611 & xx_vwind_file, xx_vwindstartdate, xx_vwindperiod,
612 & maskc, vwind, xx_vwind0, xx_vwind1, xx_vwind_dummy,
613 & xx_vwind_remo_intercept, xx_vwind_remo_slope,
614 & wvwind, myTime, myIter, myThid )
615 #endif /* ALLOW_VWIND_CONTROL */
616
617 #else
618
619 #if defined(ALLOW_UWIND_CONTROL) && defined(ALLOW_VWIND_CONTROL)
620 CALL CTRL_GET_GEN (
621 & xx_uwind_file, xx_uwindstartdate, xx_uwindperiod,
622 & maskc, tmpUE, xx_uwind0, xx_uwind1, xx_uwind_dummy,
623 & xx_uwind_remo_intercept, xx_uwind_remo_slope,
624 & wuwind, myTime, myIter, myThid )
625
626 CALL CTRL_GET_GEN (
627 & xx_vwind_file, xx_vwindstartdate, xx_vwindperiod,
628 & maskc, tmpVN, xx_vwind0, xx_vwind1, xx_vwind_dummy,
629 & xx_vwind_remo_intercept, xx_vwind_remo_slope,
630 & wvwind, myTime, myIter, myThid )
631
632 CALL ROTATE_UV2EN_RL(tmpUX,tmpVY,tmpUE,tmpVN,
633 & .FALSE.,.FALSE.,.TRUE.,1,myThid)
634
635 DO bj = myByLo(myThid),myByHi(myThid)
636 DO bi = myBxLo(myThid),mybxhi(myThid)
637 DO j = 1,sNy
638 DO i = 1,sNx
639 uwind(i,j,bi,bj)=uwind(i,j,bi,bj)+tmpUX(i,j,bi,bj)
640 vwind(i,j,bi,bj)=vwind(i,j,bi,bj)+tmpVY(i,j,bi,bj)
641 ENDDO
642 ENDDO
643 ENDDO
644 ENDDO
645 #endif
646
647 #endif /* ALLOW_ROTATE_UV_CONTROLS */
648 ENDIF
649
650 #ifdef ALLOW_ATM_MEAN_CONTROL
651 DO bj = myByLo(myThid),myByHi(myThid)
652 DO bi = myBxLo(myThid),mybxhi(myThid)
653 DO j = 1,sNy
654 DO i = 1,sNx
655 # ifdef ALLOW_ATEMP_CONTROL
656 atemp(i,j,bi,bj) =atemp(i,j,bi,bj) +xx_atemp_mean(i,j,bi,bj)
657 # endif
658 # ifdef ALLOW_AQH_CONTROL
659 aqh(i,j,bi,bj) =aqh(i,j,bi,bj) +xx_aqh_mean(i,j,bi,bj)
660 # endif
661 # ifdef ALLOW_PRECIP_CONTROL
662 precip(i,j,bi,bj)=precip(i,j,bi,bj)+xx_precip_mean(i,j,bi,bj)
663 # endif
664 # ifdef ALLOW_SWDOWN_CONTROL
665 swdown(i,j,bi,bj)=swdown(i,j,bi,bj)+xx_swdown_mean(i,j,bi,bj)
666 # endif
667 # ifdef ALLOW_UWIND_CONTROL
668 uwind(i,j,bi,bj) =uwind(i,j,bi,bj) +xx_uwind_mean(i,j,bi,bj)
669 # endif
670 # ifdef ALLOW_VWIND_CONTROL
671 vwind(i,j,bi,bj) =vwind(i,j,bi,bj) +xx_vwind_mean(i,j,bi,bj)
672 # endif
673 ENDDO
674 ENDDO
675 ENDDO
676 ENDDO
677 #endif /* ALLOW_ATM_MEAN_CONTROL */
678
679 cdm transferred from exf_init_runoff.F
680 cdm functionality needs to be checked before turning on
681 cdm #ifdef ALLOW_RUNOFF_CONTROL
682 cdm CALL CTRL_GET_GEN (
683 cdm & xx_runoff_file, xx_runoffstartdate, xx_runoffperiod,
684 cdm & maskc, runoff, xx_runoff0, xx_runoff1, xx_runoff_dummy,
685 cdm & xx_runoff_remo_intercept, xx_runoff_remo_slope,
686 cdm & wrunoff, 0., 0., myThid )
687 cdm #endif
688
689 #endif /* ALLOW_CTRL */
690
691 #endif /* undef ALLOW_ECCO) || def ECCO_CTRL_DEPRECATED */
692
693 #if (defined (ALLOW_CTRL) && defined (ALLOW_GENTIM2D_CONTROL))
694 IF ( useCTRL.AND.ctrlUseGen ) THEN
695 DO bj = myByLo(myThid),myByHi(myThid)
696 DO bi = myBxLo(myThid),mybxhi(myThid)
697 DO j = 1,sNy
698 DO i = 1,sNx
699 DO iarr = 1, maxCtrlTim2D
700 #ifdef ALLOW_ATM_TEMP
701 IF (xx_gentim2d_file(iarr)(1:8).EQ.'xx_atemp')
702 & atemp(i,j,bi,bj)=atemp(i,j,bi,bj)+
703 & xx_gentim2d(i,j,bi,bj,iarr)
704 IF (xx_gentim2d_file(iarr)(1:6).EQ.'xx_aqh')
705 & aqh(i,j,bi,bj)=aqh(i,j,bi,bj)+
706 & xx_gentim2d(i,j,bi,bj,iarr)
707 IF (xx_gentim2d_file(iarr)(1:9).EQ.'xx_precip')
708 & precip(i,j,bi,bj)=precip(i,j,bi,bj)+
709 & xx_gentim2d(i,j,bi,bj,iarr)
710 IF (xx_gentim2d_file(iarr)(1:9).EQ.'xx_lwflux')
711 & lwflux(i,j,bi,bj)=lwflux(i,j,bi,bj)+
712 & xx_gentim2d(i,j,bi,bj,iarr)
713 #endif
714 #if defined(ALLOW_ATM_TEMP) || defined(SHORTWAVE_HEATING)
715 IF (xx_gentim2d_file(iarr)(1:9).EQ.'xx_swflux')
716 & swflux(i,j,bi,bj)=swflux(i,j,bi,bj)+
717 & xx_gentim2d(i,j,bi,bj,iarr)
718 #endif
719 #ifdef ALLOW_DOWNWARD_RADIATION
720 IF (xx_gentim2d_file(iarr)(1:9).EQ.'xx_swdown')
721 & swdown(i,j,bi,bj)=swdown(i,j,bi,bj)+
722 & xx_gentim2d(i,j,bi,bj,iarr)
723 IF (xx_gentim2d_file(iarr)(1:9).EQ.'xx_lwdown')
724 & lwdown(i,j,bi,bj)=lwdown(i,j,bi,bj)+
725 & xx_gentim2d(i,j,bi,bj,iarr)
726 #endif
727 #ifdef ALLOW_RUNOFF
728 IF (xx_gentim2d_file(iarr)(1:9).EQ.'xx_runoff')
729 & runoff(i,j,bi,bj)=runoff(i,j,bi,bj)+
730 & xx_gentim2d(i,j,bi,bj,iarr)
731 #endif
732 #ifdef EXF_READ_EVAP
733 IF (xx_gentim2d_file(iarr)(1:7).EQ.'xx_evap')
734 & evap(i,j,bi,bj)=evap(i,j,bi,bj)+
735 & xx_gentim2d(i,j,bi,bj,iarr)
736 #endif
737 #ifdef ATMOSPHERIC_LOADING
738 IF (xx_gentim2d_file(iarr)(1:12).EQ.'xx_apressure')
739 & apressure(i,j,bi,bj)=apressure(i,j,bi,bj)+
740 & xx_gentim2d(i,j,bi,bj,iarr)
741 #endif
742 #ifdef EXF_SEAICE_FRACTION
743 IF (xx_gentim2d_file(iarr)(1:11).EQ.'xx_areamask')
744 & areamask(i,j,bi,bj)=areamask(i,j,bi,bj)+
745 & xx_gentim2d(i,j,bi,bj,iarr)
746 #endif
747 #ifndef ALLOW_ROTATE_UV_CONTROLS
748 IF (xx_gentim2d_file(iarr)(1:8).EQ.'xx_uwind')
749 & uwind(i,j,bi,bj)=uwind(i,j,bi,bj)+
750 & xx_gentim2d(i,j,bi,bj,iarr)
751 IF (xx_gentim2d_file(iarr)(1:8).EQ.'xx_vwind')
752 & vwind(i,j,bi,bj)=vwind(i,j,bi,bj)+
753 & xx_gentim2d(i,j,bi,bj,iarr)
754 #else
755 IF (xx_gentim2d_file(iarr)(1:8).EQ.'xx_uwind')
756 & tmpUE(i,j,bi,bj)=tmpUE(i,j,bi,bj)+
757 & xx_gentim2d(i,j,bi,bj,iarr)
758 IF (xx_gentim2d_file(iarr)(1:8).EQ.'xx_vwind')
759 & tmpVN(i,j,bi,bj)=tmpVN(i,j,bi,bj)+
760 & xx_gentim2d(i,j,bi,bj,iarr)
761 #endif
762 ENDDO
763 ENDDO
764 ENDDO
765 ENDDO
766 ENDDO
767 #ifdef ALLOW_ROTATE_UV_CONTROLS
768 CALL ROTATE_UV2EN_RL(tmpUX,tmpVY,tmpUE,tmpVN,
769 & .FALSE.,.FALSE.,.TRUE.,1,myThid)
770
771 DO bj = myByLo(myThid),myByHi(myThid)
772 DO bi = myBxLo(myThid),mybxhi(myThid)
773 DO j = 1,sNy
774 DO i = 1,sNx
775 uwind(i,j,bi,bj)=uwind(i,j,bi,bj)+tmpUX(i,j,bi,bj)
776 vwind(i,j,bi,bj)=vwind(i,j,bi,bj)+tmpVY(i,j,bi,bj)
777 ENDDO
778 ENDDO
779 ENDDO
780 ENDDO
781 #endif /* ALLOW_ROTATE_UV_CONTROLS */
782
783 ENDIF !if (ctrlUseGen) then
784 #endif
785
786 RETURN
787 END

  ViewVC Help
Powered by ViewVC 1.1.22