/[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.57 - (show annotations) (download)
Mon Oct 20 03:13:32 2014 UTC (9 years, 7 months ago) by gforget
Branch: MAIN
Changes since 1.56: +7 -1 lines
- ECCO_OPTIONS.h is needed when including ecco_cost.h, ecco.h
- AUTODIFF_OPTIONS.h is needed when including tamc.h, tamc_keys.h
- CTRL_OPTIONS.h is needed when including ctrl.h

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

  ViewVC Help
Powered by ViewVC 1.1.22