/[MITgcm]/MITgcm/pkg/exf/exf_getffields.F
ViewVC logotype

Contents of /MITgcm/pkg/exf/exf_getffields.F

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


Revision 1.66 - (show annotations) (download)
Wed Jan 11 03:45:34 2017 UTC (7 years, 4 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint66c
Changes since 1.65: +17 -2 lines
- add saltflx :: Net upward salt flux in psu.kg/m^2/s

1 C $Header: /u/gcmpack/MITgcm/pkg/exf/exf_getffields.F,v 1.65 2016/11/29 13:21:36 gforget 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-- read forcing fields from files and temporal interpolation
74
75 C Zonal and meridional wind stress.
76 IF ( .NOT.useAtmWind ) THEN
77 CALL EXF_SET_UV(
78 I ustressfile, ustressstartdate, ustressperiod,
79 I exf_inscal_ustress,
80 I ustress_exfremo_intercept, ustress_exfremo_slope,
81 U ustress, ustress0, ustress1, ustressmask,
82 I vstressfile, vstressstartdate, vstressperiod,
83 I exf_inscal_vstress,
84 I vstress_exfremo_intercept, vstress_exfremo_slope,
85 U vstress, vstress0, vstress1, vstressmask,
86 #ifdef USE_EXF_INTERPOLATION
87 I ustress_lon0, ustress_lon_inc, ustress_lat0, ustress_lat_inc,
88 I ustress_nlon, ustress_nlat, ustress_interpMethod,
89 I vstress_lon0, vstress_lon_inc, vstress_lat0, vstress_lat_inc,
90 I vstress_nlon, vstress_nlat, vstress_interpMethod,
91 I uvInterp_stress,
92 #endif /* USE_EXF_INTERPOLATION */
93 I myTime, myIter, myThid )
94 ELSE
95 DO bj = myByLo(myThid),myByHi(myThid)
96 DO bi = myBxLo(myThid),mybxhi(myThid)
97 DO j=1-OLy,sNy+OLy
98 DO i=1-OLx,sNx+OLx
99 ustress(i,j,bi,bj) = 0. _d 0
100 vstress(i,j,bi,bj) = 0. _d 0
101 ENDDO
102 ENDDO
103 ENDDO
104 ENDDO
105 ENDIF
106
107 C-- wind speed
108 CALL EXF_SET_GEN(
109 & wspeedfile, wspeedstartdate, wspeedperiod,
110 & exf_inscal_wspeed,
111 & wspeed_exfremo_intercept, wspeed_exfremo_slope,
112 & wspeed, wspeed0, wspeed1, wspeedmask,
113 #ifdef USE_EXF_INTERPOLATION
114 & wspeed_lon0, wspeed_lon_inc,
115 & wspeed_lat0, wspeed_lat_inc,
116 & wspeed_nlon, wspeed_nlat, xC, yC, wspeed_interpMethod,
117 #endif
118 & myTime, myIter, myThid )
119
120 C Zonal and meridional wind.
121 IF ( useAtmWind ) THEN
122 CALL EXF_SET_UV(
123 I uwindfile, uwindstartdate, uwindperiod,
124 I exf_inscal_uwind,
125 I uwind_exfremo_intercept, uwind_exfremo_slope,
126 U uwind, uwind0, uwind1, uwindmask,
127 I vwindfile, vwindstartdate, vwindperiod,
128 I exf_inscal_vwind,
129 I vwind_exfremo_intercept, vwind_exfremo_slope,
130 U vwind, vwind0, vwind1, vwindmask,
131 #ifdef USE_EXF_INTERPOLATION
132 I uwind_lon0, uwind_lon_inc, uwind_lat0, uwind_lat_inc,
133 I uwind_nlon, uwind_nlat, uwind_interpMethod,
134 I vwind_lon0, vwind_lon_inc, vwind_lat0, vwind_lat_inc,
135 I vwind_nlon, vwind_nlat, vwind_interpMethod, uvInterp_wind,
136 #endif /* USE_EXF_INTERPOLATION */
137 I myTime, myIter, myThid )
138
139 if (useRelativeWind) then
140 C Subtract UVEL and VVEL from UWIND and VWIND.
141 DO bj = myByLo(myThid),myByHi(myThid)
142 DO bi = myBxLo(myThid),mybxhi(myThid)
143 DO j = 1,sNy
144 DO i = 1,sNx
145 uwind(i,j,bi,bj) = uwind(i,j,bi,bj) - 0.5 _d 0 *
146 & (uVel(i,j,1,bi,bj)+uVel(i+1,j,1,bi,bj))
147 vwind(i,j,bi,bj) = vwind(i,j,bi,bj) - 0.5 _d 0 *
148 & (vVel(i,j,1,bi,bj)+vVel(i,j+1,1,bi,bj))
149 ENDDO
150 ENDDO
151 ENDDO
152 ENDDO
153 endif
154
155 ELSE
156 DO bj = myByLo(myThid),myByHi(myThid)
157 DO bi = myBxLo(myThid),mybxhi(myThid)
158 DO j=1-OLy,sNy+OLy
159 DO i=1-OLx,sNx+OLx
160 uwind(i,j,bi,bj) = 0. _d 0
161 vwind(i,j,bi,bj) = 0. _d 0
162 ENDDO
163 ENDDO
164 ENDDO
165 ENDDO
166 ENDIF
167
168 C Atmospheric heat flux.
169 CALL EXF_SET_GEN (
170 & hfluxfile, hfluxstartdate, hfluxperiod,
171 & exf_inscal_hflux,
172 & hflux_exfremo_intercept, hflux_exfremo_slope,
173 & hflux, hflux0, hflux1, hfluxmask,
174 #ifdef USE_EXF_INTERPOLATION
175 & hflux_lon0, hflux_lon_inc, hflux_lat0, hflux_lat_inc,
176 & hflux_nlon, hflux_nlat, xC, yC, hflux_interpMethod,
177 #endif
178 & myTime, myIter, myThid )
179
180 C Freshwater flux.
181 CALL EXF_SET_GEN (
182 & sfluxfile, sfluxstartdate, sfluxperiod,
183 & exf_inscal_sflux,
184 & sflux_exfremo_intercept, sflux_exfremo_slope,
185 & sflux, sflux0, sflux1, sfluxmask,
186 #ifdef USE_EXF_INTERPOLATION
187 & sflux_lon0, sflux_lon_inc, sflux_lat0, sflux_lat_inc,
188 & sflux_nlon, sflux_nlat, xC, yC, sflux_interpMethod,
189 #endif
190 & myTime, myIter, myThid )
191
192 #ifdef ALLOW_ATM_TEMP
193
194 C Atmospheric temperature.
195 CALL EXF_SET_GEN(
196 & atempfile, atempstartdate, atempperiod,
197 & exf_inscal_atemp,
198 & atemp_exfremo_intercept, atemp_exfremo_slope,
199 & atemp, atemp0, atemp1, atempmask,
200 #ifdef USE_EXF_INTERPOLATION
201 & atemp_lon0, atemp_lon_inc, atemp_lat0, atemp_lat_inc,
202 & atemp_nlon, atemp_nlat, xC, yC, atemp_interpMethod,
203 #endif
204 & myTime, myIter, myThid )
205 DO bj = myByLo(myThid),myByHi(myThid)
206 DO bi = myBxLo(myThid),mybxhi(myThid)
207 DO j = 1,sNy
208 DO i = 1,sNx
209 atemp(i,j,bi,bj) = atemp(i,j,bi,bj) + exf_offset_atemp
210 ENDDO
211 ENDDO
212 ENDDO
213 ENDDO
214
215 C Atmospheric humidity.
216 CALL EXF_SET_GEN(
217 & aqhfile, aqhstartdate, aqhperiod,
218 & exf_inscal_aqh,
219 & aqh_exfremo_intercept, aqh_exfremo_slope,
220 & aqh, aqh0, aqh1, aqhmask,
221 #ifdef USE_EXF_INTERPOLATION
222 & aqh_lon0, aqh_lon_inc, aqh_lat0, aqh_lat_inc,
223 & aqh_nlon, aqh_nlat, xC, yC, aqh_interpMethod,
224 #endif
225 & myTime, myIter, myThid )
226
227 C Net long wave radiative flux.
228 CALL EXF_SET_GEN(
229 & lwfluxfile, lwfluxstartdate, lwfluxperiod,
230 & exf_inscal_lwflux,
231 & lwflux_exfremo_intercept, lwflux_exfremo_slope,
232 & lwflux, lwflux0, lwflux1, lwfluxmask,
233 #ifdef USE_EXF_INTERPOLATION
234 & lwflux_lon0, lwflux_lon_inc, lwflux_lat0, lwflux_lat_inc,
235 & lwflux_nlon, lwflux_nlat, xC, yC, lwflux_interpMethod,
236 #endif
237 & myTime, myIter, myThid )
238
239 C Precipitation.
240 CALL EXF_SET_GEN(
241 & precipfile, precipstartdate, precipperiod,
242 & exf_inscal_precip,
243 & precip_exfremo_intercept, precip_exfremo_slope,
244 & precip, precip0, precip1, precipmask,
245 #ifdef USE_EXF_INTERPOLATION
246 & precip_lon0, precip_lon_inc, precip_lat0, precip_lat_inc,
247 & precip_nlon, precip_nlat, xC, yC, precip_interpMethod,
248 #endif
249 & myTime, myIter, myThid )
250
251 C Snow.
252 CALL EXF_SET_GEN(
253 & snowprecipfile, snowprecipstartdate, snowprecipperiod,
254 & exf_inscal_snowprecip,
255 & snowprecip_exfremo_intercept, snowprecip_exfremo_slope,
256 & snowprecip, snowprecip0, snowprecip1, snowprecipmask,
257 #ifdef USE_EXF_INTERPOLATION
258 & snowprecip_lon0, snowprecip_lon_inc,
259 & snowprecip_lat0, snowprecip_lat_inc,
260 & snowprecip_nlon, snowprecip_nlat, xC, yC,
261 & snowprecip_interpMethod,
262 #endif
263 & myTime, myIter, myThid )
264 C Take care of case where total precip is not defined
265 IF ( snowPrecipFile .NE. ' ' ) THEN
266 DO bj = myByLo(myThid),myByHi(myThid)
267 DO bi = myBxLo(myThid),mybxhi(myThid)
268 DO j = 1,sNy
269 DO i = 1,sNx
270 precip(i,j,bi,bj) =
271 & max( precip(i,j,bi,bj), snowPrecip(i,j,bi,bj) )
272 ENDDO
273 ENDDO
274 ENDDO
275 ENDDO
276 ENDIF
277
278 #endif /* ALLOW_ATM_TEMP */
279
280 #if defined(ALLOW_ATM_TEMP) || defined(SHORTWAVE_HEATING)
281 C Net short wave radiative flux.
282 CALL EXF_SET_GEN (
283 & swfluxfile, swfluxstartdate, swfluxperiod,
284 & exf_inscal_swflux,
285 & swflux_exfremo_intercept, swflux_exfremo_slope,
286 & swflux, swflux0, swflux1, swfluxmask,
287 #ifdef USE_EXF_INTERPOLATION
288 & swflux_lon0, swflux_lon_inc, swflux_lat0, swflux_lat_inc,
289 & swflux_nlon, swflux_nlat, xC, yC, swflux_interpMethod,
290 #endif
291 & myTime, myIter, myThid )
292 #endif
293
294 #ifdef EXF_READ_EVAP
295 C Evaporation
296 CALL EXF_SET_GEN (
297 & evapfile, evapstartdate, evapperiod,
298 & exf_inscal_evap,
299 & evap_exfremo_intercept, evap_exfremo_slope,
300 & evap, evap0, evap1, evapmask,
301 #ifdef USE_EXF_INTERPOLATION
302 & evap_lon0, evap_lon_inc, evap_lat0, evap_lat_inc,
303 & evap_nlon, evap_nlat, xC, yC, evap_interpMethod,
304 #endif
305 & myTime, myIter, myThid )
306 #endif
307
308 #ifdef ALLOW_DOWNWARD_RADIATION
309
310 C Downward shortwave radiation.
311 CALL EXF_SET_GEN (
312 & swdownfile, swdownstartdate, swdownperiod,
313 & exf_inscal_swdown,
314 & swdown_exfremo_intercept, swdown_exfremo_slope,
315 & swdown, swdown0, swdown1, swdownmask,
316 #ifdef USE_EXF_INTERPOLATION
317 & swdown_lon0, swdown_lon_inc, swdown_lat0, swdown_lat_inc,
318 & swdown_nlon, swdown_nlat, xC, yC, swdown_interpMethod,
319 #endif
320 & myTime, myIter, myThid )
321
322 C Downward longwave radiation.
323 CALL EXF_SET_GEN (
324 & lwdownfile, lwdownstartdate, lwdownperiod,
325 & exf_inscal_lwdown,
326 & lwdown_exfremo_intercept, lwdown_exfremo_slope,
327 & lwdown, lwdown0, lwdown1, lwdownmask,
328 #ifdef USE_EXF_INTERPOLATION
329 & lwdown_lon0, lwdown_lon_inc, lwdown_lat0, lwdown_lat_inc,
330 & lwdown_nlon, lwdown_nlat, xC, yC, lwdown_interpMethod,
331 #endif
332 & myTime, myIter, myThid )
333
334 #endif
335
336 #ifdef ATMOSPHERIC_LOADING
337 C Atmos. pressure forcing
338 CALL EXF_SET_GEN (
339 & apressurefile, apressurestartdate, apressureperiod,
340 & exf_inscal_apressure,
341 & apressure_exfremo_intercept, apressure_exfremo_slope,
342 & apressure, apressure0, apressure1, apressuremask,
343 #ifdef USE_EXF_INTERPOLATION
344 & apressure_lon0, apressure_lon_inc,
345 & apressure_lat0, apressure_lat_inc,
346 & apressure_nlon,apressure_nlat,xC,yC, apressure_interpMethod,
347 #endif
348 & myTime, myIter, myThid )
349 #endif
350
351 #ifdef EXF_SEAICE_FRACTION
352 C fractional ice-covered area MASK
353 CALL EXF_SET_GEN (
354 & areamaskfile, areamaskstartdate, areamaskperiod,
355 & exf_inscal_areamask,
356 & areamask_exfremo_intercept, areamask_exfremo_slope,
357 & areamask, areamask0, areamask1, areamaskmask,
358 #ifdef USE_EXF_INTERPOLATION
359 & areamask_lon0, areamask_lon_inc,
360 & areamask_lat0, areamask_lat_inc,
361 & areamask_nlon, areamask_nlat, xC, yC, areamask_interpMethod,
362 #endif
363 & myTime, myIter, myThid )
364 #endif
365
366 #ifdef ALLOW_RUNOFF
367 C Runoff
368 CALL EXF_SET_GEN (
369 & runofffile, runoffstartdate, runoffperiod,
370 & exf_inscal_runoff,
371 & runoff_exfremo_intercept, runoff_exfremo_slope,
372 & runoff, runoff0, runoff1, runoffmask,
373 #ifdef USE_EXF_INTERPOLATION
374 & runoff_lon0, runoff_lon_inc, runoff_lat0, runoff_lat_inc,
375 & runoff_nlon, runoff_nlat, xC, yC, runoff_interpMethod,
376 #endif
377 & myTime, myIter, myThid )
378 #endif /* ALLOW_RUNOFF */
379
380 #ifdef ALLOW_RUNOFTEMP
381 C Runoff temperature
382 CALL EXF_SET_GEN (
383 & runoftempfile, runoffstartdate, runoffperiod,
384 & exf_inscal_runoftemp,
385 & runoftemp_exfremo_intercept, runoftemp_exfremo_slope,
386 & runoftemp, runoftemp0, runoftemp1, runoffmask,
387 #ifdef USE_EXF_INTERPOLATION
388 & runoff_lon0, runoff_lon_inc, runoff_lat0, runoff_lat_inc,
389 & runoff_nlon, runoff_nlat, xC, yC, runoff_interpMethod,
390 #endif
391 & myTime, myIter, myThid )
392 #endif /* ALLOW_RUNOFTEMP */
393
394 #ifdef ALLOW_SALTFLX
395 C Salt flux
396 CALL EXF_SET_GEN (
397 & saltflxfile, saltflxstartdate, saltflxperiod,
398 & exf_inscal_saltflx,
399 & saltflx_exfremo_intercept, saltflx_exfremo_slope,
400 & saltflx, saltflx0, saltflx1, saltflxmask,
401 #ifdef USE_EXF_INTERPOLATION
402 & saltflx_lon0, saltflx_lon_inc,
403 & saltflx_lat0, saltflx_lat_inc,
404 & saltflx_nlon, saltflx_nlat, xC, yC, saltflx_interpMethod,
405 #endif
406 & myTime, myIter, myThid )
407 #endif
408
409 #ifdef ALLOW_ROTATE_UV_CONTROLS
410 if ( useCTRL ) then
411 DO bj = myByLo(myThid),myByHi(myThid)
412 DO bi = myBxLo(myThid),mybxhi(myThid)
413 DO j = 1-OLy,sNy+OLy
414 DO i = 1-OLx,sNx+OLx
415 tmpUE(i,j,bi,bj) = 0. _d 0
416 tmpVN(i,j,bi,bj) = 0. _d 0
417 tmpUX(i,j,bi,bj) = 0. _d 0
418 tmpVY(i,j,bi,bj) = 0. _d 0
419 ENDDO
420 ENDDO
421 ENDDO
422 ENDDO
423 endif
424 #endif
425
426 # if (!defined (ALLOW_ECCO) || defined (ECCO_CTRL_DEPRECATED))
427
428 C-- Control variables for atmos. state
429 #ifdef ALLOW_CTRL
430 if (.NOT.ctrlUseGen) then
431
432 #ifdef ALLOW_ATEMP_CONTROL
433 CALL CTRL_GET_GEN (
434 & xx_atemp_file, xx_atempstartdate, xx_atempperiod,
435 & maskc, atemp, xx_atemp0, xx_atemp1, xx_atemp_dummy,
436 & xx_atemp_remo_intercept, xx_atemp_remo_slope,
437 & watemp, myTime, myIter, myThid )
438 #endif
439
440 #ifdef ALLOW_AQH_CONTROL
441 CALL CTRL_GET_GEN (
442 & xx_aqh_file, xx_aqhstartdate, xx_aqhperiod,
443 & maskc, aqh, xx_aqh0, xx_aqh1, xx_aqh_dummy,
444 & xx_aqh_remo_intercept, xx_aqh_remo_slope,
445 & waqh, myTime, myIter, myThid )
446 #endif
447
448 #ifdef ALLOW_PRECIP_CONTROL
449 CALL CTRL_GET_GEN (
450 & xx_precip_file, xx_precipstartdate, xx_precipperiod,
451 & maskc, precip, xx_precip0, xx_precip1, xx_precip_dummy,
452 & xx_precip_remo_intercept, xx_precip_remo_slope,
453 & wprecip, myTime, myIter, myThid )
454 #endif
455
456 #ifdef ALLOW_SWDOWN_CONTROL
457 CALL CTRL_GET_GEN (
458 & xx_swdown_file, xx_swdownstartdate, xx_swdownperiod,
459 & maskc, swdown, xx_swdown0, xx_swdown1, xx_swdown_dummy,
460 & xx_swdown_remo_intercept, xx_swdown_remo_slope,
461 & wswdown, myTime, myIter, myThid )
462 #endif
463
464 #ifdef ALLOW_LWDOWN_CONTROL
465 CALL CTRL_GET_GEN (
466 & xx_lwdown_file, xx_lwdownstartdate, xx_lwdownperiod,
467 & maskc, lwdown, xx_lwdown0, xx_lwdown1, xx_lwdown_dummy,
468 & xx_lwdown_remo_intercept, xx_lwdown_remo_slope,
469 & wlwdown, myTime, myIter, myThid )
470 #endif
471
472 endif !if (.NOT.ctrlUseGen) then
473
474 #ifdef ALLOW_SWFLUX_CONTROL
475 CALL CTRL_GET_GEN (
476 & xx_swflux_file, xx_swfluxstartdate, xx_swfluxperiod,
477 & maskc, swflux, xx_swflux0, xx_swflux1, xx_swflux_dummy,
478 & xx_swflux_remo_intercept, xx_swflux_remo_slope,
479 & wswflux, myTime, myIter, myThid )
480 #endif
481
482 #ifdef ALLOW_LWFLUX_CONTROL
483 CALL CTRL_GET_GEN (
484 & xx_lwflux_file, xx_lwfluxstartdate, xx_lwfluxperiod,
485 & maskc, lwflux, xx_lwflux0, xx_lwflux1, xx_lwflux_dummy,
486 & xx_lwflux_remo_intercept, xx_lwflux_remo_slope,
487 & wswflux, myTime, myIter, myThid )
488 #endif
489
490 #ifdef ALLOW_EVAP_CONTROL
491 CALL CTRL_GET_GEN (
492 & xx_evap_file, xx_evapstartdate, xx_evapperiod,
493 & maskc, evap, xx_evap0, xx_evap1, xx_evap_dummy,
494 & xx_evap_remo_intercept, xx_evap_remo_slope,
495 & wevap, myTime, myIter, myThid )
496 #endif
497
498 #ifdef ALLOW_SNOWPRECIP_CONTROL
499 CALL CTRL_GET_GEN (
500 & xx_snowprecip_file, xx_snowprecipstartdate,
501 & xx_snowprecipperiod,
502 & maskc, snowprecip, xx_snowprecip0, xx_snowprecip1,
503 & xx_snowprecip_dummy,
504 & xx_snowprecip_remo_intercept, xx_snowprecip_remo_slope,
505 & wsnowprecip, myTime, myIter, myThid )
506 #endif
507
508 #ifdef ALLOW_APRESSURE_CONTROL
509 CALL CTRL_GET_GEN (
510 & xx_apressure_file, xx_apressurestartdate,
511 & xx_apressureperiod,
512 & maskc, apressure, xx_apressure0, xx_apressure1,
513 & xx_apressure_dummy,
514 & xx_apressure_remo_intercept, xx_apressure_remo_slope,
515 & wapressure, myTime, myIter, myThid )
516 #endif
517
518 IF ( useAtmWind ) THEN
519 #ifndef ALLOW_ROTATE_UV_CONTROLS
520
521 #ifdef ALLOW_UWIND_CONTROL
522 CALL CTRL_GET_GEN (
523 & xx_uwind_file, xx_uwindstartdate, xx_uwindperiod,
524 & maskc, uwind, xx_uwind0, xx_uwind1, xx_uwind_dummy,
525 & xx_uwind_remo_intercept, xx_uwind_remo_slope,
526 & wuwind, myTime, myIter, myThid )
527 #endif /* ALLOW_UWIND_CONTROL */
528
529 #ifdef ALLOW_VWIND_CONTROL
530 CALL CTRL_GET_GEN (
531 & xx_vwind_file, xx_vwindstartdate, xx_vwindperiod,
532 & maskc, vwind, xx_vwind0, xx_vwind1, xx_vwind_dummy,
533 & xx_vwind_remo_intercept, xx_vwind_remo_slope,
534 & wvwind, myTime, myIter, myThid )
535 #endif /* ALLOW_VWIND_CONTROL */
536
537 #else
538
539 #if defined(ALLOW_UWIND_CONTROL) && defined(ALLOW_VWIND_CONTROL)
540 CALL CTRL_GET_GEN (
541 & xx_uwind_file, xx_uwindstartdate, xx_uwindperiod,
542 & maskc, tmpUE, xx_uwind0, xx_uwind1, xx_uwind_dummy,
543 & xx_uwind_remo_intercept, xx_uwind_remo_slope,
544 & wuwind, myTime, myIter, myThid )
545
546 CALL CTRL_GET_GEN (
547 & xx_vwind_file, xx_vwindstartdate, xx_vwindperiod,
548 & maskc, tmpVN, xx_vwind0, xx_vwind1, xx_vwind_dummy,
549 & xx_vwind_remo_intercept, xx_vwind_remo_slope,
550 & wvwind, myTime, myIter, myThid )
551
552 CALL ROTATE_UV2EN_RL(tmpUX,tmpVY,tmpUE,tmpVN,
553 & .FALSE.,.FALSE.,.TRUE.,1,myThid)
554
555 DO bj = myByLo(myThid),myByHi(myThid)
556 DO bi = myBxLo(myThid),mybxhi(myThid)
557 DO j = 1,sNy
558 DO i = 1,sNx
559 uwind(i,j,bi,bj)=uwind(i,j,bi,bj)+tmpUX(i,j,bi,bj)
560 vwind(i,j,bi,bj)=vwind(i,j,bi,bj)+tmpVY(i,j,bi,bj)
561 ENDDO
562 ENDDO
563 ENDDO
564 ENDDO
565 #endif
566
567 #endif /* ALLOW_ROTATE_UV_CONTROLS */
568 ENDIF
569
570 #ifdef ALLOW_ATM_MEAN_CONTROL
571 DO bj = myByLo(myThid),myByHi(myThid)
572 DO bi = myBxLo(myThid),mybxhi(myThid)
573 DO j = 1,sNy
574 DO i = 1,sNx
575 # ifdef ALLOW_ATEMP_CONTROL
576 atemp(i,j,bi,bj) =atemp(i,j,bi,bj) +xx_atemp_mean(i,j,bi,bj)
577 # endif
578 # ifdef ALLOW_AQH_CONTROL
579 aqh(i,j,bi,bj) =aqh(i,j,bi,bj) +xx_aqh_mean(i,j,bi,bj)
580 # endif
581 # ifdef ALLOW_PRECIP_CONTROL
582 precip(i,j,bi,bj)=precip(i,j,bi,bj)+xx_precip_mean(i,j,bi,bj)
583 # endif
584 # ifdef ALLOW_SWDOWN_CONTROL
585 swdown(i,j,bi,bj)=swdown(i,j,bi,bj)+xx_swdown_mean(i,j,bi,bj)
586 # endif
587 # ifdef ALLOW_UWIND_CONTROL
588 uwind(i,j,bi,bj) =uwind(i,j,bi,bj) +xx_uwind_mean(i,j,bi,bj)
589 # endif
590 # ifdef ALLOW_VWIND_CONTROL
591 vwind(i,j,bi,bj) =vwind(i,j,bi,bj) +xx_vwind_mean(i,j,bi,bj)
592 # endif
593 ENDDO
594 ENDDO
595 ENDDO
596 ENDDO
597 #endif /* ALLOW_ATM_MEAN_CONTROL */
598
599 cdm transferred from exf_init_runoff.F
600 cdm functionality needs to be checked before turning on
601 cdm #ifdef ALLOW_RUNOFF_CONTROL
602 cdm CALL CTRL_GET_GEN (
603 cdm & xx_runoff_file, xx_runoffstartdate, xx_runoffperiod,
604 cdm & maskc, runoff, xx_runoff0, xx_runoff1, xx_runoff_dummy,
605 cdm & xx_runoff_remo_intercept, xx_runoff_remo_slope,
606 cdm & wrunoff, 0., 0., myThid )
607 cdm #endif
608
609 #endif /* ALLOW_CTRL */
610
611 #endif /* undef ALLOW_ECCO) || def ECCO_CTRL_DEPRECATED */
612
613 #if (defined (ALLOW_CTRL) && defined (ALLOW_GENTIM2D_CONTROL))
614 if ( useCTRL.AND.ctrlUseGen ) then
615 DO bj = myByLo(myThid),myByHi(myThid)
616 DO bi = myBxLo(myThid),mybxhi(myThid)
617 DO j = 1,sNy
618 DO i = 1,sNx
619 do iarr = 1, maxCtrlTim2D
620 #ifdef ALLOW_ATM_TEMP
621 if (xx_gentim2d_file(iarr)(1:8).EQ.'xx_atemp')
622 & atemp(i,j,bi,bj)=atemp(i,j,bi,bj)+
623 & xx_gentim2d(i,j,bi,bj,iarr)
624 if (xx_gentim2d_file(iarr)(1:6).EQ.'xx_aqh')
625 & aqh(i,j,bi,bj)=aqh(i,j,bi,bj)+
626 & xx_gentim2d(i,j,bi,bj,iarr)
627 if (xx_gentim2d_file(iarr)(1:9).EQ.'xx_precip')
628 & precip(i,j,bi,bj)=precip(i,j,bi,bj)+
629 & xx_gentim2d(i,j,bi,bj,iarr)
630 #endif
631 #ifdef ALLOW_DOWNWARD_RADIATION
632 if (xx_gentim2d_file(iarr)(1:9).EQ.'xx_swdown')
633 & swdown(i,j,bi,bj)=swdown(i,j,bi,bj)+
634 & xx_gentim2d(i,j,bi,bj,iarr)
635 if (xx_gentim2d_file(iarr)(1:9).EQ.'xx_lwdown')
636 & lwdown(i,j,bi,bj)=lwdown(i,j,bi,bj)+
637 & xx_gentim2d(i,j,bi,bj,iarr)
638 #endif
639 if (xx_gentim2d_file(iarr)(1:9).EQ.'xx_runoff')
640 & runoff(i,j,bi,bj)=runoff(i,j,bi,bj)+
641 & xx_gentim2d(i,j,bi,bj,iarr)
642 #ifndef ALLOW_ROTATE_UV_CONTROLS
643 if (xx_gentim2d_file(iarr)(1:8).EQ.'xx_uwind')
644 & uwind(i,j,bi,bj)=uwind(i,j,bi,bj)+
645 & xx_gentim2d(i,j,bi,bj,iarr)
646 if (xx_gentim2d_file(iarr)(1:8).EQ.'xx_vwind')
647 & vwind(i,j,bi,bj)=vwind(i,j,bi,bj)+
648 & xx_gentim2d(i,j,bi,bj,iarr)
649 #else
650 if (xx_gentim2d_file(iarr)(1:8).EQ.'xx_uwind')
651 & tmpUE(i,j,bi,bj)=tmpUE(i,j,bi,bj)+
652 & xx_gentim2d(i,j,bi,bj,iarr)
653 if (xx_gentim2d_file(iarr)(1:8).EQ.'xx_vwind')
654 & tmpVN(i,j,bi,bj)=tmpVN(i,j,bi,bj)+
655 & xx_gentim2d(i,j,bi,bj,iarr)
656
657 #endif
658 enddo
659 ENDDO
660 ENDDO
661 ENDDO
662 ENDDO
663 #ifdef ALLOW_ROTATE_UV_CONTROLS
664 CALL ROTATE_UV2EN_RL(tmpUX,tmpVY,tmpUE,tmpVN,
665 & .FALSE.,.FALSE.,.TRUE.,1,myThid)
666
667 DO bj = myByLo(myThid),myByHi(myThid)
668 DO bi = myBxLo(myThid),mybxhi(myThid)
669 DO j = 1,sNy
670 DO i = 1,sNx
671 uwind(i,j,bi,bj)=uwind(i,j,bi,bj)+tmpUX(i,j,bi,bj)
672 vwind(i,j,bi,bj)=vwind(i,j,bi,bj)+tmpVY(i,j,bi,bj)
673 ENDDO
674 ENDDO
675 ENDDO
676 ENDDO
677 #endif /* ALLOW_ROTATE_UV_CONTROLS */
678
679 endif !if (ctrlUseGen) then
680 #endif
681
682 RETURN
683 END

  ViewVC Help
Powered by ViewVC 1.1.22