5.13.2-dev0
mHM
The mesoscale Hydrological Model
Loading...
Searching...
No Matches
mo_meteo_handler Module Reference

Class for the meteo handler. More...

Data Types

module  meteo_handler_type
 This is a handler for the meteorological forcings. More...
 

Functions/Subroutines

subroutine clean_up (self)
 clean up
 
subroutine config (self, file_namelist, unamelist, optimize, domainmeta, processmatrix, timestep, couple_cfg)
 configure the meteo_handler_type class from the mhm namelist
 
logical function single_read (self, idomain)
 whether meteo data should be read completely at the begining
 
subroutine init_level2 (self, level0, level1)
 Initialize meteo data and level-2 grid.
 
subroutine update_timestep (self, tt, time, idomain, level1, simper)
 update the current time-step of the meteo_handler_type class
 
subroutine prepare_data (self, tt, idomain, level1, simper)
 Prepare meteorological forcings data for a given variable.
 
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_temp (self, temp_calc)
 get surface temperature for the current timestep and domain
 
subroutine get_prec (self, prec_calc)
 get precipitation for the current timestep and domain
 
subroutine get_ssrd (self, ssrd_calc)
 get surface short-wave (solar) radiation downwards for the current timestep and domain
 
subroutine get_strd (self, strd_calc)
 get surface long-wave (thermal) radiation downwards for the current timestep and domain
 
subroutine get_tann (self, tann_calc)
 get annual mean surface temperature for the current timestep and domain
 
subroutine set_meteo (self, year, month, day, hour, pre, temp, pet, tmin, tmax, netrad, absvappress, windspeed, ssrd, strd, tann)
 set meteo_data from coupling
 

Detailed Description

Class for the meteo handler.

Handler for meteorological forcings in mHM. Is independent of global variables and provides 3 methods to access forcings:

  • get_corrected_pet : get the modified pet for mHM for the current timestep
  • get_temp : get the temporal disaggregated temperature for the current timestep
  • get_prec : get the temporal disaggregated percipitation for the current timestep
Version
0.1
Authors
Sebastian Mueller
Date
Mar 2023

Function/Subroutine Documentation

◆ clean_up()

subroutine mo_meteo_handler::clean_up ( class(meteo_handler_type), intent(inout)  self)

clean up

Definition at line 206 of file mo_meteo_handler.f90.

◆ config()

subroutine mo_meteo_handler::config ( class(meteo_handler_type), intent(inout)  self,
character(*), intent(in)  file_namelist,
integer, intent(in)  unamelist,
logical, intent(in)  optimize,
type(domain_meta), intent(in)  domainmeta,
integer(i4), dimension(nprocesses, 3), intent(in)  processmatrix,
integer(i4), intent(in)  timestep,
type(couple_cfg_type), intent(in)  couple_cfg 
)
private

configure the meteo_handler_type class from the mhm namelist

Parameters
[in]file_namelistmhm namelist file
[in]unamelistunit to open namelist file
[in]optimizeOptimization flag
[in]domainmetadomain general description
[in]processmatrixInfo about which process runs in which option
[in]timestep[h] simulation time step (= TS) in [h]
[in]couple_cfgcoupling configuration class

Definition at line 243 of file mo_meteo_handler.f90.

References mo_common_constants::maxnodomains, mo_common_constants::nodata_i4, and mo_common_variables::nprocesses.

◆ get_corrected_pet()

subroutine mo_meteo_handler::get_corrected_pet ( class(meteo_handler_type), intent(inout)  self,
real(dp), dimension(:), intent(inout)  pet_calc,
real(dp), dimension(:), intent(in)  petlaicorfactorl1,
real(dp), dimension(:), intent(in)  fasp,
real(dp), dimension(:), intent(in)  harsamcoeff,
real(dp), dimension(:), intent(in)  latitude,
real(dp), dimension(:), intent(in)  prietayalpha,
real(dp), dimension(:), intent(in)  aeroresist,
real(dp), dimension(:), intent(in)  surfresist 
)
private

get corrected PET for the current timestep and domain

Parameters
[in,out]pet_calc[mm TS-1] estimated PET (if PET is input = corrected values (fAsp*PET))
[in]petlaicorfactorl1PET correction factor based on LAI at level 1
[in]fasp[1] PET correction for Aspect at level 1
[in]harsamcoeff[1] PET Hargreaves Samani coefficient at level 1
[in]latitudelatitude on level 1
[in]prietayalpha[1] PET Priestley Taylor coefficient at level 1
[in]aeroresist[s m-1] PET aerodynamical resitance at level 1
[in]surfresist[s m-1] PET bulk surface resitance at level 1

Definition at line 915 of file mo_meteo_handler.f90.

References mo_mhm_constants::harsamconst, mo_pet::pet_hargreaves(), mo_pet::pet_penman(), mo_pet::pet_priestly(), mo_meteo_temporal_tools::temporal_disagg_flux_daynight(), and mo_meteo_temporal_tools::temporal_disagg_meteo_weights().

Here is the call graph for this function:

◆ get_prec()

subroutine mo_meteo_handler::get_prec ( class(meteo_handler_type), intent(inout)  self,
real(dp), dimension(:), intent(inout)  prec_calc 
)
private

get precipitation for the current timestep and domain

Parameters
[in,out]prec_calc[mm TS-1] precipitation for current time step

Definition at line 1208 of file mo_meteo_handler.f90.

References mo_meteo_temporal_tools::temporal_disagg_flux_daynight(), and mo_meteo_temporal_tools::temporal_disagg_meteo_weights().

Here is the call graph for this function:

◆ get_ssrd()

subroutine mo_meteo_handler::get_ssrd ( class(meteo_handler_type), intent(inout)  self,
real(dp), dimension(:), intent(inout)  ssrd_calc 
)
private

get surface short-wave (solar) radiation downwards for the current timestep and domain

Parameters
[in,out]ssrd_calc[W m2] surface short-wave (solar) radiation downwards for current time step

Definition at line 1288 of file mo_meteo_handler.f90.

References mo_meteo_temporal_tools::temporal_disagg_state_daynight().

Here is the call graph for this function:

◆ get_strd()

subroutine mo_meteo_handler::get_strd ( class(meteo_handler_type), intent(inout)  self,
real(dp), dimension(:), intent(inout)  strd_calc 
)
private

get surface long-wave (thermal) radiation downwards for the current timestep and domain

Parameters
[in,out]strd_calc[W m2] surface long-wave (thermal) radiation downwards for current time step

Definition at line 1359 of file mo_meteo_handler.f90.

References mo_meteo_temporal_tools::temporal_disagg_state_daynight().

Here is the call graph for this function:

◆ get_tann()

subroutine mo_meteo_handler::get_tann ( class(meteo_handler_type), intent(inout)  self,
real(dp), dimension(:), intent(inout)  tann_calc 
)
private

get annual mean surface temperature for the current timestep and domain

Parameters
[in,out]tann_calc[degC] annual mean air temperature

Definition at line 1430 of file mo_meteo_handler.f90.

◆ get_temp()

subroutine mo_meteo_handler::get_temp ( class(meteo_handler_type), intent(inout)  self,
real(dp), dimension(:), intent(inout)  temp_calc 
)
private

get surface temperature for the current timestep and domain

Parameters
[in,out]temp_calc[degC] temperature for current time step

Definition at line 1126 of file mo_meteo_handler.f90.

References mo_meteo_temporal_tools::temporal_disagg_meteo_weights(), and mo_meteo_temporal_tools::temporal_disagg_state_daynight().

Here is the call graph for this function:

◆ init_level2()

subroutine mo_meteo_handler::init_level2 ( class(meteo_handler_type), intent(inout)  self,
type(grid), dimension(:), intent(in)  level0,
type(grid), dimension(:), intent(in)  level1 
)
private

Initialize meteo data and level-2 grid.

Parameters
[in]level0grid information at level-0
[in]level1grid information at level-1 if all meteo data is coupled

Definition at line 505 of file mo_meteo_handler.f90.

References mo_grid::init_lowres_level(), mo_read_spatial_data::read_header_ascii(), and mo_grid::set_domain_indices().

Here is the call graph for this function:

◆ prepare_data()

subroutine mo_meteo_handler::prepare_data ( class(meteo_handler_type), intent(inout)  self,
integer(i4), intent(in)  tt,
integer(i4), intent(in)  idomain,
type(grid), dimension(:), intent(in)  level1,
type(period), dimension(:), intent(in)  simper 
)
private

Prepare meteorological forcings data for a given variable.

Prepare meteorological forcings data for a given variable. Internally this subroutine calls another routine meteo_wrapper for different meterological variables

Changelog
  • Matthias Zink, Jun 2013
    • addded NetCDf reader
  • Rohini Kumar, Aug 2013
    • name changed "inputFormat" to inputFormat_meteo_forcings
  • Matthias Zink, Feb 2014
    • added read in for different PET processes (process 5)
  • Stephan Thober, Jun 2014
    • add chunk_config for chunk read, copied L2 initialization to mo_startup
  • Stephan Thober, Nov 2016
    • moved processMatrix to common variables
  • Stephan Thober, Jan 2017
    • added subroutine for meteo_weights
  • Robert Schweppe Jun 2018
    • refactoring and reformatting
  • Sebastian Müller Mar 2023
    • converted routine to meteo-handler method
Authors
Rohini Kumar
Date
Jan 2013
Parameters
[in]ttcurrent timestep
[in]idomainDomain number
[in]level1grid information at hydrologic level
[in]simperwarmPer + evalPer

Definition at line 645 of file mo_meteo_handler.f90.

References mo_meteo_helper::chunk_config(), mo_meteo_helper::meteo_forcings_wrapper(), and mo_meteo_helper::meteo_weights_wrapper().

Here is the call graph for this function:

◆ set_meteo()

subroutine mo_meteo_handler::set_meteo ( class(meteo_handler_type), intent(inout)  self,
integer(i4), intent(in)  year,
integer(i4), intent(in)  month,
integer(i4), intent(in)  day,
integer(i4), intent(in), optional  hour,
real(dp), dimension(:), intent(in), optional  pre,
real(dp), dimension(:), intent(in), optional  temp,
real(dp), dimension(:), intent(in), optional  pet,
real(dp), dimension(:), intent(in), optional  tmin,
real(dp), dimension(:), intent(in), optional  tmax,
real(dp), dimension(:), intent(in), optional  netrad,
real(dp), dimension(:), intent(in), optional  absvappress,
real(dp), dimension(:), intent(in), optional  windspeed,
real(dp), dimension(:), intent(in), optional  ssrd,
real(dp), dimension(:), intent(in), optional  strd,
real(dp), dimension(:), intent(in), optional  tann 
)
private

set meteo_data from coupling

Definition at line 1462 of file mo_meteo_handler.f90.

◆ single_read()

logical function mo_meteo_handler::single_read ( class(meteo_handler_type), intent(in)  self,
integer(i4), intent(in)  idomain 
)
private

whether meteo data should be read completely at the begining

Returns
True if meteo data is retrieved with a single read
Parameters
[in]idomaincurrent domain

Definition at line 497 of file mo_meteo_handler.f90.

References single_read().

Referenced by single_read().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ update_timestep()

subroutine mo_meteo_handler::update_timestep ( class(meteo_handler_type), intent(inout)  self,
integer(i4), intent(in)  tt,
real(dp), intent(in)  time,
integer(i4), intent(in)  idomain,
type(grid), dimension(:), intent(in)  level1,
type(period), dimension(:), intent(in)  simper 
)
private

update the current time-step of the meteo_handler_type class

Parameters
[in]ttcurrent time step
[in]timecurrent decimal Julian day
[in]idomaincurrent domain
[in]level1grid information at hydrologic level
[in]simperwarmPer + evalPer

Definition at line 575 of file mo_meteo_handler.f90.