5.13.2-dev0
mHM
The mesoscale Hydrological Model
Loading...
Searching...
No Matches
mo_mpr_eval.f90
Go to the documentation of this file.
1!> \file mo_mpr_eval.f90
2!> \brief \copybrief mo_mpr_eval
3!> \details \copydetails mo_mpr_eval
4
5!> \brief Runs MPR
6!> \details Runs MPR and writes to global effective parameters
7!> \authors Robert Schweppe
8!> \date Feb 2018
9!> \copyright Copyright 2005-\today, the mHM Developers, Luis Samaniego, Sabine Attinger: All rights reserved.
10!! mHM is released under the LGPLv3+ license \license_note
11!> \ingroup f_mpr
13
14 USE mo_kind, ONLY : i4, dp
15
16 IMPLICIT NONE
17
18 PRIVATE
19
20 PUBLIC :: mpr_eval
21
22 ! ------------------------------------------------------------------
23
24CONTAINS
25
26 ! ------------------------------------------------------------------
27
28 ! NAME
29 ! mpr_eval
30
31 ! PURPOSE
32 !> \brief Runs MPR and writes to global effective parameters
33
34 !> \details Runs MPR and writes to global effective parameters
35
36 ! INTENT(IN), OPTIONAL
37 !> \param[in] "real(dp), dimension(:), optional :: parameterset" a set of global parameter (gamma) to run mHM,
38 !> DIMENSION [no. of global_Parameters]
39
40 ! HISTORY
41 !> \authors Juliane Mai, Rohini Kumar
42
43 !> \date Feb 2013
44
45 ! Modifications:
46 ! R. Kumar Jun 2013 - restart_flag_states_read is passed to mhm call for the soil moisture initalisation
47 ! R. Kumar Jun 2013 - frac_sealed_city_area is added
48 ! R. Kumar & S. Thober Aug 2013 - code change to incorporate output timestep during writing of the netcdf file
49 ! R. Kumar Aug 2013 - added iFlag_LAI_data_format to handle LAI options, and changed within the code made accordingly
50 ! R. Kumar, J. Mai Sep 2013 - Splitting allocation and initialization of arrays
51 ! R. Kumar Nov 2013 - update intent variables in documentation
52 ! L. Samaniego Nov 2013 - relational statements == to .eq., etc.
53 ! M. Zink Feb 2014 - added PET calculation: Hargreaves-Samani (Process 5)
54 ! M. Zink Mar 2014 - added inflow from upstream areas
55 ! Stephan Thober Jun 2014 - added chunk read for meteorological input
56 ! Stephan Thober Jun 2014 - updated flag for read_restart
57 ! M. Cuntz & J. Mai Nov 2014 - LAI input from daily, monthly or yearly files
58 ! Matthias Zink Dec 2014 - adopted inflow gauges to ignore headwater cells
59 ! Stephan Thober Aug 2015 - moved writing of daily discharge to mo_write_routing, included routing related variables from mRM
60 ! David Schaefer Aug 2015 - changed to new netcdf-writing scheme
61 ! Stephan Thober Sep 2015 - updated mrm_routing call
62 ! O. Rakovec, R. Kumar Oct 2015 - added optional output for domain averaged TWS
63 ! Rohini Kumar Mar 2016 - changes for handling multiple soil database options
64 ! Stephan Thober Nov 2016 - added two options for routing
65 ! Rohini Kuamr Dec 2016 - option to handle monthly mean gridded fields of LAI
66 ! Stephan Thober Jan 2017 - added prescribed weights for tavg and pet
67 ! Zink M. Demirel C. Mar 2017 - Added Jarvis soil water stress function at SM process(3)
68 ! Robert Schweppe Dec 2017 - extracted call to mpr from inside mhm
69
70 SUBROUTINE mpr_eval(parameterset, opti_domain_indices)
71
73 use mo_message, only : message
81 ! neutron count
83
84 use mo_multi_param_reg, only : mpr
85 use mo_string_utils, only : num2str
86 use mo_timer, only : timer_clear, timer_get, timer_start, timer_stop
87
88 implicit none
89
90 ! a set of global parameter (gamma) to run mHM, DIMENSION [no. of global_Parameters]
91 real(dp), dimension(:), intent(in), optional :: parameterset
92
93 integer(i4), dimension(:), optional, intent(in) :: opti_domain_indices
94
95 ! Counters
96 integer(i4) :: idomain, ii, ndomains, itimer
97
98 ! start and end index at level 0 for current domain
99 integer(i4) :: s0, e0
100
101 ! start and end index at level 1 for current domain
102 integer(i4) :: s1, e1
103
104
105 !-------------------------------------------------------------------
106 ! NOW call MPR
107 !-------------------------------------------------------------------
108 call message(' Executing MPR ...')
109 itimer = 10
110 call timer_start(itimer)
111
112 !----------------------------------------
113 ! loop over domains
114 !----------------------------------------
115 if (present(opti_domain_indices)) then
116 ndomains = size(opti_domain_indices)
117 else
118 ndomains = domainmeta%nDomains
119 end if
120 do ii = 1, ndomains
121 if (present(opti_domain_indices)) then
122 idomain = opti_domain_indices(ii)
123 else
124 idomain = ii
125 end if
126
127 ! get domain information
128 s0 = level0(domainmeta%L0DataFrom(idomain))%iStart
129 e0 = level0(domainmeta%L0DataFrom(idomain))%iEnd
130 s1 = level1(idomain)%iStart
131 e1 = level1(idomain)%iEnd
132
133 call mpr(level0(domainmeta%L0DataFrom(idomain))%mask, l0_geounit(s0 : e0), &
134 l0_soilid(s0 : e0, :), l0_asp(s0 : e0), l0_gridded_lai(s0 : e0, :), &
135 l0_lcover(s0 : e0, :), l0_slope_emp(s0 : e0), &
136 pack(level0(domainmeta%L0DataFrom(idomain))%y, level0(domainmeta%L0DataFrom(idomain))%mask), &
137 level0(domainmeta%L0DataFrom(idomain))%Id, &
138 l0_l1_remap(idomain)%upper_bound, l0_l1_remap(idomain)%lower_bound, &
139 l0_l1_remap(idomain)%left_bound, l0_l1_remap(idomain)%right_bound, &
140 l0_l1_remap(idomain)%n_subcells, &
141 l1_fsealed(s1 : e1, :, :), &
142 l1_alpha(s1 : e1, :, :), l1_degdayinc(s1 : e1, :, :), l1_degdaymax(s1 : e1, :, :), &
143 l1_degdaynopre(s1 : e1, :, :), l1_fasp(s1 : e1, :, :), l1_harsamcoeff(s1 : e1, :, :), &
144 l1_prietayalpha(s1 : e1, :, :), l1_aeroresist(s1 : e1, :, :), l1_surfresist(s1 : e1, :, :), &
145 l1_froots(s1 : e1, :, :), l1_kfastflow(s1 : e1, :, :), &
146 l1_kslowflow(s1 : e1, :, :), l1_kbaseflow(s1 : e1, :, :), l1_kperco(s1 : e1, :, :), &
147 l1_karstloss(s1 : e1, :, :), l1_soilmoistfc(s1 : e1, :, :), l1_soilmoistsat(s1 : e1, :, :), &
148 l1_soilmoistexp(s1 : e1, :, :), l1_jarvis_thresh_c1(s1 : e1, :, :), &
149 l1_tempthresh(s1 : e1, :, :), l1_unsatthresh(s1 : e1, :, :), l1_sealedthresh(s1 : e1, :, :), &
150 l1_wiltingpoint(s1 : e1, :, :), l1_maxinter(s1 : e1, :, :), l1_petlaicorfactor(s1 : e1, :, :), &
151 l1_no_count(s1:e1,:,:), l1_bulkdens(s1:e1,:,:), l1_latticewater(s1:e1,:,:), l1_cosmicl3(s1:e1,:,:), &
152 parameterset )
153
154 end do
155 call timer_stop(itimer)
156 call message(' in ', trim(num2str(timer_get(itimer), '(F9.3)')), ' seconds.')
157 call timer_clear(itimer)
158
159 END SUBROUTINE mpr_eval
160
161END MODULE mo_mpr_eval
Provides structures needed by mHM, mRM and/or mpr.
type(domain_meta), public domainmeta
integer(i4), dimension(:, :), allocatable, public l0_lcover
type(grid), dimension(:), allocatable, target, public level1
type(grid), dimension(:), allocatable, target, public level0
type(gridremapper), dimension(:), allocatable, public l0_l1_remap
Runs MPR.
subroutine, public mpr_eval(parameterset, opti_domain_indices)
Runs MPR and writes to global effective parameters.
Global variables for mpr only.
real(dp), dimension(:, :, :), allocatable, public l1_degdaymax
real(dp), dimension(:, :, :), allocatable, public l1_soilmoistexp
real(dp), dimension(:, :, :), allocatable, public l1_harsamcoeff
real(dp), dimension(:, :, :), allocatable, public l1_karstloss
real(dp), dimension(:, :, :), allocatable, public l1_unsatthresh
real(dp), dimension(:, :, :), allocatable, public l1_kperco
real(dp), dimension(:), allocatable, public l0_slope_emp
real(dp), dimension(:, :, :), allocatable, public l1_alpha
real(dp), dimension(:, :, :), allocatable, public l1_surfresist
integer(i4), dimension(:), allocatable, public l0_geounit
real(dp), dimension(:,:,:), allocatable, public l1_cosmicl3
real(dp), dimension(:), allocatable, public l0_asp
real(dp), dimension(:, :, :), allocatable, public l1_degdayinc
real(dp), dimension(:, :, :), allocatable, public l1_petlaicorfactor
real(dp), dimension(:, :, :), allocatable, public l1_fasp
real(dp), dimension(:, :, :), allocatable, public l1_kbaseflow
real(dp), dimension(:, :, :), allocatable, public l1_soilmoistfc
real(dp), dimension(:, :, :), allocatable, public l1_maxinter
real(dp), dimension(:,:,:), allocatable, public l1_bulkdens
real(dp), dimension(:, :, :), allocatable, public l1_degdaynopre
real(dp), dimension(:, :, :), allocatable, public l1_wiltingpoint
integer(i4), dimension(:, :), allocatable, public l0_soilid
real(dp), dimension(:, :), allocatable, public l0_gridded_lai
real(dp), dimension(:, :, :), allocatable, public l1_prietayalpha
real(dp), dimension(:, :, :), allocatable, public l1_fsealed
real(dp), dimension(:, :, :), allocatable, public l1_froots
real(dp), dimension(:, :, :), allocatable, public l1_jarvis_thresh_c1
real(dp), dimension(:, :, :), allocatable, public l1_tempthresh
real(dp), dimension(:,:,:), allocatable, public l1_no_count
real(dp), dimension(:, :, :), allocatable, public l1_kfastflow
real(dp), dimension(:, :, :), allocatable, public l1_kslowflow
real(dp), dimension(:, :, :), allocatable, public l1_aeroresist
real(dp), dimension(:, :, :), allocatable, public l1_sealedthresh
real(dp), dimension(:, :, :), allocatable, public l1_soilmoistsat
real(dp), dimension(:,:,:), allocatable, public l1_latticewater
Multiscale parameter regionalization (MPR).
subroutine, public mpr(mask0, geounit0, soilid0, asp0, gridded_lai0, lcover0, slope_emp0, y0, id0, upper_bound1, lower_bound1, left_bound1, right_bound1, n_subcells1, fsealed1, alpha1, degdayinc1, degdaymax1, degdaynopre1, fasp1, harsamcoeff1, prietayalpha1, aeroresist1, surfresist1, froots1, kfastflow1, kslowflow1, kbaseflow1, kperco1, karstloss1, soilmoistfc1, soilmoistsat1, soilmoistexp1, jarvis_thresh_c1, tempthresh1, unsatthresh1, sealedthresh1, wiltingpoint1, maxinter1, petlaicorfactor, no_count1, bulkdens1, latticewater1, cosmicl31, parameterset)
Regionalizing and Upscaling process parameters.