/[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.48 - (show annotations) (download)
Tue Jul 31 16:08:16 2012 UTC (11 years, 10 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint63r
Changes since 1.47: +2 -1 lines
Attempt at adding CTRL_SIZE.h

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

  ViewVC Help
Powered by ViewVC 1.1.22