5.13.2-dev0
mHM
The mesoscale Hydrological Model
Loading...
Searching...
No Matches
mo_meteo_handler.f90
Go to the documentation of this file.
1!> \dir meteo
2!> \brief \copybrief f_meteo
3!> \details \copydetails f_meteo
4
5!> \defgroup f_meteo Meteo - Fortran modules
6!> \brief Core modules to deal with meteorological forcings.
7!> \details These modules provide the meteo handler, the spatial and temporal remapping algorithms and helper routines.
8
9!> \file mo_meteo_handler.f90
10!> \brief \copybrief mo_meteo_handler
11!> \details \copydetails mo_meteo_handler
12
13!> \brief Class for the meteo handler
14!> \details Handler for meteorological forcings in mHM.
15!! Is independent of global variables and provides 3 methods to access forcings:
16!! - get_corrected_pet : get the modified pet for mHM for the current timestep
17!! - get_temp : get the temporal disaggregated temperature for the current timestep
18!! - get_prec : get the temporal disaggregated percipitation for the current timestep
19!!
20!> \version 0.1
21!> \authors Sebastian Mueller
22!> \date Mar 2023
23!> \copyright Copyright 2005-\today, the mHM Developers, Luis Samaniego, Sabine Attinger: All rights reserved.
24!! mHM is released under the LGPLv3+ license \license_note
25!> \ingroup f_meteo
27
28 USE mo_kind, ONLY : i4, dp
29 USE mo_constants, ONLY : yearmonths
30 use mo_common_types, only: grid, period
31 use mo_message, only : message, error_message
33 use mo_sentinel, only : set_sentinel, check_sentinel
34 use mo_datetime, only : datetime, timedelta, zero_delta, one_hour, one_day
35
36 implicit none
37
38 private
39
40 !> \class meteo_handler_type
41 !> \brief This is a handler for the meteorological forcings
42 !> \details This class provides all procedures to handle meteorological forcings from file or interfaces.
43 type, public :: meteo_handler_type
44 ! config settings
45 character(256) :: dir_nml_name = 'directories_mHM' !< namelist name in mhm namelist
46 character(256) :: weight_nml_name = 'nightDayRatio' !< namelist name in mhm namelist
47 !> Input nCols and nRows of binary meteo and LAI files are in header file
48 CHARACTER(256) :: file_meteo_header = 'header.txt'
49 !> Unit for meteo header file
50 INTEGER :: umeteo_header = 50
51 !> .FALSE. to only warn about bound (lower, upper) violations in meteo files, default = .TRUE. - raise an error
52 logical :: bound_error
53 integer(i4) :: pet_case !< process case for PET (processCase(5))
54 integer(i4) :: riv_temp_case !< process case for river temperature (processCase(11))
55 type(period), public :: readper !< start and end dates of read period
56 integer(i4), public :: ntstepforcingday !< Number of forcing intervals per day
57 !> flag wether forcings are given at hourly timestep
58 logical, public :: is_hourly_forcing
59 integer(i4), public :: timestep !< [h] simulation time step (= TS) in [h]
60 integer(i4), public :: ntstepday !< Number of time intervals per day
61 real(dp), public :: ntstepday_dp !< Number of time intervals per day as real
62 ! -------------------------------------------------------------------
63 ! level 2 description
64 ! -------------------------------------------------------------------
65 type(grid), dimension(:), allocatable, public :: level2 !< Reference of the metereological variables
66 ! -------------------------------------------------------------------
67 ! INPUT variables for configuration of mHM
68 ! -------------------------------------------------------------------
69 integer(i4), dimension(:), allocatable, public :: timestep_model_inputs !< frequency for reading meteo input
70 logical, public :: read_meteo_weights !< read weights for tavg and pet
71 character(256), public :: inputformat_meteo_forcings !< format of meteo input data (nc)
72 ! ------------------------------------------------------------------
73 ! DIRECTORIES
74 ! ------------------------------------------------------------------
75 ! has the dimension of nDomains
76 character(256), dimension(:), allocatable, public :: dir_meteo_header !< Directory where the meteo header file is located
77 character(256), dimension(:), allocatable, public :: dirprecipitation !< Directory where precipitation files are located
78 character(256), dimension(:), allocatable, public :: dirtemperature !< Directory where temperature files are located
79 character(256), dimension(:), allocatable, public :: dirmintemperature !< Directory where minimum temp. files are located
80 character(256), dimension(:), allocatable, public :: dirmaxtemperature !< Directory where maximum temp. files are located
81 character(256), dimension(:), allocatable, public :: dirnetradiation !< Directory where abs. vap. pressure files are located
82 character(256), dimension(:), allocatable, public :: dirabsvappressure !< Directory where abs. vap. pressure files are located
83 character(256), dimension(:), allocatable, public :: dirwindspeed !< Directory where windspeed files are located
84 character(256), dimension(:), allocatable, public :: dirreferenceet !< Directory where reference-ET files are located
85 ! riv-temp releated
86 character(256), dimension(:), allocatable, public :: dirradiation !< Directory where short/long-wave rad. files are located
87 ! -------------------------------------------------------------------
88 ! L1 DOMAIN description
89 ! -------------------------------------------------------------------
90 ! Forcings
91 ! dim1 = number grid cells L1
92 ! dim2 = number of meteorological time steps
93 ! dim2 = month of year
94 real(dp), public, dimension(:, :, :), allocatable :: l1_temp_weights !< hourly temperature weights for daily values
95 real(dp), public, dimension(:, :, :), allocatable :: l1_pet_weights !< hourly pet weights for daily values
96 real(dp), public, dimension(:, :, :), allocatable :: l1_pre_weights !< hourly pre weights for daily values
97 real(dp), public, dimension(:, :), allocatable :: l1_pre !< [mm] Precipitation
98 real(dp), public, dimension(:, :), allocatable :: l1_temp !< [degC] Air temperature
99 real(dp), public, dimension(:, :), allocatable :: l1_pet !< [mm TS-1] Potential evapotranspiration
100 real(dp), public, dimension(:, :), allocatable :: l1_tmin !< [degC] minimum daily air temperature
101 real(dp), public, dimension(:, :), allocatable :: l1_tmax !< [degC] maximum daily air temperature
102 real(dp), public, dimension(:, :), allocatable :: l1_netrad !< [W m2] net radiation
103 real(dp), public, dimension(:, :), allocatable :: l1_absvappress !< [Pa] absolute vapour pressure
104 real(dp), public, dimension(:, :), allocatable :: l1_windspeed !< [m s-1] windspeed
105 ! riv-temp related
106 real(dp), public, dimension(:, :), allocatable :: l1_ssrd !< [W m2] short wave radiation
107 real(dp), public, dimension(:, :), allocatable :: l1_strd !< [W m2] long wave radiation
108 real(dp), public, dimension(:, :), allocatable :: l1_tann !< [degC] annual mean air temperature
109 ! -------------------------------------------------------------------
110 ! Monthly day/night variation of Meteorological variables
111 ! for temporal disaggregation
112 ! -------------------------------------------------------------------
113 ! dim1 = number of months in a year
114 real(dp), public, dimension(int(YearMonths, i4)) :: fday_prec !< [-] Day ratio precipitation < 1
115 real(dp), public, dimension(int(YearMonths, i4)) :: fnight_prec !< [-] Night ratio precipitation < 1
116 real(dp), public, dimension(int(YearMonths, i4)) :: fday_pet !< [-] Day ratio PET < 1
117 real(dp), public, dimension(int(YearMonths, i4)) :: fnight_pet !< [-] Night ratio PET < 1
118 real(dp), public, dimension(int(YearMonths, i4)) :: fday_temp !< [-] Day factor mean temp
119 real(dp), public, dimension(int(YearMonths, i4)) :: fnight_temp !< [-] Night factor mean temp
120 real(dp), public, dimension(int(YearMonths, i4)) :: fday_ssrd !< [-] Day factor short-wave rad.
121 real(dp), public, dimension(int(YearMonths, i4)) :: fnight_ssrd !< [-] Night factor short-wave rad.
122 real(dp), public, dimension(int(YearMonths, i4)) :: fday_strd !< [-] Day factor long-wave rad.
123 real(dp), public, dimension(int(YearMonths, i4)) :: fnight_strd !< [-] Night factor long-wave rad.
124
125 ! -------------------------------------------------------------------
126 ! current mHM array indizes
127 ! -------------------------------------------------------------------
128 !> start index of meteo variables
129 integer(i4) :: s_meteo
130 !> end index of meteo variables
131 integer(i4) :: e_meteo
132 !> index of meteo time-step
133 integer(i4) :: imeteots
134 !> start index of level-1 variables
135 integer(i4) :: s1
136 !> end index of level-1 variables
137 integer(i4) :: e1
138 !> number of domains
139 integer(i4) :: ndomains
140 integer(i4), dimension(:), allocatable :: indices !< indices
141 integer(i4), dimension(:), allocatable :: l0datafrom !< index of associated level-0 domain
142 !> current julian time
143 real(dp) :: time
144
145 ! -------------------------------------------------------------------
146 ! coupling settings
147 ! -------------------------------------------------------------------
148 type(couple_cfg_type), public :: couple_cfg !< coupling configuration class
149 type(timedelta) :: couple_step_delta !< timedelta for the coupling meteo time-step
150 type(datetime) :: couple_pre_time !< current time from coupling for pre
151 type(datetime) :: couple_temp_time !< current time from coupling for temp
152 type(datetime) :: couple_pet_time !< current time from coupling for pet
153 type(datetime) :: couple_tmin_time !< current time from coupling for tmin
154 type(datetime) :: couple_tmax_time !< current time from coupling for tmax
155 type(datetime) :: couple_netrad_time !< current time from coupling for netrad
156 type(datetime) :: couple_absvappress_time !< current time from coupling for absvappress
157 type(datetime) :: couple_windspeed_time !< current time from coupling for windspeed
158 type(datetime) :: couple_ssrd_time !< current time from coupling for ssrd
159 type(datetime) :: couple_strd_time !< current time from coupling for strd
160 type(datetime) :: couple_tann_time !< current time from coupling for tann
161 logical, public :: couple_pre !< coupling config for pre
162 logical, public :: couple_temp !< coupling config for temp
163 logical, public :: couple_pet !< coupling config for pet
164 logical, public :: couple_tmin !< coupling config for tmin
165 logical, public :: couple_tmax !< coupling config for tmax
166 logical, public :: couple_netrad !< coupling config for netrad
167 logical, public :: couple_absvappress !< coupling config for absvappress
168 logical, public :: couple_windspeed !< coupling config for windspeed
169 logical, public :: couple_ssrd !< coupling config for ssrd
170 logical, public :: couple_strd !< coupling config for strd
171 logical, public :: couple_tann !< coupling config for tann
172 logical, public :: couple_all !< flag to indicated that all meteo-data is coming from the coupler
173 logical, public :: couple_is_hourly !< flag to indicated hourly data from coupler
174 contains
175 !> \copydoc mo_meteo_handler::clean_up
176 procedure :: clean_up !< \see mo_meteo_handler::clean_up
177 !> \copydoc mo_meteo_handler::config
178 procedure :: config !< \see mo_meteo_handler::config
179 !> \copydoc mo_meteo_handler::single_read
180 procedure :: single_read !< \see mo_meteo_handler::single_read
181 !> \copydoc mo_meteo_handler::init_level2
182 procedure :: init_level2 !< \see mo_meteo_handler::init_level2
183 !> \copydoc mo_meteo_handler::prepare_data
184 procedure :: prepare_data !< \see mo_meteo_handler::prepare_data
185 !> \copydoc mo_meteo_handler::update_timestep
186 procedure :: update_timestep !< \see mo_meteo_handler::update_timestep
187 !> \copydoc mo_meteo_handler::get_corrected_pet
188 procedure :: get_corrected_pet !< \see mo_meteo_handler::get_corrected_pet
189 !> \copydoc mo_meteo_handler::get_temp
190 procedure :: get_temp !< \see mo_meteo_handler::get_temp
191 !> \copydoc mo_meteo_handler::get_prec
192 procedure :: get_prec !< \see mo_meteo_handler::get_prec
193 !> \copydoc mo_meteo_handler::get_ssrd
194 procedure :: get_ssrd !< \see mo_meteo_handler::get_ssrd
195 !> \copydoc mo_meteo_handler::get_strd
196 procedure :: get_strd !< \see mo_meteo_handler::get_strd
197 !> \copydoc mo_meteo_handler::get_tann
198 procedure :: get_tann !< \see mo_meteo_handler::get_tann
199 !> \copydoc mo_meteo_handler::set_meteo
200 procedure :: set_meteo !< \see mo_meteo_handler::set_meteo
201 end type meteo_handler_type
202
203contains
204
205 !> \brief clean up
206 subroutine clean_up(self)
207 implicit none
208
209 class(meteo_handler_type), intent(inout) :: self
210
211 if ( allocated(self%indices) ) deallocate(self%indices)
212 if ( allocated(self%L0DataFrom) ) deallocate(self%L0DataFrom)
213 if ( allocated(self%timeStep_model_inputs) ) deallocate(self%timeStep_model_inputs)
214 if ( allocated(self%dir_meteo_header) ) deallocate(self%dir_meteo_header)
215 if ( allocated(self%dirPrecipitation) ) deallocate(self%dirPrecipitation)
216 if ( allocated(self%dirTemperature) ) deallocate(self%dirTemperature)
217 if ( allocated(self%dirMinTemperature) ) deallocate(self%dirMinTemperature)
218 if ( allocated(self%dirMaxTemperature) ) deallocate(self%dirMaxTemperature)
219 if ( allocated(self%dirNetRadiation) ) deallocate(self%dirNetRadiation)
220 if ( allocated(self%dirabsVapPressure) ) deallocate(self%dirabsVapPressure)
221 if ( allocated(self%dirwindspeed) ) deallocate(self%dirwindspeed)
222 if ( allocated(self%dirReferenceET) ) deallocate(self%dirReferenceET)
223 if ( allocated(self%dirRadiation) ) deallocate(self%dirRadiation)
224 if ( allocated(self%level2) ) deallocate(self%level2)
225 if ( allocated(self%L1_temp_weights) ) deallocate(self%L1_temp_weights)
226 if ( allocated(self%L1_pet_weights) ) deallocate(self%L1_pet_weights)
227 if ( allocated(self%L1_pre_weights) ) deallocate(self%L1_pre_weights)
228 if ( allocated(self%L1_pre) ) deallocate(self%L1_pre)
229 if ( allocated(self%L1_temp) ) deallocate(self%L1_temp)
230 if ( allocated(self%L1_pet) ) deallocate(self%L1_pet)
231 if ( allocated(self%L1_tmin) ) deallocate(self%L1_tmin)
232 if ( allocated(self%L1_tmax) ) deallocate(self%L1_tmax)
233 if ( allocated(self%L1_netrad) ) deallocate(self%L1_netrad)
234 if ( allocated(self%L1_absvappress) ) deallocate(self%L1_absvappress)
235 if ( allocated(self%L1_windspeed) ) deallocate(self%L1_windspeed)
236 if ( allocated(self%L1_ssrd) ) deallocate(self%L1_ssrd)
237 if ( allocated(self%L1_strd) ) deallocate(self%L1_strd)
238 if ( allocated(self%L1_tann) ) deallocate(self%L1_tann)
239
240 end subroutine clean_up
241
242 !> \brief configure the \ref meteo_handler_type class from the mhm namelist
243 subroutine config(self, file_namelist, unamelist, optimize, domainMeta, processMatrix, timeStep, couple_cfg)
244
246 use mo_common_types, only : domain_meta
247 use mo_nml, only : close_nml, open_nml, position_nml
249
250 implicit none
251
252 class(meteo_handler_type), intent(inout) :: self
253 character(*), intent(in) :: file_namelist !< mhm namelist file
254 integer, intent(in) :: unamelist !< unit to open namelist file
255 logical, intent(in) :: optimize !< Optimization flag
256 type(domain_meta), intent(in) :: domainMeta !< domain general description
257 integer(i4), dimension(nProcesses, 3), intent(in) :: processMatrix !< Info about which process runs in which option
258 integer(i4), intent(in) :: timeStep !< [h] simulation time step (= TS) in [h]
259 type(couple_cfg_type), intent(in) :: couple_cfg !< coupling configuration class
260
261 integer(i4), dimension(maxNoDomains) :: time_step_model_inputs
262 character(256), dimension(maxNoDomains) :: dir_meteo_header
263 character(256), dimension(maxNoDomains) :: dir_Precipitation
264 character(256), dimension(maxNoDomains) :: dir_Temperature
265 character(256), dimension(maxNoDomains) :: dir_MinTemperature
266 character(256), dimension(maxNoDomains) :: dir_MaxTemperature
267 character(256), dimension(maxNoDomains) :: dir_NetRadiation
268 character(256), dimension(maxNoDomains) :: dir_windspeed
269 character(256), dimension(maxNoDomains) :: dir_absVapPressure
270 character(256), dimension(maxNoDomains) :: dir_ReferenceET
271 character(256), dimension(maxNoDomains) :: dir_Radiation ! riv-temp related
272 character(256) :: inputFormat_meteo_forcings
273
274 logical :: read_meteo_weights, bound_error
275 real(dp), dimension(int(YearMonths, i4)) :: fnight_prec
276 real(dp), dimension(int(YearMonths, i4)) :: fnight_pet
277 real(dp), dimension(int(YearMonths, i4)) :: fnight_temp
278 real(dp), dimension(int(YearMonths, i4)) :: fnight_ssrd
279 real(dp), dimension(int(YearMonths, i4)) :: fnight_strd
280
281 integer(i4) :: domainID, iDomain
282
283 ! namelist directories
284 namelist /directories_mhm/ &
285 inputformat_meteo_forcings, &
286 bound_error, &
287 dir_meteo_header, &
288 dir_precipitation, &
289 dir_temperature, &
290 dir_referenceet, &
291 dir_mintemperature, &
292 dir_maxtemperature, &
293 dir_absvappressure, &
294 dir_windspeed, &
295 dir_netradiation, &
296 dir_radiation, &
297 time_step_model_inputs
298
299 ! namelist for night-day ratio of precipitation, referenceET and temperature
300 namelist /nightdayratio/ &
301 read_meteo_weights, &
302 fnight_prec, &
303 fnight_pet, &
304 fnight_temp, &
305 fnight_ssrd, &
306 fnight_strd
307
308 ! store coupling config
309 self%couple_cfg = couple_cfg
310
311 ! store needed domain meta infos
312 self%nDomains = domainmeta%nDomains
313 allocate(self%indices(self%nDomains))
314 allocate(self%L0DataFrom(self%nDomains))
315 self%indices(:) = domainmeta%indices(:)
316 self%L0DataFrom(:) = domainmeta%L0DataFrom(:)
317
318 ! # init of number of forcing timesteps, will be set when reading forcings
319 self%nTStepForcingDay = nodata_i4
320
321 ! store important process cases
322 self%pet_case = processmatrix(5,1)
323 self%riv_temp_case = processmatrix(11,1)
324 ! store time-stepping info
325 self%timeStep = timestep
326 self%nTStepDay = 24_i4 / timestep ! # of time steps per day
327 self%nTstepDay_dp = real(self%nTStepDay, dp)
328
329 ! allocate variables
330 allocate(self%dir_meteo_header(self%nDomains))
331 allocate(self%dirPrecipitation(self%nDomains))
332 allocate(self%dirTemperature(self%nDomains))
333 allocate(self%dirwindspeed(self%nDomains))
334 allocate(self%dirabsVapPressure(self%nDomains))
335 allocate(self%dirReferenceET(self%nDomains))
336 allocate(self%dirMinTemperature(self%nDomains))
337 allocate(self%dirMaxTemperature(self%nDomains))
338 allocate(self%dirNetRadiation(self%nDomains))
339 allocate(self%dirRadiation(self%nDomains))
340 ! allocate time periods
341 allocate(self%timestep_model_inputs(self%nDomains))
342
343 ! open the namelist file
344 call open_nml(file_namelist, unamelist, quiet=.true.)
345
346 !===============================================================
347 ! Read namelist main directories
348 !===============================================================
349 call set_sentinel(dir_meteo_header) ! set sentinal to check reading
350 inputformat_meteo_forcings = "nc"
351 bound_error = .true.
352 call position_nml(self%dir_nml_name, unamelist)
353 read(unamelist, nml = directories_mhm)
354
355 self%bound_error = bound_error
356 self%inputFormat_meteo_forcings = inputformat_meteo_forcings
357
358 do idomain = 1, self%nDomains
359 domainid = self%indices(idomain)
360 self%timestep_model_inputs(idomain) = time_step_model_inputs(domainid)
361 self%dirPrecipitation(idomain) = dir_precipitation(domainid)
362 self%dirTemperature(idomain) = dir_temperature(domainid)
363 self%dirReferenceET(idomain) = dir_referenceet(domainid)
364 self%dirMinTemperature(idomain) = dir_mintemperature(domainid)
365 self%dirMaxTemperature(idomain) = dir_maxtemperature(domainid)
366 self%dirNetRadiation(idomain) = dir_netradiation(domainid)
367 self%dirwindspeed(idomain) = dir_windspeed(domainid)
368 self%dirabsVapPressure(idomain) = dir_absvappressure(domainid)
369 ! riv-temp related
370 self%dirRadiation(idomain) = dir_radiation(domainid)
371 ! meteo header directory (if not given, use precipitation dir)
372 if (check_sentinel(dir_meteo_header(domainid))) then
373 self%dir_meteo_header(idomain) = self%dirPrecipitation(idomain)
374 else
375 self%dir_meteo_header(idomain) = dir_meteo_header(domainid)
376 end if
377 end do
378
379 ! consistency check for timestep_model_inputs
380 if (any(self%timestep_model_inputs .ne. 0) .and. .not. all(self%timestep_model_inputs .ne. 0)) then
381 call error_message('***ERROR: timestep_model_inputs either have to be all zero or all non-zero')
382 end if
383 ! check for optimzation and timestep_model_inputs options
384 if (optimize .and. (any(self%timestep_model_inputs .ne. 0))) then
385 call error_message('***ERROR: optimize and chunk read is switched on! (set timestep_model_inputs to zero)')
386 end if
387
388 !===============================================================
389 ! Read night-day ratios and pan evaporation
390 !===============================================================
391 ! default values for long/shortwave rad.
392 fnight_ssrd = 0.0_dp
393 fnight_strd = 0.45_dp
394
395 ! namelist for night-day ratio of precipitation, referenceET and temperature
396 call position_nml(self%weight_nml_name, unamelist)
397 read(unamelist, nml = nightdayratio)
398 self%read_meteo_weights = read_meteo_weights
399 self%fnight_prec = fnight_prec
400 self%fnight_pet = fnight_pet
401 self%fnight_temp = fnight_temp
402 self%fnight_ssrd = fnight_ssrd
403 self%fnight_strd = fnight_strd
404 self%fday_prec = 1.0_dp - self%fnight_prec
405 self%fday_pet = 1.0_dp - self%fnight_pet
406 self%fday_temp = -1.0_dp * self%fnight_temp
407 self%fday_ssrd = 1.0_dp - self%fnight_ssrd
408 self%fday_strd = 1.0_dp - self%fnight_strd
409
410 ! TODO-RIV-TEMP:
411 ! - add short- and long-wave raidiation weights (nc files)
412
413 ! closing the namelist file
414 call close_nml(unamelist)
415
416 ! check coupling configuration matching process cases
417 self%couple_is_hourly = .false.
418 self%couple_all = .false.
419 self%couple_pre = self%couple_cfg%active() .and. self%couple_cfg%meteo_expect_pre
420 self%couple_temp = self%couple_cfg%active() .and. self%couple_cfg%meteo_expect_temp
421 self%couple_pet = self%couple_cfg%active() .and. self%couple_cfg%meteo_expect_pet
422 self%couple_tmin = self%couple_cfg%active() .and. self%couple_cfg%meteo_expect_tmin
423 self%couple_tmax = self%couple_cfg%active() .and. self%couple_cfg%meteo_expect_tmax
424 self%couple_netrad = self%couple_cfg%active() .and. self%couple_cfg%meteo_expect_netrad
425 self%couple_absvappress = self%couple_cfg%active() .and. self%couple_cfg%meteo_expect_absvappress
426 self%couple_windspeed = self%couple_cfg%active() .and. self%couple_cfg%meteo_expect_windspeed
427 self%couple_ssrd = self%couple_cfg%active() .and. self%couple_cfg%meteo_expect_ssrd
428 self%couple_strd = self%couple_cfg%active() .and. self%couple_cfg%meteo_expect_strd
429 self%couple_tann = self%couple_cfg%active() .and. self%couple_cfg%meteo_expect_tann
430 if (self%couple_cfg%active()) then
431 self%couple_step_delta = timedelta(hours=self%couple_cfg%meteo_timestep)
432 self%couple_is_hourly = self%couple_step_delta == one_hour()
433 ! default init values for coupling times: 0001-01-01
434 if (self%couple_cfg%meteo_expect_pre) self%couple_pre_time = datetime()
435 if (self%couple_cfg%meteo_expect_temp) self%couple_temp_time = datetime()
436 ! PET releated
437 if (self%couple_cfg%meteo_expect_pet) self%couple_pet_time = datetime()
438 if (self%couple_cfg%meteo_expect_tmin) self%couple_tmin_time = datetime()
439 if (self%couple_cfg%meteo_expect_tmax) self%couple_tmax_time = datetime()
440 if (self%couple_cfg%meteo_expect_netrad) self%couple_netrad_time = datetime()
441 if (self%couple_cfg%meteo_expect_absvappress) self%couple_absvappress_time = datetime()
442 if (self%couple_cfg%meteo_expect_windspeed) self%couple_windspeed_time = datetime()
443 ! RIV-TEMP releated
444 if (self%couple_cfg%meteo_expect_ssrd) self%couple_ssrd_time = datetime()
445 if (self%couple_cfg%meteo_expect_strd) self%couple_strd_time = datetime()
446 if (self%couple_cfg%meteo_expect_tann) self%couple_tann_time = datetime()
447 ! PET related meteo
448 self%couple_all = self%couple_cfg%meteo_expect_pre .and. self%couple_cfg%meteo_expect_temp
449 select case (self%pet_case)
450 case(-1 : 0) ! pet is input
451 self%couple_all = self%couple_all .and. self%couple_cfg%meteo_expect_pet
452 if (self%couple_cfg%meteo_expect_tmin) call error_message("Coupling: tmin expected but not needed for PET.")
453 if (self%couple_cfg%meteo_expect_tmax) call error_message("Coupling: tmax expected but not needed for PET.")
454 if (self%couple_cfg%meteo_expect_netrad) call error_message("Coupling: netrad expected but not needed for PET.")
455 if (self%couple_cfg%meteo_expect_absvappress) call error_message("Coupling: absvappress expected but not needed for PET.")
456 if (self%couple_cfg%meteo_expect_windspeed) call error_message("Coupling: windspeed expected but not needed for PET.")
457
458 case(1) ! Hargreaves-Samani formulation (input: minimum and maximum Temperature)
459 self%couple_all = self%couple_all .and. self%couple_cfg%meteo_expect_tmin
460 self%couple_all = self%couple_all .and. self%couple_cfg%meteo_expect_tmax
461 if (self%couple_cfg%meteo_expect_pet) call error_message("Coupling: pet expected but not needed for PET.")
462 if (self%couple_cfg%meteo_expect_netrad) call error_message("Coupling: netrad expected but not needed for PET.")
463 if (self%couple_cfg%meteo_expect_absvappress) call error_message("Coupling: absvappress expected but not needed for PET.")
464 if (self%couple_cfg%meteo_expect_windspeed) call error_message("Coupling: windspeed expected but not needed for PET.")
465
466 case(2) ! Priestley-Taylor formulation (input: net radiation)
467 self%couple_all = self%couple_all .and. self%couple_cfg%meteo_expect_netrad
468 if (self%couple_cfg%meteo_expect_pet) call error_message("Coupling: pet expected but not needed for PET.")
469 if (self%couple_cfg%meteo_expect_tmin) call error_message("Coupling: tmin expected but not needed for PET.")
470 if (self%couple_cfg%meteo_expect_tmax) call error_message("Coupling: tmax expected but not needed for PET.")
471 if (self%couple_cfg%meteo_expect_absvappress) call error_message("Coupling: absvappress expected but not needed for PET.")
472 if (self%couple_cfg%meteo_expect_windspeed) call error_message("Coupling: windspeed expected but not needed for PET.")
473
474 case(3) ! Penman-Monteith formulation (input: net radiationm absulute vapour pressure, windspeed)
475 self%couple_all = self%couple_all .and. self%couple_cfg%meteo_expect_netrad
476 self%couple_all = self%couple_all .and. self%couple_cfg%meteo_expect_absvappress
477 self%couple_all = self%couple_all .and. self%couple_cfg%meteo_expect_windspeed
478 if (self%couple_cfg%meteo_expect_pet) call error_message("Coupling: pet expected but not needed for PET.")
479 if (self%couple_cfg%meteo_expect_tmin) call error_message("Coupling: tmin expected but not needed for PET.")
480 if (self%couple_cfg%meteo_expect_tmax) call error_message("Coupling: tmax expected but not needed for PET.")
481 end select
482 ! river temperature related meteo
483 if ( self%riv_temp_case == 0 ) then
484 if (self%couple_cfg%meteo_expect_ssrd) call error_message("Coupling: ssrd expected but river temperature not activated.")
485 if (self%couple_cfg%meteo_expect_strd) call error_message("Coupling: strd expected but river temperature not activated.")
486 if (self%couple_cfg%meteo_expect_tann) call error_message("Coupling: tann expected but river temperature not activated.")
487 else
488 self%couple_all = self%couple_all .and. self%couple_cfg%meteo_expect_ssrd
489 self%couple_all = self%couple_all .and. self%couple_cfg%meteo_expect_strd
490 self%couple_all = self%couple_all .and. self%couple_cfg%meteo_expect_tann
491 end if
492 end if
493 end subroutine config
494
495 !> \brief whether meteo data should be read completely at the begining
496 !> \return True if meteo data is retrieved with a single read
497 logical function single_read(self, iDomain)
498 implicit none
499 class(meteo_handler_type), intent(in) :: self
500 integer(i4), intent(in) :: idomain !< current domain
501 single_read = self%timeStep_model_inputs(idomain) == 0_i4
502 end function single_read
503
504 !> \brief Initialize meteo data and level-2 grid
505 subroutine init_level2(self, level0, level1)
506
507 use mo_grid, only : set_domain_indices
508 use mo_common_types, only : grid
509 use mo_grid, only : init_lowres_level
511 use mo_string_utils, only : num2str
512
513 implicit none
514
515 class(meteo_handler_type), intent(inout) :: self
516 !> grid information at level-0
517 type(grid), dimension(:), intent(in) :: level0
518 !> grid information at level-1 if all meteo data is coupled
519 type(grid), dimension(:), intent(in) :: level1
520
521 ! header info
522 integer(i4) :: nrows2, ncols2
523 real(dp) :: xllcorner2, yllcorner2
524 real(dp) :: cellsize2, nodata_dummy
525 ! file name
526 character(256) :: fName
527 ! looping variable
528 integer(i4) :: iDomain
529
530 ! create level-2 info
531 allocate(self%level2(self%nDomains))
532
533 ! we don't need level 2 if all meteo data comes from the coupler
534 if (self%couple_all) then
535 self%level2(:) = level1(:)
536 return
537 end if
538
539 do idomain = 1, self%nDomains
540 ! read header
541 fname = trim(adjustl(self%dir_meteo_header(idomain))) // trim(adjustl(self%file_meteo_header))
542 call read_header_ascii(trim(fname), self%umeteo_header, &
543 nrows2, ncols2, xllcorner2, yllcorner2, cellsize2, nodata_dummy)
544 ! check grid compatibility
545 call init_lowres_level(level0(self%L0DataFrom(idomain)), cellsize2, self%level2(idomain))
546 ! check
547 if ((ncols2 .ne. self%level2(idomain)%ncols) .or. &
548 (nrows2 .ne. self%level2(idomain)%nrows) .or. &
549 (abs(xllcorner2 - self%level2(idomain)%xllcorner) .gt. tiny(1.0_dp)) .or. &
550 (abs(yllcorner2 - self%level2(idomain)%yllcorner) .gt. tiny(1.0_dp)) .or. &
551 (abs(cellsize2 - self%level2(idomain)%cellsize) .gt. tiny(1.0_dp))) then
552 call error_message(' ***ERROR: subroutine L2_variable_init: size mismatch in grid file for level2 in domain ', &
553 trim(adjustl(num2str(idomain))), '!', raise=.false.)
554 call error_message(' Expected to have following properties (based on L0):', raise=.false.)
555 call error_message('... rows: ', trim(adjustl(num2str(self%level2(idomain)%nrows))), ', ', raise=.false.)
556 call error_message('... cols: ', trim(adjustl(num2str(self%level2(idomain)%ncols))), ', ', raise=.false.)
557 call error_message('... cellsize: ', trim(adjustl(num2str(self%level2(idomain)%cellsize))), ', ', raise=.false.)
558 call error_message('... xllcorner:', trim(adjustl(num2str(self%level2(idomain)%xllcorner))), ', ', raise=.false.)
559 call error_message('... yllcorner:', trim(adjustl(num2str(self%level2(idomain)%yllcorner))), ', ', raise=.false.)
560 call error_message(' Provided (in precipitation file):', raise=.false.)
561 call error_message('... rows: ', trim(adjustl(num2str(nrows2))), ', ', raise=.false.)
562 call error_message('... cols: ', trim(adjustl(num2str(ncols2))), ', ', raise=.false.)
563 call error_message('... cellsize: ', trim(adjustl(num2str(cellsize2))), ', ', raise=.false.)
564 call error_message('... xllcorner:', trim(adjustl(num2str(xllcorner2))), ', ', raise=.false.)
565 call error_message('... yllcorner:', trim(adjustl(num2str(yllcorner2))), ', ')
566 end if
567 end do
568
569 ! set indices
570 call set_domain_indices(self%level2)
571
572 end subroutine init_level2
573
574 !> \brief update the current time-step of the \ref meteo_handler_type class
575 subroutine update_timestep(self, tt, time, iDomain, level1, simPer)
576
577 implicit none
578
579 class(meteo_handler_type), intent(inout) :: self
580 integer(i4), intent(in) :: tt !< current time step
581 real(dp), intent(in) :: time !< current decimal Julian day
582 integer(i4), intent(in) :: iDomain !< current domain
583 !> grid information at hydrologic level
584 type(grid), dimension(:), intent(in) :: level1
585 !> warmPer + evalPer
586 type(period), dimension(:), intent(in) :: simPer
587
588 ! store current indizes and time
589 ! only needed for the "get_<var>" methods
590 self%s1 = level1(idomain)%iStart
591 self%e1 = level1(idomain)%iEnd
592 self%time = time
593
594 ! time increment is done right after call to mrm (and initially before looping)
595 if (self%single_read(idomain) .or. self%couple_all) then
596 ! whole meteorology is already read or all meteo is coupled
597
598 ! set start and end of meteo position
599 self%s_meteo = level1(idomain)%iStart
600 self%e_meteo = level1(idomain)%iEnd
601
602 ! time step for meteorological variable (daily values)
603 ! iMeteoTS = ceiling(real(tt, dp) / real(nTstepDay, dp))
604 if (self%couple_all) then
605 self%iMeteoTS = 1_i4
606 else
607 self%iMeteoTS = ceiling(real(tt, dp) / real(nint( 24._dp / real(self%nTstepForcingDay, dp)), dp))
608 end if
609 else
610 ! read chunk of meteorological forcings data (reading, upscaling/downscaling)
611 call self%prepare_data(tt, idomain, level1, simper)
612 ! set start and end of meteo position
613 self%s_meteo = 1
614 self%e_meteo = level1(idomain)%iEnd - level1(idomain)%iStart + 1
615 ! time step for meteorological variable (daily values)
616 self%iMeteoTS = ceiling(real(tt, dp) / real(nint( 24._dp / real(self%nTstepForcingDay, dp)), dp)) &
617 - (self%readPer%julStart - simper(idomain)%julStart)
618 end if
619
620 end subroutine update_timestep
621
622 !> \brief Prepare meteorological forcings data for a given variable
623 !> \details Prepare meteorological forcings data for a given variable.
624 !! Internally this subroutine calls another routine meteo_wrapper
625 !! for different meterological variables
626 !> \changelog
627 !! - Matthias Zink, Jun 2013
628 !! - addded NetCDf reader
629 !! - Rohini Kumar, Aug 2013
630 !! - name changed "inputFormat" to inputFormat_meteo_forcings
631 !! - Matthias Zink, Feb 2014
632 !! - added read in for different PET processes (process 5)
633 !! - Stephan Thober, Jun 2014
634 !! - add chunk_config for chunk read, copied L2 initialization to mo_startup
635 !! - Stephan Thober, Nov 2016
636 !! - moved processMatrix to common variables
637 !! - Stephan Thober, Jan 2017
638 !! - added subroutine for meteo_weights
639 !! - Robert Schweppe Jun 2018
640 !! - refactoring and reformatting
641 !! - Sebastian Müller Mar 2023
642 !! - converted routine to meteo-handler method
643 !> \authors Rohini Kumar
644 !> \date Jan 2013
645 subroutine prepare_data(self, tt, iDomain, level1, simPer)
646
647 use mo_string_utils, only : num2str
648 use mo_timer, only : timer_get, timer_start, timer_stop
650
651 implicit none
652
653 class(meteo_handler_type), intent(inout) :: self
654 integer(i4), intent(in) :: tt !< current timestep
655 integer(i4), intent(in) :: iDomain !< Domain number
656 !> grid information at hydrologic level
657 type(grid), dimension(:), intent(in) :: level1
658 !> warmPer + evalPer
659 type(period), dimension(:), intent(in) :: simPer
660
661 ! indicate whether data should be read
662 logical :: read_flag
663 integer(i4) :: domainID ! current domain ID
664
665 ! allocate arrays only once if they are coupled
666 if (self%couple_pre .and. .not. allocated(self%L1_pre)) allocate(self%L1_pre(level1(idomain)%nCells, 1))
667 if (self%couple_temp .and. .not. allocated(self%L1_temp)) allocate(self%L1_temp(level1(idomain)%nCells, 1))
668 if (self%couple_pet .and. .not. allocated(self%L1_pet)) allocate(self%L1_pet(level1(idomain)%nCells, 1))
669 if (self%couple_tmin .and. .not. allocated(self%L1_tmin)) allocate(self%L1_tmin(level1(idomain)%nCells, 1))
670 if (self%couple_tmax .and. .not. allocated(self%L1_tmax)) allocate(self%L1_tmax(level1(idomain)%nCells, 1))
671 if (self%couple_netrad .and. .not. allocated(self%L1_netrad)) allocate(self%L1_netrad(level1(idomain)%nCells, 1))
672 if (self%couple_absvappress .and. .not. allocated(self%L1_absvappress)) allocate(self%L1_absvappress(level1(idomain)%nCells, 1))
673 if (self%couple_windspeed .and. .not. allocated(self%L1_windspeed)) allocate(self%L1_windspeed(level1(idomain)%nCells, 1))
674 if (self%couple_ssrd .and. .not. allocated(self%L1_ssrd)) allocate(self%L1_ssrd(level1(idomain)%nCells, 1))
675 if (self%couple_strd .and. .not. allocated(self%L1_strd)) allocate(self%L1_strd(level1(idomain)%nCells, 1))
676 if (self%couple_tann .and. .not. allocated(self%L1_tann)) allocate(self%L1_tann(level1(idomain)%nCells, 1))
677
678 domainid = self%indices(idomain)
679
680 ! configuration of chunk_read
681 call chunk_config(idomain, tt, self%nTstepDay, simper, self%timestep, self%timeStep_model_inputs, read_flag, self%readPer)
682
683 ! only read, if read_flag is true
684 if (read_flag) then
685
686 ! read weights for hourly disaggregation
687 if (tt .eq. 1) then
688 ! TODO-RIV-TEMP: No NC files for weights for radiation at the moment
689 if (self%single_read(idomain)) call message(' read meteo weights for tavg ...')
690 call meteo_weights_wrapper(idomain, self%read_meteo_weights, self%dirTemperature(idomain), &
691 self%L1_temp_weights, level1=level1, level2=self%level2, ncvarname = 'tavg_weight')
692
693 if (self%single_read(idomain)) call message(' read meteo weights for pet ...')
694 call meteo_weights_wrapper(idomain, self%read_meteo_weights, self%dirReferenceET(idomain), &
695 self%L1_pet_weights, level1=level1, level2=self%level2, ncvarname = 'pet_weight')
696
697 if (self%single_read(idomain)) call message(' read meteo weights for pre ...')
698 call meteo_weights_wrapper(idomain, self%read_meteo_weights, self%dirPrecipitation(idomain), &
699 self%L1_pre_weights, level1=level1, level2=self%level2, ncvarname = 'pre_weight')
700 end if
701
702 ! free L1 variables if chunk read is activated
703 if (self%timeStep_model_inputs(idomain) .ne. 0) then
704 if (.not. self%couple_pre .and. allocated(self%L1_pre)) deallocate(self%L1_pre)
705 if (.not. self%couple_temp .and. allocated(self%L1_temp)) deallocate(self%L1_temp)
706 if (.not. self%couple_pet .and. allocated(self%L1_pet)) deallocate(self%L1_pet)
707 if (.not. self%couple_tmin .and. allocated(self%L1_tmin)) deallocate(self%L1_tmin)
708 if (.not. self%couple_tmax .and. allocated(self%L1_tmax)) deallocate(self%L1_tmax)
709 if (.not. self%couple_netrad .and. allocated(self%L1_netrad)) deallocate(self%L1_netrad)
710 if (.not. self%couple_absvappress .and. allocated(self%L1_absvappress)) deallocate(self%L1_absvappress)
711 if (.not. self%couple_windspeed .and. allocated(self%L1_windspeed)) deallocate(self%L1_windspeed)
712 end if
713
714 ! Domain characteristics and read meteo header
715 if (self%single_read(idomain)) then
716 call message(' Reading meteorological forcings for Domain: ', trim(adjustl(num2str(domainid))), ' ...')
717 call timer_start(1)
718 end if
719
720 ! precipitation
721 if (.not. self%couple_pre) then
722 if (self%single_read(idomain)) call message(' read precipitation ...')
723 ! upper bound: 1825 mm/d in La Réunion 7-8 Jan 1966
724 call meteo_forcings_wrapper(idomain, self%dirPrecipitation(idomain), self%inputFormat_meteo_forcings, &
725 dataout1=self%L1_pre, &
726 readper=self%readPer, ntstepforcingday=self%nTstepForcingDay, level1=level1, level2=self%level2, &
727 lower = 0.0_dp, upper = 2000._dp, ncvarname = 'pre', bound_error=self%bound_error)
728 end if
729
730 ! temperature
731 if (.not. self%couple_temp) then
732 if (self%single_read(idomain)) call message(' read temperature ...')
733 call meteo_forcings_wrapper(idomain, self%dirTemperature(idomain), self%inputFormat_meteo_forcings, &
734 dataout1=self%L1_temp, &
735 readper=self%readPer, ntstepforcingday=self%nTstepForcingDay, level1=level1, level2=self%level2, &
736 lower = -100._dp, upper = 100._dp, ncvarname = 'tavg', bound_error=self%bound_error)
737 end if
738
739 ! read input for PET (process 5) depending on specified option
740 ! 0 - input, 1 - Hargreaves-Samani, 2 - Priestley-Taylor, 3 - Penman-Monteith
741 select case (self%pet_case)
742 case(-1 : 0) ! pet is input
743 if (.not. self%couple_pet) then
744 if (self%single_read(idomain)) call message(' read pet ...')
745 call meteo_forcings_wrapper(idomain, self%dirReferenceET(idomain), self%inputFormat_meteo_forcings, &
746 dataout1=self%L1_pet, &
747 readper=self%readPer, ntstepforcingday=self%nTstepForcingDay, level1=level1, level2=self%level2, &
748 lower = 0.0_dp, upper = 1000._dp, ncvarname = 'pet', bound_error=self%bound_error)
749 end if
750 ! allocate PET and dummies for mhm_call
751 if ((idomain.eq.self%nDomains) .OR. (self%timeStep_model_inputs(idomain) .NE. 0)) then
752 allocate(self%L1_tmin(1, 1))
753 allocate(self%L1_tmax(1, 1))
754 allocate(self%L1_netrad(1, 1))
755 allocate(self%L1_absvappress(1, 1))
756 allocate(self%L1_windspeed(1, 1))
757 end if
758
759 case(1) ! Hargreaves-Samani formulation (input: minimum and maximum Temperature)
760 if (.not. self%couple_tmin) then
761 if (self%single_read(idomain)) call message(' read min. temperature ...')
762 call meteo_forcings_wrapper(idomain, self%dirMinTemperature(idomain), self%inputFormat_meteo_forcings, &
763 dataout1=self%L1_tmin, &
764 readper=self%readPer, ntstepforcingday=self%nTstepForcingDay, level1=level1, level2=self%level2, &
765 lower = -100.0_dp, upper = 100._dp, ncvarname = 'tmin', bound_error=self%bound_error)
766 end if
767 if (.not. self%couple_tmax) then
768 if (self%single_read(idomain)) call message(' read max. temperature ...')
769 call meteo_forcings_wrapper(idomain, self%dirMaxTemperature(idomain), self%inputFormat_meteo_forcings, &
770 dataout1=self%L1_tmax, &
771 readper=self%readPer, ntstepforcingday=self%nTstepForcingDay, level1=level1, level2=self%level2, &
772 lower = -100.0_dp, upper = 100._dp, ncvarname = 'tmax', bound_error=self%bound_error)
773 end if
774 ! allocate PET and dummies for mhm_call
775 if ((idomain .eq. self%nDomains) .OR. (self%timeStep_model_inputs(idomain) .NE. 0)) then
776 allocate(self%L1_pet (size(self%L1_tmax, dim = 1), size(self%L1_tmax, dim = 2)))
777 allocate(self%L1_netrad(1, 1))
778 allocate(self%L1_absvappress(1, 1))
779 allocate(self%L1_windspeed(1, 1))
780 end if
781
782 case(2) ! Priestley-Taylor formulation (input: net radiation)
783 if (.not. self%couple_netrad) then
784 if (self%single_read(idomain)) call message(' read net radiation ...')
785 call meteo_forcings_wrapper(idomain, self%dirNetRadiation(idomain), self%inputFormat_meteo_forcings, &
786 dataout1=self%L1_netrad, &
787 readper=self%readPer, ntstepforcingday=self%nTstepForcingDay, level1=level1, level2=self%level2, &
788 lower = -500.0_dp, upper = 1500._dp, ncvarname = 'net_rad', bound_error=self%bound_error)
789 end if
790 ! allocate PET and dummies for mhm_call
791 if ((idomain .eq. self%nDomains) .OR. (self%timeStep_model_inputs(idomain) .NE. 0)) then
792 allocate(self%L1_pet (size(self%L1_netrad, dim = 1), size(self%L1_netrad, dim = 2)))
793 allocate(self%L1_tmin(1, 1))
794 allocate(self%L1_tmax(1, 1))
795 allocate(self%L1_absvappress(1, 1))
796 allocate(self%L1_windspeed(1, 1))
797 end if
798
799 case(3) ! Penman-Monteith formulation (input: net radiationm absulute vapour pressure, windspeed)
800 if (.not. self%couple_netrad) then
801 if (self%single_read(idomain)) call message(' read net radiation ...')
802 call meteo_forcings_wrapper(idomain, self%dirNetRadiation(idomain), self%inputFormat_meteo_forcings, &
803 dataout1=self%L1_netrad, &
804 readper=self%readPer, ntstepforcingday=self%nTstepForcingDay, level1=level1, level2=self%level2, &
805 lower = -500.0_dp, upper = 1500._dp, ncvarname = 'net_rad', bound_error=self%bound_error)
806 end if
807 if (.not. self%couple_absvappress) then
808 if (self%single_read(idomain)) call message(' read absolute vapour pressure ...')
809 call meteo_forcings_wrapper(idomain, self%dirabsVapPressure(idomain), self%inputFormat_meteo_forcings, &
810 dataout1=self%L1_absvappress, &
811 readper=self%readPer, ntstepforcingday=self%nTstepForcingDay, level1=level1, level2=self%level2, &
812 lower = 0.0_dp, upper = 15000.0_dp, ncvarname = 'eabs', bound_error=self%bound_error)
813 end if
814 if (.not. self%couple_windspeed) then
815 if (self%single_read(idomain)) call message(' read windspeed ...')
816 call meteo_forcings_wrapper(idomain, self%dirwindspeed(idomain), self%inputFormat_meteo_forcings, &
817 dataout1=self%L1_windspeed, &
818 readper=self%readPer, ntstepforcingday=self%nTstepForcingDay, level1=level1, level2=self%level2, &
819 lower = 0.0_dp, upper = 250.0_dp, ncvarname = 'windspeed', bound_error=self%bound_error)
820 end if
821 ! allocate PET and dummies for mhm_call
822 if ((idomain.eq.self%nDomains) .OR. (self%timeStep_model_inputs(idomain) .NE. 0)) then
823 allocate(self%L1_pet (size(self%L1_absvappress, dim = 1), size(self%L1_absvappress, dim = 2)))
824 allocate(self%L1_tmin(1, 1))
825 allocate(self%L1_tmax(1, 1))
826 end if
827 end select
828
829 ! long/short-wave radiation and annual mean temperature for river-temperature routing
830 if ( self%riv_temp_case .ne. 0 ) then
831 ! free L1 variables if chunk read is activated
832 if (self%timeStep_model_inputs(idomain) .ne. 0) then
833 if (.not. self%couple_ssrd .and. allocated(self%L1_ssrd)) deallocate(self%L1_ssrd)
834 if (.not. self%couple_strd .and. allocated(self%L1_strd)) deallocate(self%L1_strd)
835 if (.not. self%couple_tann .and. allocated(self%L1_tann)) deallocate(self%L1_tann)
836 end if
837 if (.not. self%couple_ssrd) then
838 if (self%single_read(idomain)) call message(' read short-wave radiation ...')
840 idomain, self%dirRadiation(idomain), self%inputFormat_meteo_forcings, &
841 dataout1=self%L1_ssrd, &
842 readper=self%readPer, ntstepforcingday=self%nTstepForcingDay, level1=level1, level2=self%level2, &
843 lower = 0.0_dp, upper = 1500._dp, ncvarname = 'ssrd', bound_error=self%bound_error)
844 end if
845 if (.not. self%couple_strd) then
846 if (self%single_read(idomain)) call message(' read long-wave radiation ...')
848 idomain, self%dirRadiation(idomain), self%inputFormat_meteo_forcings, &
849 dataout1=self%L1_strd, &
850 readper=self%readPer, ntstepforcingday=self%nTstepForcingDay, level1=level1, level2=self%level2, &
851 lower = 0.0_dp, upper = 1500._dp, ncvarname = 'strd', bound_error=self%bound_error)
852 end if
853 if (.not. self%couple_tann) then
854 if (self%single_read(idomain)) call message(' read annual mean temperature ...')
856 idomain, self%dirTemperature(idomain), self%inputFormat_meteo_forcings, &
857 dataout1=self%L1_tann, &
858 readper=self%readPer, ntstepforcingday=self%nTstepForcingDay, level1=level1, level2=self%level2, &
859 lower = -100.0_dp, upper = 100._dp, ncvarname = 'tann', bound_error=self%bound_error)
860 end if
861 end if
862
863 if (self%single_read(idomain)) then
864 call timer_stop(1)
865 call message(' in ', trim(num2str(timer_get(1), '(F9.3)')), ' seconds.')
866 end if
867 end if
868
869 if (tt .eq. 1) then
870 ! set hourly flags (once at begining)
871 if (self%couple_all) then
872 self%nTstepForcingDay = int(one_day() / self%couple_step_delta, i4)
873 self%is_hourly_forcing = self%couple_is_hourly
874 else
875 self%is_hourly_forcing = (self%nTstepForcingDay .eq. 24_i4)
876 end if
877
878 ! check hourly data for PET (once at begining)
879 select case (self%pet_case)
880 case(1) ! Hargreaves-Samani formulation
881 if (self%couple_temp .and. self%couple_tmin .and. self%couple_tmax) then
882 if (self%couple_is_hourly) call error_message("Coupling: PET - Hargreaves-Samani needs daily data. Got hourly.")
883 else if (self%couple_temp .or. self%couple_tmin .or. self%couple_tmax) then
884 if (self%couple_is_hourly) call error_message("Coupling: PET - Hargreaves-Samani needs daily data. Got hourly.")
885 if (self%is_hourly_forcing) call error_message("Meteo: PET - Hargreaves-Samani needs daily data. Got hourly.")
886 else
887 if (self%is_hourly_forcing) call error_message("Meteo: PET - Hargreaves-Samani needs daily data. Got hourly.")
888 end if
889
890 case(2) ! Priestley-Taylor formulation
891 if (self%couple_temp .and. self%couple_netrad) then
892 if (self%couple_is_hourly) call error_message("Coupling: PET - Priestley-Taylor needs daily data. Got hourly.")
893 else if (self%couple_temp .or. self%couple_netrad) then
894 if (self%couple_is_hourly) call error_message("Coupling: PET - Priestley-Taylor needs daily data. Got hourly.")
895 if (self%is_hourly_forcing) call error_message("Meteo: PET - Priestley-Taylor needs daily data. Got hourly.")
896 else
897 if (self%is_hourly_forcing) call error_message("Meteo: PET - Priestley-Taylor needs daily data. Got hourly.")
898 end if
899
900 case(3) ! Penman-Monteith formulation
901 if (self%couple_temp .and. self%couple_netrad .and. self%couple_absvappress .and. self%couple_windspeed) then
902 if (self%couple_is_hourly) call error_message("Coupling: PET - Penman-Monteith needs daily data. Got hourly.")
903 else if (self%couple_temp .or. self%couple_netrad .or. self%couple_absvappress .or. self%couple_windspeed) then
904 if (self%couple_is_hourly) call error_message("Coupling: PET - Penman-Monteith needs daily data. Got hourly.")
905 if (self%is_hourly_forcing) call error_message("Meteo: PET - Penman-Monteith needs daily data. Got hourly.")
906 else
907 if (self%is_hourly_forcing) call error_message("Meteo: PET - Penman-Monteith needs daily data. Got hourly.")
908 end if
909 end select
910 end if
911
912 end subroutine prepare_data
913
914 !> \brief get corrected PET for the current timestep and domain
915 subroutine get_corrected_pet(self, pet_calc, &
916 petLAIcorFactorL1, fAsp, HarSamCoeff, latitude, PrieTayAlpha, aeroResist, surfResist)
917
918 use mo_mhm_constants, only : harsamconst
919 use mo_julian, only : date2dec, dec2date
921 use mo_string_utils, only : num2str
923
924 implicit none
925
926 class(meteo_handler_type), intent(inout) :: self
927 !> [mm TS-1] estimated PET (if PET is input = corrected values (fAsp*PET))
928 real(dp), dimension(:), intent(inout) :: pet_calc
929 !> PET correction factor based on LAI at level 1
930 real(dp), dimension(:), intent(in) :: petLAIcorFactorL1
931 !> [1] PET correction for Aspect at level 1
932 real(dp), dimension(:), intent(in) :: fAsp
933 !> [1] PET Hargreaves Samani coefficient at level 1
934 real(dp), dimension(:), intent(in) :: HarSamCoeff
935 !> latitude on level 1
936 real(dp), dimension(:), intent(in) :: latitude
937 !> [1] PET Priestley Taylor coefficient at level 1
938 real(dp), dimension(:), intent(in) :: PrieTayAlpha
939 !> [s m-1] PET aerodynamical resitance at level 1
940 real(dp), dimension(:), intent(in) :: aeroResist
941 !> [s m-1] PET bulk surface resitance at level 1
942 real(dp), dimension(:), intent(in) :: surfResist
943
944 ! pet in [mm d-1]
945 real(dp) :: pet
946 logical :: isday, is_hourly
947 integer(i4) :: year, month, day, hour
948 type(datetime) :: curr_dt
949 type(timedelta) :: meteo_time_delta
950
951 ! doy of the year [1-365 or 1-366]
952 integer(i4) :: doy
953
954 ! number of L1 cells
955 integer(i4) :: nCells1
956 ! cell index
957 integer(i4) :: k, i, s1
958 ! individual meteo time-steps for all variables due to possible coupling
959 integer(i4) :: mTS_pet, mTS_temp, mTS_tmin, mTS_tmax, mTS_rn, mTS_avp, mTS_ws
960
961 ! date and month of this timestep
962 call dec2date(self%time, yy = year, mm = month, dd = day, hh = hour)
963 curr_dt = datetime(year, month, day, hour)
964 doy = curr_dt%doy()
965
966 ncells1 = self%e1 - self%s1 + 1
967 s1 = self%s1
968
969 if (self%couple_pet) then
970 meteo_time_delta = curr_dt - self%couple_pet_time
971 ! check that the PET from the interface has the correct time-stamp
972 if (meteo_time_delta < zero_delta() .or. meteo_time_delta >= self%couple_step_delta) &
973 call error_message("meteo_handler: PET was expected from coupler, but has a wrong time-stamp.")
974 mts_pet = 1_i4
975 is_hourly = self%couple_is_hourly
976 else
977 mts_pet = self%iMeteoTS
978 is_hourly = self%is_hourly_forcing ! not needed to set with other variables
979 end if
980
981 if (self%couple_temp) then
982 meteo_time_delta = curr_dt - self%couple_temp_time
983 ! check that the temperature from the interface has the correct time-stamp
984 if (meteo_time_delta < zero_delta() .or. meteo_time_delta >= self%couple_step_delta) &
985 call error_message("meteo_handler: temperature was expected from coupler, but has a wrong time-stamp.")
986 mts_temp = 1_i4
987 is_hourly = .false.
988 else
989 mts_temp = self%iMeteoTS
990 end if
991
992 if (self%couple_tmin) then
993 meteo_time_delta = curr_dt - self%couple_tmin_time
994 ! check that the tmin from the interface has the correct time-stamp
995 if (meteo_time_delta < zero_delta() .or. meteo_time_delta >= self%couple_step_delta) &
996 call error_message("meteo_handler: minimum temperature was expected from coupler, but has a wrong time-stamp.")
997 mts_tmin = 1_i4
998 is_hourly = .false.
999 else
1000 mts_tmin = self%iMeteoTS
1001 end if
1002
1003 if (self%couple_tmax) then
1004 meteo_time_delta = curr_dt - self%couple_tmax_time
1005 ! check that the tmax from the interface has the correct time-stamp
1006 if (meteo_time_delta < zero_delta() .or. meteo_time_delta >= self%couple_step_delta) &
1007 call error_message("meteo_handler: maximum temperature was expected from coupler, but has a wrong time-stamp.")
1008 mts_tmax = 1_i4
1009 is_hourly = .false.
1010 else
1011 mts_tmax = self%iMeteoTS
1012 end if
1013
1014 if (self%couple_netrad) then
1015 meteo_time_delta = curr_dt - self%couple_netrad_time
1016 ! check that the netrad from the interface has the correct time-stamp
1017 if (meteo_time_delta < zero_delta() .or. meteo_time_delta >= self%couple_step_delta) &
1018 call error_message("meteo_handler: net-radiation was expected from coupler, but has a wrong time-stamp.")
1019 mts_rn = 1_i4
1020 is_hourly = .false.
1021 else
1022 mts_rn = self%iMeteoTS
1023 end if
1024
1025 if (self%couple_absvappress) then
1026 meteo_time_delta = curr_dt - self%couple_absvappress_time
1027 ! check that the absvappress from the interface has the correct time-stamp
1028 if (meteo_time_delta < zero_delta() .or. meteo_time_delta >= self%couple_step_delta) &
1029 call error_message("meteo_handler: Abs. vapour pressure was expected from coupler, but has a wrong time-stamp.")
1030 mts_avp = 1_i4
1031 is_hourly = .false.
1032 else
1033 mts_avp = self%iMeteoTS
1034 end if
1035
1036 if (self%couple_windspeed) then
1037 meteo_time_delta = curr_dt - self%couple_windspeed_time
1038 ! check that the windspeed from the interface has the correct time-stamp
1039 if (meteo_time_delta < zero_delta() .or. meteo_time_delta >= self%couple_step_delta) &
1040 call error_message("meteo_handler: windspeed was expected from coupler, but has a wrong time-stamp.")
1041 mts_ws = 1_i4
1042 is_hourly = .false.
1043 else
1044 mts_ws = self%iMeteoTS
1045 end if
1046
1047 ! flag for day or night depending on hours of the day
1048 isday = (hour .gt. 6) .AND. (hour .le. 18)
1049
1050 !$OMP parallel default(shared) &
1051 !$OMP private(k, pet, i)
1052 !$OMP do SCHEDULE(STATIC)
1053 do k = 1, ncells1
1054
1055 ! correct index on concatenated arrays
1056 i = self%s_meteo - 1 + k
1057
1058 ! PET calculation
1059 select case (self%pet_case)
1060 case(-1) ! PET is input ! correct pet for every day only once at the first time step
1061 pet = petlaicorfactorl1(k) * self%L1_pet(i, mts_pet)
1062
1063 case(0) ! PET is input ! correct pet for every day only once at the first time step
1064 pet = fasp(k) * self%L1_pet(i, mts_pet)
1065
1066 case(1) ! Hargreaves-Samani
1067 ! estimate day of the year (doy) for approximation of the extraterrestrial radiation
1068 if (self%L1_tmax(i, mts_tmax) .lt. self%L1_tmin(i, mts_tmin)) &
1069 call message('WARNING: tmax smaller than tmin at doy ', &
1070 num2str(doy), ' in year ', num2str(year), ' at cell', num2str(k), '!')
1071 pet = fasp(k) * pet_hargreaves( &
1072 harsamcoeff=harsamcoeff(k), &
1074 tavg=self%L1_temp(i, mts_temp), &
1075 tmax=self%L1_tmax(i, mts_tmax), &
1076 tmin=self%L1_tmin(i, mts_tmin), &
1077 latitude=latitude(k), &
1078 doy=doy)
1079
1080 case(2) ! Priestley-Taylor
1081 ! Priestley Taylor is not defined for values netrad < 0.0_dp
1082 pet = pet_priestly( &
1083 prietayparam=prietayalpha(k), &
1084 rn=max(self%L1_netrad(i, mts_rn), 0.0_dp), &
1085 tavg=self%L1_temp(i, mts_temp))
1086
1087 case(3) ! Penman-Monteith
1088 pet = pet_penman( &
1089 net_rad=max(self%L1_netrad(i, mts_rn), 0.0_dp), &
1090 tavg=self%L1_temp(i, mts_temp), &
1091 act_vap_pressure=self%L1_absvappress(i, mts_avp) / 1000.0_dp, &
1092 aerodyn_resistance=aeroresist(k) / self%L1_windspeed(i, mts_ws), &
1093 bulksurface_resistance=surfresist(k), &
1094 a_s=1.0_dp, &
1095 a_sh=1.0_dp)
1096 end select
1097
1098 ! temporal disaggreagtion of forcing variables
1099 if (self%is_hourly_forcing) then
1100 pet_calc(k) = pet
1101 else
1102 if (self%read_meteo_weights) then
1103 ! all meteo forcings are disaggregated with given weights
1105 meteo_val_day=pet, &
1106 meteo_val_weights=self%L1_pet_weights(s1 - 1 + k, month, hour + 1), &
1107 meteo_val=pet_calc(k))
1108 else
1109 ! all meteo forcings are disaggregated with day-night correction values
1111 isday=isday, &
1112 ntimesteps_day=self%nTstepDay_dp, &
1113 meteo_val_day=pet, &
1114 fday_meteo_val=self%fday_pet(month), &
1115 fnight_meteo_val=self%fnight_pet(month), &
1116 meteo_val=pet_calc(k))
1117 end if
1118 end if
1119 end do
1120 !$OMP end do
1121 !$OMP end parallel
1122
1123 end subroutine get_corrected_pet
1124
1125 !> \brief get surface temperature for the current timestep and domain
1126 subroutine get_temp(self, temp_calc)
1127
1128 use mo_julian, only : dec2date
1130 use mo_constants, only : t0_dp ! 273.15 - Celcius <-> Kelvin [K]
1131
1132 implicit none
1133
1134 class(meteo_handler_type), intent(inout) :: self
1135 !> [degC] temperature for current time step
1136 real(dp), dimension(:), intent(inout) :: temp_calc
1137
1138 logical :: isday, is_hourly
1139 integer(i4) :: year, month, day, hour
1140 type(datetime) :: curr_dt
1141 type(timedelta) :: meteo_time_delta
1142
1143 ! number of L1 cells
1144 integer(i4) :: nCells1
1145 ! cell index
1146 integer(i4) :: k, i, s1, mTS
1147
1148 ! date and month of this timestep
1149 call dec2date(self%time, yy=year, mm=month, dd=day, hh=hour)
1150
1151 ncells1 = self%e1 - self%s1 + 1
1152 s1 = self%s1
1153 if (self%couple_temp) then
1154 curr_dt = datetime(year, month, day, hour)
1155 meteo_time_delta = curr_dt - self%couple_temp_time
1156 ! check that the temperature from the interface has the correct time-stamp
1157 if (meteo_time_delta < zero_delta() .or. meteo_time_delta >= self%couple_step_delta) &
1158 call error_message("meteo_handler: temperature was expected from coupler, but has a wrong time-stamp.")
1159 mts = 1_i4
1160 is_hourly = self%couple_is_hourly
1161 else
1162 mts = self%iMeteoTS
1163 is_hourly = self%is_hourly_forcing
1164 end if
1165
1166 ! shortcut hourly data
1167 if (is_hourly) then
1168 temp_calc(:) = self%L1_temp(self%s_meteo : self%e_meteo, mts)
1169 return
1170 end if
1171
1172 ! flag for day or night depending on hours of the day
1173 isday = (hour .gt. 6) .AND. (hour .le. 18)
1174
1175 !$OMP parallel default(shared) &
1176 !$OMP private(k, i)
1177 !$OMP do SCHEDULE(STATIC)
1178 do k = 1, ncells1
1179
1180 ! correct index on concatenated arrays
1181 i = self%s_meteo - 1 + k
1182
1183 if (self%read_meteo_weights) then
1184 ! all meteo forcings are disaggregated with given weights
1186 meteo_val_day=self%L1_temp(i, mts), &
1187 meteo_val_weights=self%L1_temp_weights(s1 - 1 + k, month, hour + 1), &
1188 meteo_val=temp_calc(k), &
1189 weights_correction=t0_dp)
1190 else
1191 ! all meteo forcings are disaggregated with day-night correction values
1193 isday=isday, &
1194 ntimesteps_day=self%nTstepDay_dp, &
1195 meteo_val_day=self%L1_temp(i, mts), &
1196 fday_meteo_val=self%fday_temp(month), &
1197 fnight_meteo_val=self%fnight_temp(month), &
1198 meteo_val=temp_calc(k), &
1199 add_correction=.true.)
1200 end if
1201 end do
1202 !$OMP end do
1203 !$OMP end parallel
1204
1205 end subroutine get_temp
1206
1207 !> \brief get precipitation for the current timestep and domain
1208 subroutine get_prec(self, prec_calc)
1209
1210 use mo_julian, only : dec2date
1212
1213 implicit none
1214
1215 class(meteo_handler_type), intent(inout) :: self
1216 !> [mm TS-1] precipitation for current time step
1217 real(dp), dimension(:), intent(inout) :: prec_calc
1218
1219 logical :: isday, is_hourly
1220 integer(i4) :: year, month, day, hour
1221 type(datetime) :: curr_dt
1222 type(timedelta) :: meteo_time_delta
1223
1224 ! number of L1 cells
1225 integer(i4) :: nCells1
1226 ! cell index
1227 integer(i4) :: k, i, s1, mTS
1228
1229 ! date and month of this timestep
1230 call dec2date(self%time, yy=year, mm=month, dd=day, hh=hour)
1231
1232 ncells1 = self%e1 - self%s1 + 1
1233 s1 = self%s1
1234 if (self%couple_pre) then
1235 curr_dt = datetime(year, month, day, hour)
1236 meteo_time_delta = curr_dt - self%couple_pre_time
1237 ! check that the precipitation from the interface has the correct time-stamp
1238 if (meteo_time_delta < zero_delta() .or. meteo_time_delta >= self%couple_step_delta) &
1239 call error_message("meteo_handler: precipitation was expected from coupler, but has a wrong time-stamp.")
1240 mts = 1_i4
1241 is_hourly = self%couple_is_hourly
1242 else
1243 mts = self%iMeteoTS
1244 is_hourly = self%is_hourly_forcing
1245 end if
1246
1247 ! shortcut hourly data
1248 if (is_hourly) then
1249 prec_calc(:) = self%L1_pre(self%s_meteo : self%e_meteo, mts)
1250 return
1251 end if
1252
1253 ! flag for day or night depending on hours of the day
1254 isday = (hour .gt. 6) .AND. (hour .le. 18)
1255
1256 !$OMP parallel default(shared) &
1257 !$OMP private(k, i)
1258 !$OMP do SCHEDULE(STATIC)
1259 do k = 1, ncells1
1260
1261 ! correct index on concatenated arrays
1262 i = self%s_meteo - 1 + k
1263
1264 ! temporal disaggreagtion of forcing variables
1265 if (self%read_meteo_weights) then
1266 ! all meteo forcings are disaggregated with given weights
1268 meteo_val_day=self%L1_pre(i, mts), &
1269 meteo_val_weights=self%L1_pre_weights(s1 - 1 + k, month, hour + 1), &
1270 meteo_val=prec_calc(k))
1271 else
1272 ! all meteo forcings are disaggregated with day-night correction values
1274 isday=isday, &
1275 ntimesteps_day=self%nTstepDay_dp, &
1276 meteo_val_day=self%L1_pre(i, mts), &
1277 fday_meteo_val=self%fday_prec(month), &
1278 fnight_meteo_val=self%fnight_prec(month), &
1279 meteo_val=prec_calc(k))
1280 end if
1281 end do
1282 !$OMP end do
1283 !$OMP end parallel
1284
1285 end subroutine get_prec
1286
1287 !> \brief get surface short-wave (solar) radiation downwards for the current timestep and domain
1288 subroutine get_ssrd(self, ssrd_calc)
1289
1290 use mo_julian, only : dec2date
1292
1293 implicit none
1294
1295 class(meteo_handler_type), intent(inout) :: self
1296 !> [W m2] surface short-wave (solar) radiation downwards for current time step
1297 real(dp), dimension(:), intent(inout) :: ssrd_calc
1298
1299 logical :: isday, is_hourly
1300 integer(i4) :: year, month, day, hour
1301 type(datetime) :: curr_dt
1302 type(timedelta) :: meteo_time_delta
1303
1304 ! number of L1 cells
1305 integer(i4) :: nCells1
1306 ! cell index
1307 integer(i4) :: k, i, s1, mTS
1308
1309 ! date and month of this timestep
1310 call dec2date(self%time, yy=year, mm=month, dd=day, hh=hour)
1311
1312 ncells1 = self%e1 - self%s1 + 1
1313 s1 = self%s1
1314 if (self%couple_ssrd) then
1315 curr_dt = datetime(year, month, day, hour)
1316 meteo_time_delta = curr_dt - self%couple_ssrd_time
1317 ! check that the ssrd from the interface has the correct time-stamp
1318 if (meteo_time_delta < zero_delta() .or. meteo_time_delta >= self%couple_step_delta) &
1319 call error_message("meteo_handler: ssrd was expected from coupler, but has a wrong time-stamp.")
1320 mts = 1_i4
1321 is_hourly = self%couple_is_hourly
1322 else
1323 mts = self%iMeteoTS
1324 is_hourly = self%is_hourly_forcing
1325 end if
1326
1327 ! shortcut hourly data
1328 if (is_hourly) then
1329 ssrd_calc(:) = self%L1_ssrd(self%s_meteo : self%e_meteo, mts)
1330 return
1331 end if
1332
1333 ! flag for day or night depending on hours of the day
1334 isday = (hour .gt. 6) .AND. (hour .le. 18)
1335
1336 !$OMP parallel default(shared) &
1337 !$OMP private(k, i)
1338 !$OMP do SCHEDULE(STATIC)
1339 do k = 1, ncells1
1340 ! correct index on concatenated arrays
1341 i = self%s_meteo - 1 + k
1342 ! TODO-RIV-TEMP: add weights for ssrd
1343 ! temporal disaggreagtion of forcing variables
1345 isday=isday, &
1346 ntimesteps_day=self%nTstepDay_dp, &
1347 meteo_val_day=self%L1_ssrd(i, mts), &
1348 fday_meteo_val=self%fday_ssrd(month), &
1349 fnight_meteo_val=self%fnight_ssrd(month), &
1350 meteo_val=ssrd_calc(k) &
1351 )
1352 end do
1353 !$OMP end do
1354 !$OMP end parallel
1355
1356 end subroutine get_ssrd
1357
1358 !> \brief get surface long-wave (thermal) radiation downwards for the current timestep and domain
1359 subroutine get_strd(self, strd_calc)
1360
1361 use mo_julian, only : dec2date
1363
1364 implicit none
1365
1366 class(meteo_handler_type), intent(inout) :: self
1367 !> [W m2] surface long-wave (thermal) radiation downwards for current time step
1368 real(dp), dimension(:), intent(inout) :: strd_calc
1369
1370 logical :: isday, is_hourly
1371 integer(i4) :: year, month, day, hour
1372 type(datetime) :: curr_dt
1373 type(timedelta) :: meteo_time_delta
1374
1375 ! number of L1 cells
1376 integer(i4) :: nCells1
1377 ! cell index
1378 integer(i4) :: k, i, s1, mTS
1379
1380 ! date and month of this timestep
1381 call dec2date(self%time, yy=year, mm=month, dd=day, hh=hour)
1382
1383 ncells1 = self%e1 - self%s1 + 1
1384 s1 = self%s1
1385 if (self%couple_strd) then
1386 curr_dt = datetime(year, month, day, hour)
1387 meteo_time_delta = curr_dt - self%couple_strd_time
1388 ! check that the strd from the interface has the correct time-stamp
1389 if (meteo_time_delta < zero_delta() .or. meteo_time_delta >= self%couple_step_delta) &
1390 call error_message("meteo_handler: strd was expected from coupler, but has a wrong time-stamp.")
1391 mts = 1_i4
1392 is_hourly = self%couple_is_hourly
1393 else
1394 mts = self%iMeteoTS
1395 is_hourly = self%is_hourly_forcing
1396 end if
1397
1398 ! shortcut hourly data
1399 if (is_hourly) then
1400 strd_calc(:) = self%L1_strd(self%s_meteo : self%e_meteo, mts)
1401 return
1402 end if
1403
1404 ! flag for day or night depending on hours of the day
1405 isday = (hour .gt. 6) .AND. (hour .le. 18)
1406
1407 !$OMP parallel default(shared) &
1408 !$OMP private(k, i)
1409 !$OMP do SCHEDULE(STATIC)
1410 do k = 1, ncells1
1411 ! correct index on concatenated arrays
1412 i = self%s_meteo - 1 + k
1413 ! TODO-RIV-TEMP: add weights for strd
1414 ! temporal disaggreagtion of forcing variables
1416 isday=isday, &
1417 ntimesteps_day=self%nTstepDay_dp, &
1418 meteo_val_day=self%L1_strd(i, mts), &
1419 fday_meteo_val=self%fday_strd(month), &
1420 fnight_meteo_val=self%fnight_strd(month), &
1421 meteo_val=strd_calc(k) &
1422 )
1423 end do
1424 !$OMP end do
1425 !$OMP end parallel
1426
1427 end subroutine get_strd
1428
1429 !> \brief get annual mean surface temperature for the current timestep and domain
1430 subroutine get_tann(self, tann_calc)
1431 use mo_julian, only : dec2date
1432
1433 implicit none
1434
1435 class(meteo_handler_type), intent(inout) :: self
1436 !> [degC] annual mean air temperature
1437 real(dp), dimension(:), intent(inout) :: tann_calc
1438
1439 type(datetime) :: curr_dt
1440 type(timedelta) :: meteo_time_delta
1441 integer(i4) :: year, month, day, hour, mTS
1442
1443 if (self%couple_tann) then
1444 ! date and month of this timestep
1445 call dec2date(self%time, yy=year, mm=month, dd=day, hh=hour)
1446 curr_dt = datetime(year, month, day, hour)
1447 meteo_time_delta = curr_dt - self%couple_tann_time
1448 ! check that the tann from the interface has the correct time-stamp
1449 if (meteo_time_delta < zero_delta() .or. meteo_time_delta >= self%couple_step_delta) &
1450 call error_message("meteo_handler: tann was expected from coupler, but has a wrong time-stamp.")
1451 mts = 1_i4
1452 else
1453 mts = self%iMeteoTS
1454 end if
1455
1456 ! annual temperature is not disaggregated
1457 tann_calc(:) = self%L1_tann(self%s_meteo : self%e_meteo, mts)
1458
1459 end subroutine get_tann
1460
1461 !> \brief set meteo_data from coupling
1462 subroutine set_meteo(self, &
1463 year, month, day, hour, &
1464 pre, &
1465 temp, &
1466 pet, &
1467 tmin, &
1468 tmax, &
1469 netrad, &
1470 absvappress, &
1471 windspeed, &
1472 ssrd, &
1473 strd, &
1474 tann)
1475 implicit none
1476 class(meteo_handler_type), intent(inout) :: self
1477 integer(i4), intent(in) :: year
1478 integer(i4), intent(in) :: month
1479 integer(i4), intent(in) :: day
1480 integer(i4), intent(in), optional :: hour
1481 real(dp), dimension(:), optional, intent(in) :: pre
1482 real(dp), dimension(:), optional, intent(in) :: temp
1483 real(dp), dimension(:), optional, intent(in) :: pet
1484 real(dp), dimension(:), optional, intent(in) :: tmin
1485 real(dp), dimension(:), optional, intent(in) :: tmax
1486 real(dp), dimension(:), optional, intent(in) :: netrad
1487 real(dp), dimension(:), optional, intent(in) :: absvappress
1488 real(dp), dimension(:), optional, intent(in) :: windspeed
1489 real(dp), dimension(:), optional, intent(in) :: ssrd
1490 real(dp), dimension(:), optional, intent(in) :: strd
1491 real(dp), dimension(:), optional, intent(in) :: tann
1492
1493 integer(i4) :: hour_
1494 type(datetime) :: input_time
1495
1496 if (.not. self%couple_cfg%active()) &
1497 call error_message("meteo_handler%set_meteo: coupling was not activated.")
1498
1499 ! determine input time
1500 hour_ = -1_i4
1501 if (present(hour)) then
1502 hour_ = hour
1503 else if (self%couple_cfg%meteo_timestep == 24_i4) then
1504 hour_ = 0_i4
1505 end if
1506 if (hour_ == -1_i4) &
1507 call error_message("meteo_handler%set_meteo: hour for the meteo date needs to be given if the timestep is not daily.")
1508 input_time = datetime(year, month, day, hour_)
1509
1510 ! fix input time, if reference point is at the end of the time interval
1511 if (self%couple_cfg%meteo_time_ref_endpoint) input_time = input_time - self%couple_step_delta
1512
1513 ! check if input time matches the required time step
1514 if (mod(input_time%hour, self%couple_cfg%meteo_timestep) /= 0) &
1515 call error_message("meteo_handler%set_meteo: given time doesn't match couple timestep: ", input_time%str())
1516
1517 ! precipitation
1518 if (present(pre)) then
1519 if (.not. self%couple_pre) &
1520 call error_message("meteo_handler%set_meteo: precipitation was not set to be coupled.")
1521 self%couple_pre_time = input_time
1522 self%L1_pre(:, 1_i4) = pre(:)
1523 end if
1524
1525 ! temperature
1526 if (present(temp)) then
1527 if (.not. self%couple_temp) &
1528 call error_message("meteo_handler%set_meteo: avg. temperature was not set to be coupled.")
1529 self%couple_temp_time = input_time
1530 self%L1_temp(:, 1_i4) = temp(:)
1531 end if
1532
1533 ! PET
1534 if (present(pet)) then
1535 if (.not. self%couple_pet) &
1536 call error_message("meteo_handler%set_meteo: PET was not set to be coupled.")
1537 self%couple_pet_time = input_time
1538 self%L1_pet(:, 1_i4) = pet(:)
1539 end if
1540
1541 ! tmin
1542 if (present(tmin)) then
1543 if (.not. self%couple_tmin) &
1544 call error_message("meteo_handler%set_meteo: tmin was not set to be coupled.")
1545 self%couple_tmin_time = input_time
1546 self%L1_tmin(:, 1_i4) = tmin(:)
1547 end if
1548
1549 ! tmax
1550 if (present(tmax)) then
1551 if (.not. self%couple_tmax) &
1552 call error_message("meteo_handler%set_meteo: tmax was not set to be coupled.")
1553 self%couple_tmax_time = input_time
1554 self%L1_tmax(:, 1_i4) = tmax(:)
1555 end if
1556
1557 ! netrad
1558 if (present(netrad)) then
1559 if (.not. self%couple_netrad) &
1560 call error_message("meteo_handler%set_meteo: netrad was not set to be coupled.")
1561 self%couple_netrad_time = input_time
1562 self%L1_netrad(:, 1_i4) = netrad(:)
1563 end if
1564
1565 ! absvappress
1566 if (present(absvappress)) then
1567 if (.not. self%couple_absvappress) &
1568 call error_message("meteo_handler%set_meteo: absvappress was not set to be coupled.")
1569 self%couple_absvappress_time = input_time
1570 self%L1_absvappress(:, 1_i4) = absvappress(:)
1571 end if
1572
1573 ! windspeed
1574 if (present(windspeed)) then
1575 if (.not. self%couple_windspeed) &
1576 call error_message("meteo_handler%set_meteo: windspeed was not set to be coupled.")
1577 self%couple_windspeed_time = input_time
1578 self%L1_windspeed(:, 1_i4) = windspeed(:)
1579 end if
1580
1581 ! ssrd
1582 if (present(ssrd)) then
1583 if (.not. self%couple_ssrd) &
1584 call error_message("meteo_handler%set_meteo: ssrd was not set to be coupled.")
1585 self%couple_ssrd_time = input_time
1586 self%L1_ssrd(:, 1_i4) = ssrd(:)
1587 end if
1588
1589 ! strd
1590 if (present(strd)) then
1591 if (.not. self%couple_strd) &
1592 call error_message("meteo_handler%set_meteo: strd was not set to be coupled.")
1593 self%couple_strd_time = input_time
1594 self%L1_strd(:, 1_i4) = strd(:)
1595 end if
1596
1597 ! tann
1598 if (present(tann)) then
1599 if (.not. self%couple_tann) &
1600 call error_message("meteo_handler%set_meteo: tann was not set to be coupled.")
1601 self%couple_tann_time = input_time
1602 self%L1_tann(:, 1_i4) = tann(:)
1603 end if
1604
1605 end subroutine set_meteo
1606
1607end module mo_meteo_handler
Provides constants commonly used by mHM, mRM and MPR.
integer(i4), parameter, public maxnodomains
integer(i4), parameter, public nodata_i4
Provides common types needed by mHM, mRM and/or mpr.
Provides structures needed by mHM, mRM and/or mpr.
integer(i4), parameter, public nprocesses
Types to specify the coupling configuration of mHM.
subroutine clean_up(self)
clean up
gridding tools
Definition mo_grid.f90:12
subroutine, public set_domain_indices(grids, indices)
TODO: add description.
Definition mo_grid.f90:205
subroutine, public init_lowres_level(highres, target_resolution, lowres, highres_lowres_remap)
Level-1 variable initialization.
Definition mo_grid.f90:59
Class for the meteo handler.
subroutine set_meteo(self, year, month, day, hour, pre, temp, pet, tmin, tmax, netrad, absvappress, windspeed, ssrd, strd, tann)
set meteo_data from coupling
subroutine get_temp(self, temp_calc)
get surface temperature for the current timestep and domain
subroutine update_timestep(self, tt, time, idomain, level1, simper)
update the current time-step of the meteo_handler_type class
subroutine config(self, file_namelist, unamelist, optimize, domainmeta, processmatrix, timestep, couple_cfg)
configure the meteo_handler_type class from the mhm namelist
subroutine get_ssrd(self, ssrd_calc)
get surface short-wave (solar) radiation downwards for the current timestep and domain
subroutine get_corrected_pet(self, pet_calc, petlaicorfactorl1, fasp, harsamcoeff, latitude, prietayalpha, aeroresist, surfresist)
get corrected PET for the current timestep and domain
subroutine get_tann(self, tann_calc)
get annual mean surface temperature for the current timestep and domain
subroutine get_prec(self, prec_calc)
get precipitation for the current timestep and domain
subroutine init_level2(self, level0, level1)
Initialize meteo data and level-2 grid.
logical function single_read(self, idomain)
whether meteo data should be read completely at the begining
subroutine prepare_data(self, tt, idomain, level1, simper)
Prepare meteorological forcings data for a given variable.
subroutine get_strd(self, strd_calc)
get surface long-wave (thermal) radiation downwards for the current timestep and domain
Prepare meteorological forcings data for mHM.
subroutine, public chunk_config(idomain, tt, ntstepday, simper, timestep, timestep_model_inputs, read_flag, readper)
determines the start date, end date, and read_flag given Domain id and current timestep
subroutine, public meteo_weights_wrapper(idomain, read_meteo_weights, datapath, dataout1, level1, level2, lower, upper, ncvarname)
Prepare weights for meteorological forcings data for mHM at Level-1.
subroutine, public meteo_forcings_wrapper(idomain, datapath, inputformat, dataout1, readper, ntstepforcingday, level1, level2, lower, upper, ncvarname, bound_error)
Prepare meteorological forcings data for mHM at Level-1.
Temporal disaggregation of daily input values.
elemental subroutine, public temporal_disagg_state_daynight(isday, ntimesteps_day, meteo_val_day, fday_meteo_val, fnight_meteo_val, meteo_val, add_correction)
Temporally distribute daily mean state forcings onto time step.
elemental subroutine, public temporal_disagg_meteo_weights(meteo_val_day, meteo_val_weights, meteo_val, weights_correction)
Temporally distribute daily mean forcings onto time step.
elemental subroutine, public temporal_disagg_flux_daynight(isday, ntimesteps_day, meteo_val_day, fday_meteo_val, fnight_meteo_val, meteo_val)
Temporally distribute daily mean forcings onto time step.
Provides mHM specific constants.
real(dp), parameter, public harsamconst
Hargreaves-Samani ref.
Module for calculating reference/potential evapotranspiration [mm d-1].
Definition mo_pet.f90:15
elemental pure real(dp) function, public pet_penman(net_rad, tavg, act_vap_pressure, aerodyn_resistance, bulksurface_resistance, a_s, a_sh)
Reference Evapotranspiration after Penman-Monteith.
Definition mo_pet.f90:237
elemental pure real(dp) function, public pet_hargreaves(harsamcoeff, harsamconst, tavg, tmax, tmin, latitude, doy)
Reference Evapotranspiration after Hargreaves.
Definition mo_pet.f90:74
elemental pure real(dp) function, public pet_priestly(prietayparam, rn, tavg)
Reference Evapotranspiration after Priestly-Taylor.
Definition mo_pet.f90:150
Reads spatial input data.
subroutine, public read_header_ascii(filename, fileunit, header_ncols, header_nrows, header_xllcorner, header_yllcorner, header_cellsize, header_nodata)
Reads header lines of ASCII files.
DOMAIN general description.
This is a container to hold all coupling configurations for mHM.
This is a handler for the meteorological forcings.