5.13.3-dev0
mHM
The mesoscale Hydrological Model
Loading...
Searching...
No Matches
mo_mpr_read_config.f90
Go to the documentation of this file.
1!> \file mo_mpr_read_config.f90
2!> \brief \copybrief mo_mpr_read_config
3!> \details \copydetails mo_mpr_read_config
4
5!> \brief read mpr config
6!> \details This module contains all mpr subroutines related to reading the mpr configuration from file.
7!> \changelog
8!! - Robert Schweppe Dec 2017
9!! - adapted for MPR
10!! - Robert Schweppe Jun 2018
11!! - refactoring and reformatting
12!! - M. Cuneyd Demirel, Simon Stisen Jun 2020
13!! - added Feddes and FC dependency on root fraction coefficient processCase(3) = 4
14!! - Rohini Kumar Oct 2021
15!! - Added Neutron count module to mHM integrate into develop branch (5.11.2)
16!> \authors Stephan Thober
17!> \date Aug 2015
18!> \copyright Copyright 2005-\today, the mHM Developers, Luis Samaniego, Sabine Attinger: All rights reserved.
19!! mHM is released under the LGPLv3+ license \license_note
20!> \ingroup f_mpr
22
23 use mo_kind, only : i4, dp
24
25 implicit none
26
27 public :: mpr_read_config
28
29contains
30
31 ! ------------------------------------------------------------------
32
33 ! NAME
34 ! mpr_read_config
35
36 ! PURPOSE
37 !> \brief Read the general config of mpr
38
39 !> \details Depending on the variable mrm_coupling_config, the
40 !> mRM config is either read from mrm.nml and parameters from
41 !> mrm_parameter.nml or copied from mHM.
42
43 ! INTENT(IN)
44 !> \param[in] "character(*) :: file_namelist"
45 !> \param[in] "character(*) :: file_namelist_param"
46
47 ! HISTORY
48 !> \authors Stephan Thober
49
50 !> \date Aug 2015
51
52 ! Modifications:
53 ! Stephan Thober Sep 2015 - removed stop condition when routing resolution is smaller than hydrologic resolution
54 ! Stephan Thober Oct 2015 - added NLoutputResults namelist, fileLatLon to directories_general namelist, and readLatLon flag
55 ! Robert Schweppe Dec 2017 - adapted for MPR
56 ! Rohini Kumar Oct 2021 - Added Neutron count module to mHM integrate into develop branch (5.11.2)
57
58 subroutine mpr_read_config(file_namelist, file_namelist_param)
59
60 use mo_namelists, only : &
66 nml_snow1, &
73 nml_pet0, &
74 nml_pet1, &
75 nml_pet2, &
76 nml_pet3, &
82 use mo_append, only : append
86 use mo_message, only : message, error_message
87 use mo_message, only : message, error_message
88 use mo_mpr_constants, only : maxgeounit, &
93 use mo_string_utils, only : num2str
94 use mo_utils, only : eq
95
96 implicit none
97
98 character(*), intent(in) :: file_namelist
99
100 character(*), intent(in) :: file_namelist_param
101
102 integer(i4) :: ii
103
104 ! depth of the single horizons
105 real(dp), dimension(maxNoSoilHorizons) :: soil_depth
106
107 ! directory of gridded LAI data
108 ! used when timeStep_LAI_input<0
109 character(256), dimension(maxNoDomains) :: dir_gridded_lai
110
111 character(256) :: dummy
112
113 ! space holder for routing parameters
114 real(dp), dimension(5, nColPars) :: dummy_2d_dp
115
116 ! space holder for routing parameters
117 real(dp), dimension(1, nColPars) :: dummy_2d_dp_2
118
119 real(dp), dimension(nColPars) :: canopyinterceptionfactor
120
121 real(dp), dimension(nColPars) :: snowtreshholdtemperature
122 real(dp), dimension(nColPars) :: degreedayfactor_forest
123 real(dp), dimension(nColPars) :: degreedayfactor_impervious
124 real(dp), dimension(nColPars) :: degreedayfactor_pervious
125 real(dp), dimension(nColPars) :: increasedegreedayfactorbyprecip
126 real(dp), dimension(nColPars) :: maxdegreedayfactor_forest
127 real(dp), dimension(nColPars) :: maxdegreedayfactor_impervious
128 real(dp), dimension(nColPars) :: maxdegreedayfactor_pervious
129
130 real(dp), dimension(nColPars) :: orgmattercontent_forest
131 real(dp), dimension(nColPars) :: orgmattercontent_impervious
132 real(dp), dimension(nColPars) :: orgmattercontent_pervious
133 real(dp), dimension(nColPars) :: ptf_lower66_5_constant
134 real(dp), dimension(nColPars) :: ptf_lower66_5_clay
135 real(dp), dimension(nColPars) :: ptf_lower66_5_db
136 real(dp), dimension(nColPars) :: ptf_higher66_5_constant
137 real(dp), dimension(nColPars) :: ptf_higher66_5_clay
138 real(dp), dimension(nColPars) :: ptf_higher66_5_db
139 real(dp), dimension(nColPars) :: ptf_ks_constant
140 real(dp), dimension(nColPars) :: ptf_ks_sand
141 real(dp), dimension(nColPars) :: ptf_ks_clay
142 real(dp), dimension(nColPars) :: ptf_ks_curveslope
143 real(dp), dimension(nColPars) :: rootfractioncoefficient_forest
144 real(dp), dimension(nColPars) :: rootfractioncoefficient_impervious
145 real(dp), dimension(nColPars) :: rootfractioncoefficient_pervious
146 real(dp), dimension(nColPars) :: infiltrationshapefactor
147 real(dp), dimension(nColPars) :: jarvis_sm_threshold_c1
148
149 real(dp), dimension(nColPars) :: fcmin_glob
150 real(dp), dimension(nColPars) :: fcdelta_glob
151 real(dp), dimension(nColPars) :: rootfractioncoefficient_sand
152 real(dp), dimension(nColPars) :: rootfractioncoefficient_clay
153
154 real(dp), dimension(nColPars) :: imperviousstoragecapacity
155
156 real(dp), dimension(nColPars) :: pet_a_forest
157 real(dp), dimension(nColPars) :: pet_a_impervious
158 real(dp), dimension(nColPars) :: pet_a_pervious
159 real(dp), dimension(nColPars) :: pet_b
160 real(dp), dimension(nColPars) :: pet_c
161
162 real(dp), dimension(nColPars) :: mincorrectionfactorpet
163 real(dp), dimension(nColPars) :: maxcorrectionfactorpet
164 real(dp), dimension(nColPars) :: aspecttresholdpet
165 real(dp), dimension(nColPars) :: hargreavessamanicoeff
166
167 real(dp), dimension(nColPars) :: priestleytaylorcoeff
168 real(dp), dimension(nColPars) :: priestleytaylorlaicorr
169
170 real(dp), dimension(nColPars) :: canopyheigth_forest
171 real(dp), dimension(nColPars) :: canopyheigth_impervious
172 real(dp), dimension(nColPars) :: canopyheigth_pervious
173 real(dp), dimension(nColPars) :: displacementheight_coeff
174 real(dp), dimension(nColPars) :: roughnesslength_momentum_coeff
175 real(dp), dimension(nColPars) :: roughnesslength_heat_coeff
176 real(dp), dimension(nColPars) :: stomatal_resistance
177
178 real(dp), dimension(nColPars) :: interflowstoragecapacityfactor
179 real(dp), dimension(nColPars) :: interflowrecession_slope
180 real(dp), dimension(nColPars) :: fastinterflowrecession_forest
181 real(dp), dimension(nColPars) :: slowinterflowrecession_ks
182 real(dp), dimension(nColPars) :: exponentslowinterflow
183
184 real(dp), dimension(nColPars) :: rechargecoefficient
185 real(dp), dimension(nColPars) :: rechargefactor_karstic
186 real(dp), dimension(nColPars) :: gain_loss_gwreservoir_karstic
187
188 real(dp), dimension(maxGeoUnit, nColPars) :: geoparam
189
190 real(dp), dimension(nColPars) :: desilets_n0
191 real(dp), dimension(nColPars) :: desilets_lw0
192 real(dp), dimension(nColPars) :: desilets_lw1
193
194
195 real(dp), dimension(nColPars) :: cosmic_n0
196 real(dp), dimension(nColPars) :: cosmic_n1
197 real(dp), dimension(nColPars) :: cosmic_n2
198 real(dp), dimension(nColPars) :: cosmic_alpha0
199 real(dp), dimension(nColPars) :: cosmic_alpha1
200 real(dp), dimension(nColPars) :: cosmic_l30
201 real(dp), dimension(nColPars) :: cosmic_l31
202 real(dp), dimension(nColPars) :: cosmic_lw0
203 real(dp), dimension(nColPars) :: cosmic_lw1
204
205 integer(i4) :: idomain, domainid
206
207 !===============================================================
208 ! INITIALIZATION
209 !===============================================================
210 dummy_2d_dp = nodata_dp
211 dummy_2d_dp_2 = nodata_dp
212
213 !===============================================================
214 ! Read namelist for LCover
215 !===============================================================
216 call nml_lcover_mpr%read(file_namelist)
217 fracsealed_cityarea = nml_lcover_mpr%fracSealed_cityArea
218
219 !===============================================================
220 ! Read soil layering information
221 !===============================================================
222 call nml_soildata%read(file_namelist)
223 iflag_soildb = nml_soildata%iFlag_soilDB
224 tillagedepth = nml_soildata%tillageDepth
225 nsoilhorizons_mhm = nml_soildata%nSoilHorizons_mHM
226 soil_depth = nml_soildata%soil_Depth
227
229 horizondepth_mhm(:) = 0.0_dp
230 ! last layer is reset to 0 in MPR in case of iFlag_soilDB is 0
232
233 ! counter checks -- soil horizons
235 call error_message('***ERROR: Number of soil horizons is resticted to ', trim(num2str(maxnosoilhorizons)), '!')
236
237 ! the default is the HorizonDepths are all set up to last
238 ! as is the default for option-1 where horizon specific information are taken into consideration
239 if(iflag_soildb .eq. 0) then
240 ! classical mhm soil database
242 else if(iflag_soildb .ne. 1) then
243 call error_message('***ERROR: iFlag_soilDB option given does not exist. Only 0 and 1 is taken at the moment.')
244 call error_message('***ERROR: iFlag_soilDB option given does not exist. Only 0 and 1 is taken at the moment.')
245 end if
246
247 ! some consistency checks for the specification of the tillage depth
248 if(iflag_soildb .eq. 1) then
249 if(count(abs(horizondepth_mhm(:) - tillagedepth) .lt. eps_dp) .eq. 0) &
250 call error_message('***ERROR: Soil tillage depth must conform with one of the specified horizon (lower) depth.')
251 end if
252
253 !===============================================================
254 ! Read LAI related information
255 !===============================================================
256 call nml_lai_data_information%read(file_namelist)
257 inputformat_gridded_lai = nml_lai_data_information%inputFormat_gridded_LAI
259
260 if (timestep_lai_input .ne. 0) then
261 !===============================================================
262 ! Read namelist for main directories
263 !===============================================================
264 call nml_directories_mpr%read(file_namelist)
265 dir_gridded_lai = nml_directories_mpr%dir_gridded_LAI
266
267 allocate(dirgridded_lai(domainmeta%nDomains))
268 do idomain = 1, domainmeta%nDomains
269 domainid = domainmeta%indices(idomain)
270 dirgridded_lai(idomain) = dir_gridded_lai(domainid)
271 end do
272
273 if (timestep_lai_input .GT. 1) &
274 call error_message('***ERROR: option for selected timeStep_LAI_input not coded yet')
275 end if
276
277 !===============================================================
278 ! Read namelist global parameters
279 !===============================================================
280 ! decide which parameters to read depending on specified processes
281
282 ! Process 1 - interception
283 select case (processmatrix(1, 1))
284 ! 1 - maximum Interception
285 case(1)
286 call nml_interception1%read(file_namelist_param)
287 canopyinterceptionfactor = nml_interception1%canopyInterceptionFactor
288
289 processmatrix(1, 2) = 1_i4
290 processmatrix(1, 3) = 1_i4
291 call append(global_parameters, reshape(canopyinterceptionfactor, (/1, ncolpars/)))
292
293 call append(global_parameters_name, (/ &
294 'canopyInterceptionFactor'/))
295
296 ! check if parameter are in range
297 if (.not. in_bound(global_parameters)) &
298 call error_message('***ERROR: parameter in namelist "interception1" out of bound in ', trim(adjustl(file_namelist_param)))
299
300 case DEFAULT
301 call error_message('***ERROR: Process description for process "interception" does not exist!')
302 stop
303 end select
304
305 ! Process 2 - snow
306 select case (processmatrix(2, 1))
307 ! 1 - degree-day approach
308 case(1)
309 call nml_snow1%read(file_namelist_param)
310 snowtreshholdtemperature = nml_snow1%snowTreshholdTemperature
311 degreedayfactor_forest = nml_snow1%degreeDayFactor_forest
312 degreedayfactor_impervious = nml_snow1%degreeDayFactor_impervious
313 degreedayfactor_pervious = nml_snow1%degreeDayFactor_pervious
314 increasedegreedayfactorbyprecip = nml_snow1%increaseDegreeDayFactorByPrecip
315 maxdegreedayfactor_forest = nml_snow1%maxDegreeDayFactor_forest
316 maxdegreedayfactor_impervious = nml_snow1%maxDegreeDayFactor_impervious
317 maxdegreedayfactor_pervious = nml_snow1%maxDegreeDayFactor_pervious
318
319 processmatrix(2, 2) = 8_i4
320 processmatrix(2, 3) = sum(processmatrix(1 : 2, 2))
321 call append(global_parameters, reshape(snowtreshholdtemperature, (/1, ncolpars/)))
322 call append(global_parameters, reshape(degreedayfactor_forest, (/1, ncolpars/)))
323 call append(global_parameters, reshape(degreedayfactor_impervious, (/1, ncolpars/)))
324 call append(global_parameters, reshape(degreedayfactor_pervious, (/1, ncolpars/)))
325 call append(global_parameters, reshape(increasedegreedayfactorbyprecip, (/1, ncolpars/)))
326 call append(global_parameters, reshape(maxdegreedayfactor_forest, (/1, ncolpars/)))
327 call append(global_parameters, reshape(maxdegreedayfactor_impervious, (/1, ncolpars/)))
328 call append(global_parameters, reshape(maxdegreedayfactor_pervious, (/1, ncolpars/)))
329
330 call append(global_parameters_name, (/ &
331 'snowTreshholdTemperature ', &
332 'degreeDayFactor_forest ', &
333 'degreeDayFactor_impervious ', &
334 'degreeDayFactor_pervious ', &
335 'increaseDegreeDayFactorByPrecip', &
336 'maxDegreeDayFactor_forest ', &
337 'maxDegreeDayFactor_impervious ', &
338 'maxDegreeDayFactor_pervious '/))
339
340 ! check if parameter are in range
341 if (.not. in_bound(global_parameters)) &
342 call error_message('***ERROR: parameter in namelist "snow1" out of bound in ', trim(adjustl(file_namelist_param)))
343
344 case DEFAULT
345 call error_message('***ERROR: Process description for process "snow" does not exist!')
346 end select
347
348 ! Process 3 - soilmoisture
349 select case (processmatrix(3, 1))
350
351 ! 1 - Feddes equation for PET reduction, bucket approach, Brooks-Corey like
352 case(1)
353 call nml_soilmoisture1%read(file_namelist_param)
354 orgmattercontent_forest = nml_soilmoisture1%orgMatterContent_forest
355 orgmattercontent_impervious = nml_soilmoisture1%orgMatterContent_impervious
356 orgmattercontent_pervious = nml_soilmoisture1%orgMatterContent_pervious
357 ptf_lower66_5_constant = nml_soilmoisture1%PTF_lower66_5_constant
358 ptf_lower66_5_clay = nml_soilmoisture1%PTF_lower66_5_clay
359 ptf_lower66_5_db = nml_soilmoisture1%PTF_lower66_5_Db
360 ptf_higher66_5_constant = nml_soilmoisture1%PTF_higher66_5_constant
361 ptf_higher66_5_clay = nml_soilmoisture1%PTF_higher66_5_clay
362 ptf_higher66_5_db = nml_soilmoisture1%PTF_higher66_5_Db
363 ptf_ks_constant = nml_soilmoisture1%PTF_Ks_constant
364 ptf_ks_sand = nml_soilmoisture1%PTF_Ks_sand
365 ptf_ks_clay = nml_soilmoisture1%PTF_Ks_clay
366 ptf_ks_curveslope = nml_soilmoisture1%PTF_Ks_curveSlope
367 rootfractioncoefficient_forest = nml_soilmoisture1%rootFractionCoefficient_forest
368 rootfractioncoefficient_impervious = nml_soilmoisture1%rootFractionCoefficient_impervious
369 rootfractioncoefficient_pervious = nml_soilmoisture1%rootFractionCoefficient_pervious
370 infiltrationshapefactor = nml_soilmoisture1%infiltrationShapeFactor
371
372 processmatrix(3, 2) = 17_i4
373 processmatrix(3, 3) = sum(processmatrix(1 : 3, 2))
374 call append(global_parameters, reshape(orgmattercontent_forest, (/1, ncolpars/)))
375 call append(global_parameters, reshape(orgmattercontent_impervious, (/1, ncolpars/)))
376 call append(global_parameters, reshape(orgmattercontent_pervious, (/1, ncolpars/)))
377 call append(global_parameters, reshape(ptf_lower66_5_constant, (/1, ncolpars/)))
378 call append(global_parameters, reshape(ptf_lower66_5_clay, (/1, ncolpars/)))
379 call append(global_parameters, reshape(ptf_lower66_5_db, (/1, ncolpars/)))
380 call append(global_parameters, reshape(ptf_higher66_5_constant, (/1, ncolpars/)))
381 call append(global_parameters, reshape(ptf_higher66_5_clay, (/1, ncolpars/)))
382 call append(global_parameters, reshape(ptf_higher66_5_db, (/1, ncolpars/)))
383 call append(global_parameters, reshape(ptf_ks_constant, (/1, ncolpars/)))
384 call append(global_parameters, reshape(ptf_ks_sand, (/1, ncolpars/)))
385 call append(global_parameters, reshape(ptf_ks_clay, (/1, ncolpars/)))
386 call append(global_parameters, reshape(ptf_ks_curveslope, (/1, ncolpars/)))
387 call append(global_parameters, reshape(rootfractioncoefficient_forest, (/1, ncolpars/)))
388 call append(global_parameters, reshape(rootfractioncoefficient_impervious, (/1, ncolpars/)))
389 call append(global_parameters, reshape(rootfractioncoefficient_pervious, (/1, ncolpars/)))
390 call append(global_parameters, reshape(infiltrationshapefactor, (/1, ncolpars/)))
391
392 call append(global_parameters_name, (/ &
393 'orgMatterContent_forest ', &
394 'orgMatterContent_impervious ', &
395 'orgMatterContent_pervious ', &
396 'PTF_lower66_5_constant ', &
397 'PTF_lower66_5_clay ', &
398 'PTF_lower66_5_Db ', &
399 'PTF_higher66_5_constant ', &
400 'PTF_higher66_5_clay ', &
401 'PTF_higher66_5_Db ', &
402 'PTF_Ks_constant ', &
403 'PTF_Ks_sand ', &
404 'PTF_Ks_clay ', &
405 'PTF_Ks_curveSlope ', &
406 'rootFractionCoefficient_forest ', &
407 'rootFractionCoefficient_impervious', &
408 'rootFractionCoefficient_pervious ', &
409 'infiltrationShapeFactor '/))
410
411 ! check if parameter are in range
412 if (.not. in_bound(global_parameters)) &
413 call error_message('***ERROR: parameter in namelist "soilmoisture1" out of bound in ', trim(adjustl(file_namelist_param)))
414
415 ! 2- Jarvis equation for PET reduction, bucket approach, Brooks-Corey like
416 case(2)
417 call nml_soilmoisture2%read(file_namelist_param)
418 orgmattercontent_forest = nml_soilmoisture2%orgMatterContent_forest
419 orgmattercontent_impervious = nml_soilmoisture2%orgMatterContent_impervious
420 orgmattercontent_pervious = nml_soilmoisture2%orgMatterContent_pervious
421 ptf_lower66_5_constant = nml_soilmoisture2%PTF_lower66_5_constant
422 ptf_lower66_5_clay = nml_soilmoisture2%PTF_lower66_5_clay
423 ptf_lower66_5_db = nml_soilmoisture2%PTF_lower66_5_Db
424 ptf_higher66_5_constant = nml_soilmoisture2%PTF_higher66_5_constant
425 ptf_higher66_5_clay = nml_soilmoisture2%PTF_higher66_5_clay
426 ptf_higher66_5_db = nml_soilmoisture2%PTF_higher66_5_Db
427 ptf_ks_constant = nml_soilmoisture2%PTF_Ks_constant
428 ptf_ks_sand = nml_soilmoisture2%PTF_Ks_sand
429 ptf_ks_clay = nml_soilmoisture2%PTF_Ks_clay
430 ptf_ks_curveslope = nml_soilmoisture2%PTF_Ks_curveSlope
431 rootfractioncoefficient_forest = nml_soilmoisture2%rootFractionCoefficient_forest
432 rootfractioncoefficient_impervious = nml_soilmoisture2%rootFractionCoefficient_impervious
433 rootfractioncoefficient_pervious = nml_soilmoisture2%rootFractionCoefficient_pervious
434 infiltrationshapefactor = nml_soilmoisture2%infiltrationShapeFactor
435 jarvis_sm_threshold_c1 = nml_soilmoisture2%jarvis_sm_threshold_c1
436
437 processmatrix(3, 2) = 18_i4
438 processmatrix(3, 3) = sum(processmatrix(1 : 3, 2))
439 call append(global_parameters, reshape(orgmattercontent_forest, (/1, ncolpars/)))
440 call append(global_parameters, reshape(orgmattercontent_impervious, (/1, ncolpars/)))
441 call append(global_parameters, reshape(orgmattercontent_pervious, (/1, ncolpars/)))
442 call append(global_parameters, reshape(ptf_lower66_5_constant, (/1, ncolpars/)))
443 call append(global_parameters, reshape(ptf_lower66_5_clay, (/1, ncolpars/)))
444 call append(global_parameters, reshape(ptf_lower66_5_db, (/1, ncolpars/)))
445 call append(global_parameters, reshape(ptf_higher66_5_constant, (/1, ncolpars/)))
446 call append(global_parameters, reshape(ptf_higher66_5_clay, (/1, ncolpars/)))
447 call append(global_parameters, reshape(ptf_higher66_5_db, (/1, ncolpars/)))
448 call append(global_parameters, reshape(ptf_ks_constant, (/1, ncolpars/)))
449 call append(global_parameters, reshape(ptf_ks_sand, (/1, ncolpars/)))
450 call append(global_parameters, reshape(ptf_ks_clay, (/1, ncolpars/)))
451 call append(global_parameters, reshape(ptf_ks_curveslope, (/1, ncolpars/)))
452 call append(global_parameters, reshape(rootfractioncoefficient_forest, (/1, ncolpars/)))
453 call append(global_parameters, reshape(rootfractioncoefficient_impervious, (/1, ncolpars/)))
454 call append(global_parameters, reshape(rootfractioncoefficient_pervious, (/1, ncolpars/)))
455 call append(global_parameters, reshape(infiltrationshapefactor, (/1, ncolpars/)))
456 call append(global_parameters, reshape(jarvis_sm_threshold_c1, (/1, ncolpars/)))
457
458 call append(global_parameters_name, (/ &
459 'orgMatterContent_forest ', &
460 'orgMatterContent_impervious ', &
461 'orgMatterContent_pervious ', &
462 'PTF_lower66_5_constant ', &
463 'PTF_lower66_5_clay ', &
464 'PTF_lower66_5_Db ', &
465 'PTF_higher66_5_constant ', &
466 'PTF_higher66_5_clay ', &
467 'PTF_higher66_5_Db ', &
468 'PTF_Ks_constant ', &
469 'PTF_Ks_sand ', &
470 'PTF_Ks_clay ', &
471 'PTF_Ks_curveSlope ', &
472 'rootFractionCoefficient_forest ', &
473 'rootFractionCoefficient_impervious', &
474 'rootFractionCoefficient_pervious ', &
475 'infiltrationShapeFactor ', &
476 'jarvis_sm_threshold_c1 '/))
477
478 ! check if parameter are in range
479 if (.not. in_bound(global_parameters)) &
480 call error_message('***ERROR: parameter in namelist "soilmoisture2" out of bound in ', trim(adjustl(file_namelist_param)))
481
482 ! 3- Jarvis equation for ET reduction and FC dependency on root fraction coefficient
483 case(3)
484 call nml_soilmoisture3%read(file_namelist_param)
485 orgmattercontent_forest = nml_soilmoisture3%orgMatterContent_forest
486 orgmattercontent_impervious = nml_soilmoisture3%orgMatterContent_impervious
487 orgmattercontent_pervious = nml_soilmoisture3%orgMatterContent_pervious
488 ptf_lower66_5_constant = nml_soilmoisture3%PTF_lower66_5_constant
489 ptf_lower66_5_clay = nml_soilmoisture3%PTF_lower66_5_clay
490 ptf_lower66_5_db = nml_soilmoisture3%PTF_lower66_5_Db
491 ptf_higher66_5_constant = nml_soilmoisture3%PTF_higher66_5_constant
492 ptf_higher66_5_clay = nml_soilmoisture3%PTF_higher66_5_clay
493 ptf_higher66_5_db = nml_soilmoisture3%PTF_higher66_5_Db
494 ptf_ks_constant = nml_soilmoisture3%PTF_Ks_constant
495 ptf_ks_sand = nml_soilmoisture3%PTF_Ks_sand
496 ptf_ks_clay = nml_soilmoisture3%PTF_Ks_clay
497 ptf_ks_curveslope = nml_soilmoisture3%PTF_Ks_curveSlope
498 rootfractioncoefficient_forest = nml_soilmoisture3%rootFractionCoefficient_forest
499 rootfractioncoefficient_impervious = nml_soilmoisture3%rootFractionCoefficient_impervious
500 rootfractioncoefficient_pervious = nml_soilmoisture3%rootFractionCoefficient_pervious
501 infiltrationshapefactor = nml_soilmoisture3%infiltrationShapeFactor
502 rootfractioncoefficient_sand = nml_soilmoisture3%rootFractionCoefficient_sand
503 rootfractioncoefficient_clay = nml_soilmoisture3%rootFractionCoefficient_clay
504 fcmin_glob = nml_soilmoisture3%FCmin_glob
505 fcdelta_glob = nml_soilmoisture3%FCdelta_glob
506 jarvis_sm_threshold_c1 = nml_soilmoisture3%jarvis_sm_threshold_c1
507
508 processmatrix(3, 2) = 22_i4
509 processmatrix(3, 3) = sum(processmatrix(1 : 3, 2))
510 call append(global_parameters, reshape(orgmattercontent_forest, (/1, ncolpars/)))
511 call append(global_parameters, reshape(orgmattercontent_impervious, (/1, ncolpars/)))
512 call append(global_parameters, reshape(orgmattercontent_pervious, (/1, ncolpars/)))
513 call append(global_parameters, reshape(ptf_lower66_5_constant, (/1, ncolpars/)))
514 call append(global_parameters, reshape(ptf_lower66_5_clay, (/1, ncolpars/)))
515 call append(global_parameters, reshape(ptf_lower66_5_db, (/1, ncolpars/)))
516 call append(global_parameters, reshape(ptf_higher66_5_constant, (/1, ncolpars/)))
517 call append(global_parameters, reshape(ptf_higher66_5_clay, (/1, ncolpars/)))
518 call append(global_parameters, reshape(ptf_higher66_5_db, (/1, ncolpars/)))
519 call append(global_parameters, reshape(ptf_ks_constant, (/1, ncolpars/)))
520 call append(global_parameters, reshape(ptf_ks_sand, (/1, ncolpars/)))
521 call append(global_parameters, reshape(ptf_ks_clay, (/1, ncolpars/)))
522 call append(global_parameters, reshape(ptf_ks_curveslope, (/1, ncolpars/)))
523 call append(global_parameters, reshape(rootfractioncoefficient_forest, (/1, ncolpars/)))
524 call append(global_parameters, reshape(rootfractioncoefficient_impervious, (/1, ncolpars/)))
525 call append(global_parameters, reshape(rootfractioncoefficient_pervious, (/1, ncolpars/)))
526 call append(global_parameters, reshape(infiltrationshapefactor, (/1, ncolpars/)))
527 call append(global_parameters, reshape(rootfractioncoefficient_sand, (/1, ncolpars/)))
528 call append(global_parameters, reshape(rootfractioncoefficient_clay, (/1, ncolpars/)))
529 call append(global_parameters, reshape(fcmin_glob, (/1, ncolpars/)))
530 call append(global_parameters, reshape(fcdelta_glob, (/1, ncolpars/)))
531 call append(global_parameters, reshape(jarvis_sm_threshold_c1, (/1, ncolpars/)))
532
533 call append(global_parameters_name, (/ &
534 'orgMatterContent_forest ', &
535 'orgMatterContent_impervious ', &
536 'orgMatterContent_pervious ', &
537 'PTF_lower66_5_constant ', &
538 'PTF_lower66_5_clay ', &
539 'PTF_lower66_5_Db ', &
540 'PTF_higher66_5_constant ', &
541 'PTF_higher66_5_clay ', &
542 'PTF_higher66_5_Db ', &
543 'PTF_Ks_constant ', &
544 'PTF_Ks_sand ', &
545 'PTF_Ks_clay ', &
546 'PTF_Ks_curveSlope ', &
547 'rootFractionCoefficient_forest ', &
548 'rootFractionCoefficient_impervious', &
549 'rootFractionCoefficient_pervious ', &
550 'infiltrationShapeFactor ', &
551 'rootFractionCoefficient_sand ', &
552 'rootFractionCoefficient_clay ', &
553 'FCmin_glob ', &
554 'FCdelta_glob ', &
555 'jarvis_sm_threshold_c1 '/))
556
557 ! check if parameter are in range
558 if (.not. in_bound(global_parameters)) &
559 call error_message('***ERROR: parameter in namelist "soilmoisture3" out of bound in ', trim(adjustl(file_namelist_param)))
560
561 ! 4- Feddes equation for ET reduction and FC dependency on root fraction coefficient
562 case(4)
563 call nml_soilmoisture4%read(file_namelist_param)
564 orgmattercontent_forest = nml_soilmoisture4%orgMatterContent_forest
565 orgmattercontent_impervious = nml_soilmoisture4%orgMatterContent_impervious
566 orgmattercontent_pervious = nml_soilmoisture4%orgMatterContent_pervious
567 ptf_lower66_5_constant = nml_soilmoisture4%PTF_lower66_5_constant
568 ptf_lower66_5_clay = nml_soilmoisture4%PTF_lower66_5_clay
569 ptf_lower66_5_db = nml_soilmoisture4%PTF_lower66_5_Db
570 ptf_higher66_5_constant = nml_soilmoisture4%PTF_higher66_5_constant
571 ptf_higher66_5_clay = nml_soilmoisture4%PTF_higher66_5_clay
572 ptf_higher66_5_db = nml_soilmoisture4%PTF_higher66_5_Db
573 ptf_ks_constant = nml_soilmoisture4%PTF_Ks_constant
574 ptf_ks_sand = nml_soilmoisture4%PTF_Ks_sand
575 ptf_ks_clay = nml_soilmoisture4%PTF_Ks_clay
576 ptf_ks_curveslope = nml_soilmoisture4%PTF_Ks_curveSlope
577 rootfractioncoefficient_forest = nml_soilmoisture4%rootFractionCoefficient_forest
578 rootfractioncoefficient_impervious = nml_soilmoisture4%rootFractionCoefficient_impervious
579 rootfractioncoefficient_pervious = nml_soilmoisture4%rootFractionCoefficient_pervious
580 infiltrationshapefactor = nml_soilmoisture4%infiltrationShapeFactor
581 rootfractioncoefficient_sand = nml_soilmoisture4%rootFractionCoefficient_sand
582 rootfractioncoefficient_clay = nml_soilmoisture4%rootFractionCoefficient_clay
583 fcmin_glob = nml_soilmoisture4%FCmin_glob
584 fcdelta_glob = nml_soilmoisture4%FCdelta_glob
585
586 processmatrix(3, 2) = 21_i4
587 processmatrix(3, 3) = sum(processmatrix(1 : 3, 2))
588 call append(global_parameters, reshape(orgmattercontent_forest, (/1, ncolpars/)))
589 call append(global_parameters, reshape(orgmattercontent_impervious, (/1, ncolpars/)))
590 call append(global_parameters, reshape(orgmattercontent_pervious, (/1, ncolpars/)))
591 call append(global_parameters, reshape(ptf_lower66_5_constant, (/1, ncolpars/)))
592 call append(global_parameters, reshape(ptf_lower66_5_clay, (/1, ncolpars/)))
593 call append(global_parameters, reshape(ptf_lower66_5_db, (/1, ncolpars/)))
594 call append(global_parameters, reshape(ptf_higher66_5_constant, (/1, ncolpars/)))
595 call append(global_parameters, reshape(ptf_higher66_5_clay, (/1, ncolpars/)))
596 call append(global_parameters, reshape(ptf_higher66_5_db, (/1, ncolpars/)))
597 call append(global_parameters, reshape(ptf_ks_constant, (/1, ncolpars/)))
598 call append(global_parameters, reshape(ptf_ks_sand, (/1, ncolpars/)))
599 call append(global_parameters, reshape(ptf_ks_clay, (/1, ncolpars/)))
600 call append(global_parameters, reshape(ptf_ks_curveslope, (/1, ncolpars/)))
601 call append(global_parameters, reshape(rootfractioncoefficient_forest, (/1, ncolpars/)))
602 call append(global_parameters, reshape(rootfractioncoefficient_impervious, (/1, ncolpars/)))
603 call append(global_parameters, reshape(rootfractioncoefficient_pervious, (/1, ncolpars/)))
604 call append(global_parameters, reshape(infiltrationshapefactor, (/1, ncolpars/)))
605 call append(global_parameters, reshape(rootfractioncoefficient_sand, (/1, ncolpars/)))
606 call append(global_parameters, reshape(rootfractioncoefficient_clay, (/1, ncolpars/)))
607 call append(global_parameters, reshape(fcmin_glob, (/1, ncolpars/)))
608 call append(global_parameters, reshape(fcdelta_glob, (/1, ncolpars/)))
609
610 call append(global_parameters_name, (/ &
611 'orgMatterContent_forest ', &
612 'orgMatterContent_impervious ', &
613 'orgMatterContent_pervious ', &
614 'PTF_lower66_5_constant ', &
615 'PTF_lower66_5_clay ', &
616 'PTF_lower66_5_Db ', &
617 'PTF_higher66_5_constant ', &
618 'PTF_higher66_5_clay ', &
619 'PTF_higher66_5_Db ', &
620 'PTF_Ks_constant ', &
621 'PTF_Ks_sand ', &
622 'PTF_Ks_clay ', &
623 'PTF_Ks_curveSlope ', &
624 'rootFractionCoefficient_forest ', &
625 'rootFractionCoefficient_impervious', &
626 'rootFractionCoefficient_pervious ', &
627 'infiltrationShapeFactor ', &
628 'rootFractionCoefficient_sand ', &
629 'rootFractionCoefficient_clay ', &
630 'FCmin_glob ', &
631 'FCdelta_glob '/))
632
633 ! check if parameter are in range
634 if (.not. in_bound(global_parameters)) &
635 call error_message('***ERROR: parameter in namelist "soilmoisture4" out of bound in ', trim(adjustl(file_namelist_param)))
636
637 case DEFAULT
638 call error_message('***ERROR: Process description for process "soilmoisture" does not exist!')
639 call error_message('***ERROR: Process description for process "soilmoisture" does not exist!')
640 end select
641
642 ! Process 4 - sealed area directRunoff
643 select case (processmatrix(4, 1))
644 ! 1 - bucket exceedance approach
645 case(1)
646 call nml_directrunoff1%read(file_namelist_param)
647 imperviousstoragecapacity = nml_directrunoff1%imperviousStorageCapacity
648
649 processmatrix(4, 2) = 1_i4
650 processmatrix(4, 3) = sum(processmatrix(1 : 4, 2))
651 call append(global_parameters, reshape(imperviousstoragecapacity, (/1, ncolpars/)))
652
653 call append(global_parameters_name, (/'imperviousStorageCapacity'/))
654
655 ! check if parameter are in range
656 if (.not. in_bound(global_parameters)) &
657 call error_message('***ERROR: parameter in namelist "directRunoff1" out of bound in ', trim(adjustl(file_namelist_param)))
658
659 case DEFAULT
660 call error_message('***ERROR: Process description for process "directRunoff" does not exist!')
661 call error_message('***ERROR: Process description for process "directRunoff" does not exist!')
662 end select
663
664 ! Process 5 - potential evapotranspiration (PET)
665 select case (processmatrix(5, 1))
666 case(-1) ! 0 - PET is input, correct PET by LAI
667 call nml_petminus1%read(file_namelist_param)
668 pet_a_forest = nml_petminus1%PET_a_forest
669 pet_a_impervious = nml_petminus1%PET_a_impervious
670 pet_a_pervious = nml_petminus1%PET_a_pervious
671 pet_b = nml_petminus1%PET_b
672 pet_c = nml_petminus1%PET_c
673
674 processmatrix(5, 2) = 5_i4
675 processmatrix(5, 3) = sum(processmatrix(1 : 5, 2))
676 call append(global_parameters, reshape(pet_a_forest, (/1, ncolpars/)))
677 call append(global_parameters, reshape(pet_a_impervious, (/1, ncolpars/)))
678 call append(global_parameters, reshape(pet_a_pervious, (/1, ncolpars/)))
679 call append(global_parameters, reshape(pet_b, (/1, ncolpars/)))
680 call append(global_parameters, reshape(pet_c, (/1, ncolpars/)))
681
682 call append(global_parameters_name, (/ &
683 'PET_a_forest ', &
684 'PET_a_impervious ', &
685 'PET_a_pervious ', &
686 'PET_b ', &
687 'PET_c '/))
688
689 ! check if parameter are in range
690 if (.not. in_bound(global_parameters)) &
691 call error_message('***ERROR: parameter in namelist "PETminus1" out of bound n ', trim(adjustl(file_namelist_param)))
692
693 case(0) ! 0 - PET is input, correct PET by aspect
694 call nml_pet0%read(file_namelist_param)
695 mincorrectionfactorpet = nml_pet0%minCorrectionFactorPET
696 maxcorrectionfactorpet = nml_pet0%maxCorrectionFactorPET
697 aspecttresholdpet = nml_pet0%aspectTresholdPET
698
699 processmatrix(5, 2) = 3_i4
700 processmatrix(5, 3) = sum(processmatrix(1 : 5, 2))
701 call append(global_parameters, reshape(mincorrectionfactorpet, (/1, ncolpars/)))
702 call append(global_parameters, reshape(maxcorrectionfactorpet, (/1, ncolpars/)))
703 call append(global_parameters, reshape(aspecttresholdpet, (/1, ncolpars/)))
704
705 call append(global_parameters_name, (/ &
706 'minCorrectionFactorPET ', &
707 'maxCorrectionFactorPET ', &
708 'aspectTresholdPET '/))
709
710 ! check if parameter are in range
711 if (.not. in_bound(global_parameters)) &
712 call error_message('***ERROR: parameter in namelist "PET0" out of bound in ', trim(adjustl(file_namelist_param)))
713
714 case(1) ! 1 - Hargreaves-Samani method (HarSam) - additional input needed: Tmin, Tmax
715 call nml_pet1%read(file_namelist_param)
716 mincorrectionfactorpet = nml_pet1%minCorrectionFactorPET
717 maxcorrectionfactorpet = nml_pet1%maxCorrectionFactorPET
718 aspecttresholdpet = nml_pet1%aspectTresholdPET
719 hargreavessamanicoeff = nml_pet1%HargreavesSamaniCoeff
720
721 processmatrix(5, 2) = 4_i4
722 processmatrix(5, 3) = sum(processmatrix(1 : 5, 2))
723 call append(global_parameters, reshape(mincorrectionfactorpet, (/1, ncolpars/)))
724 call append(global_parameters, reshape(maxcorrectionfactorpet, (/1, ncolpars/)))
725 call append(global_parameters, reshape(aspecttresholdpet, (/1, ncolpars/)))
726 call append(global_parameters, reshape(hargreavessamanicoeff, (/1, ncolpars/)))
727 call append(global_parameters_name, (/ &
728 'minCorrectionFactorPET', &
729 'maxCorrectionFactorPET', &
730 'aspectTresholdPET ', &
731 'HargreavesSamaniCoeff '/))
732
733 ! check if parameter are in range
734 if (.not. in_bound(global_parameters)) &
735 call error_message('***ERROR: parameter in namelist "PET1" out of bound in ', trim(adjustl(file_namelist_param)))
736
737 case(2) ! 2 - Priestley-Taylor method (PrieTay) - additional input needed: net_rad
738 call nml_pet2%read(file_namelist_param)
739 priestleytaylorcoeff = nml_pet2%PriestleyTaylorCoeff
740 priestleytaylorlaicorr = nml_pet2%PriestleyTaylorLAIcorr
741
742 processmatrix(5, 2) = 2_i4
743 processmatrix(5, 3) = sum(processmatrix(1 : 5, 2))
744 call append(global_parameters, reshape(priestleytaylorcoeff, (/1, ncolpars/)))
745 call append(global_parameters, reshape(priestleytaylorlaicorr, (/1, ncolpars/)))
746 call append(global_parameters_name, (/ &
747 'PriestleyTaylorCoeff ', &
748 'PriestleyTaylorLAIcorr'/))
749
750 ! check if parameter are in range
751 if (.not. in_bound(global_parameters)) &
752 call error_message('***ERROR: parameter in namelist "PET2" out of bound in ', trim(adjustl(file_namelist_param)))
753
754 case(3) ! 3 - Penman-Monteith method - additional input needed: net_rad, abs. vapour pressue, windspeed
755 call nml_pet3%read(file_namelist_param)
756 canopyheigth_forest = nml_pet3%canopyheigth_forest
757 canopyheigth_impervious = nml_pet3%canopyheigth_impervious
758 canopyheigth_pervious = nml_pet3%canopyheigth_pervious
759 displacementheight_coeff = nml_pet3%displacementheight_coeff
760 roughnesslength_momentum_coeff = nml_pet3%roughnesslength_momentum_coeff
761 roughnesslength_heat_coeff = nml_pet3%roughnesslength_heat_coeff
762 stomatal_resistance = nml_pet3%stomatal_resistance
763
764 processmatrix(5, 2) = 7_i4
765 processmatrix(5, 3) = sum(processmatrix(1 : 5, 2))
766
767 call append(global_parameters, reshape(canopyheigth_forest, (/1, ncolpars/)))
768 call append(global_parameters, reshape(canopyheigth_impervious, (/1, ncolpars/)))
769 call append(global_parameters, reshape(canopyheigth_pervious, (/1, ncolpars/)))
770 call append(global_parameters, reshape(displacementheight_coeff, (/1, ncolpars/)))
771 call append(global_parameters, reshape(roughnesslength_momentum_coeff, (/1, ncolpars/)))
772 call append(global_parameters, reshape(roughnesslength_heat_coeff, (/1, ncolpars/)))
773 call append(global_parameters, reshape(stomatal_resistance, (/1, ncolpars/)))
774
775 call append(global_parameters_name, (/ &
776 'canopyheigth_forest ', &
777 'canopyheigth_impervious ', &
778 'canopyheigth_pervious ', &
779 'displacementheight_coeff ', &
780 'roughnesslength_momentum_coeff', &
781 'roughnesslength_heat_coeff ', &
782 'stomatal_resistance '/))
783
784 ! check if parameter are in range
785 if (.not. in_bound(global_parameters)) &
786 call error_message('***ERROR: parameter in namelist "PET3" out of bound in ', trim(adjustl(file_namelist_param)))
787
788 case DEFAULT
789 call error_message('***ERROR: Process description for process "actualET" does not exist!')
790 call error_message('***ERROR: Process description for process "actualET" does not exist!')
791 end select
792
793
794 ! Process 6 - interflow
795 select case (processmatrix(6, 1))
796 ! 1 - parallel soil reservoir approach
797 case(1)
798 call nml_interflow1%read(file_namelist_param)
799 interflowstoragecapacityfactor = nml_interflow1%interflowStorageCapacityFactor
800 interflowrecession_slope = nml_interflow1%interflowRecession_slope
801 fastinterflowrecession_forest = nml_interflow1%fastInterflowRecession_forest
802 slowinterflowrecession_ks = nml_interflow1%slowInterflowRecession_Ks
803 exponentslowinterflow = nml_interflow1%exponentSlowInterflow
804
805 processmatrix(6, 2) = 5_i4
806 processmatrix(6, 3) = sum(processmatrix(1 : 6, 2))
807 call append(global_parameters, reshape(interflowstoragecapacityfactor, (/1, ncolpars/)))
808 call append(global_parameters, reshape(interflowrecession_slope, (/1, ncolpars/)))
809 call append(global_parameters, reshape(fastinterflowrecession_forest, (/1, ncolpars/)))
810 call append(global_parameters, reshape(slowinterflowrecession_ks, (/1, ncolpars/)))
811 call append(global_parameters, reshape(exponentslowinterflow, (/1, ncolpars/)))
812
813 call append(global_parameters_name, (/ &
814 'interflowStorageCapacityFactor', &
815 'interflowRecession_slope ', &
816 'fastInterflowRecession_forest ', &
817 'slowInterflowRecession_Ks ', &
818 'exponentSlowInterflow '/))
819
820 ! check if parameter are in range
821 if (.not. in_bound(global_parameters)) &
822 call error_message('***ERROR: parameter in namelist "interflow1" out of bound in ', trim(adjustl(file_namelist_param)))
823
824 case DEFAULT
825 call error_message('***ERROR: Process description for process "interflow" does not exist!')
826 call error_message('***ERROR: Process description for process "interflow" does not exist!')
827 end select
828
829 ! Process 7 - percolation
830 select case (processmatrix(7, 1))
831 ! 1 - GW layer is assumed as bucket
832 case(1)
833 call nml_percolation1%read(file_namelist_param)
834 rechargecoefficient = nml_percolation1%rechargeCoefficient
835 rechargefactor_karstic = nml_percolation1%rechargeFactor_karstic
836 gain_loss_gwreservoir_karstic = nml_percolation1%gain_loss_GWreservoir_karstic
837
838 processmatrix(7, 2) = 3_i4
839 processmatrix(7, 3) = sum(processmatrix(1 : 7, 2))
840 call append(global_parameters, reshape(rechargecoefficient, (/1, ncolpars/)))
841 call append(global_parameters, reshape(rechargefactor_karstic, (/1, ncolpars/)))
842 call append(global_parameters, reshape(gain_loss_gwreservoir_karstic, (/1, ncolpars/)))
843
844 call append(global_parameters_name, (/ &
845 'rechargeCoefficient ', &
846 'rechargeFactor_karstic ', &
847 'gain_loss_GWreservoir_karstic'/))
848
849 ! check if parameter are in range
850 if (.not. in_bound(global_parameters)) &
851 call error_message('***ERROR: parameter in namelist "percolation1" out of bound in ', trim(adjustl(file_namelist_param)))
852
853 case DEFAULT
854 call error_message('***ERROR: Process description for process "percolation" does not exist!')
855 call error_message('***ERROR: Process description for process "percolation" does not exist!')
856 end select
857
858 ! Process 8 - routing
859 select case (processmatrix(8, 1))
860 case(0)
861 ! 0 - deactivated
862 call message()
863 call message('***CAUTION: Routing is deativated! ')
864
865 processmatrix(8, 2) = 0_i4
866 processmatrix(8, 3) = sum(processmatrix(1 : 8, 2))
867 case(1)
868 ! parameter values and names are set in mRM
869 ! 1 - Muskingum approach
870 processmatrix(8, 2) = 5_i4
871 processmatrix(8, 3) = sum(processmatrix(1 : 8, 2))
872 call append(global_parameters, dummy_2d_dp)
873 call append(global_parameters_name, (/'dummy', 'dummy', 'dummy', 'dummy', 'dummy'/))
874 case(2)
875 processmatrix(8, 2) = 1_i4
876 processmatrix(8, 3) = sum(processmatrix(1 : 8, 2))
877 call append(global_parameters, dummy_2d_dp_2)
878 call append(global_parameters_name, (/'dummy'/))
879 case(3)
880 processmatrix(8, 2) = 1_i4
881 processmatrix(8, 3) = sum(processmatrix(1 : 8, 2))
882 call append(global_parameters, dummy_2d_dp_2)
883 call append(global_parameters_name, (/'dummy'/))
884 case DEFAULT
885 call error_message('***ERROR: Process description for process "routing" does not exist!')
886 end select
887
888 !===============================================================
889 ! Geological formations
890 !===============================================================
891 dummy = dummy // '' ! only to avoid warning
892
893 ! Process 9 - geoparameter
894 select case (processmatrix(9, 1))
895 case(1)
896 ! read in global parameters (NOT REGIONALIZED, i.e. these are <beta> and not <gamma>) for each geological formation used
897 call nml_geoparameter%read(file_namelist_param)
898 geoparam = nml_geoparameter%GeoParam
899
900 ! search number of geological parameters
901 do ii = 1, size(geoparam, 1) ! no while loop to avoid risk of endless loop
902 if (eq(geoparam(ii, 1), nodata_dp)) then
903 ngeounits = ii - 1
904 exit
905 end if
906 end do
907
908 ! for geology parameters
910 processmatrix(9, 3) = sum(processmatrix(1 : 9, 2))
911
912 call append(global_parameters, geoparam(1 : ngeounits, :))
913
914 ! create names
915 do ii = 1, ngeounits
916 dummy = 'GeoParam(' // trim(adjustl(num2str(ii))) // ',:)'
917 call append(global_parameters_name, (/ trim(dummy) /))
918 end do
919
920 ! check if parameter are in range
921 if (.not. in_bound(global_parameters)) &
922 call error_message('***ERROR: parameter in namelist "geoparameter" out of bound in ', trim(adjustl(file_namelist_param)))
923
924 case DEFAULT
925 call error_message('***ERROR: Process description for process "geoparameter" does not exist!')
926 end select
927
928 !===============================================================
929 ! NEUTRON COUNT
930 !===============================================================
931 ! Process 10 - neutrons
932 ! 0 - deactivated
933 ! 1 - inverse N0 based on Desilets et al. 2010
934 ! 2 - COSMIC forward operator by Shuttlworth et al. 2013
935 select case (processmatrix(10, 1))
936 case(0)
937 ! 0 - deactivated
938 call message()
939 call message('***SELECTION: Neutron count routine is deativated! ')
940
941
942 case(1)
943 ! 1 - inverse N0 based on Desilets et al. 2010
944 call nml_neutrons1%read(file_namelist_param)
945 desilets_n0 = nml_neutrons1%Desilets_N0
946 desilets_lw0 = nml_neutrons1%Desilets_LW0
947 desilets_lw1 = nml_neutrons1%Desilets_LW1
948
949 processmatrix(10,2) = 3_i4
950 processmatrix(10,3) = sum(processmatrix(1:10, 2))
951 call append(global_parameters, reshape(desilets_n0, (/1, ncolpars/)))
952 call append(global_parameters, reshape(desilets_lw0, (/1, ncolpars/)))
953 call append(global_parameters, reshape(desilets_lw1, (/1, ncolpars/)))
954
955 call append(global_parameters_name, (/ &
956 'Desilets_N0 ', &
957 'Desilets_LW0 ', &
958 'Desilets_LW1 '/))
959
960 ! check if parameter are in range
961 if (.not. in_bound(global_parameters)) &
962 call error_message('***ERROR: parameter in namelist "neutrons1" out of bound in ', trim(adjustl(file_namelist_param)))
963
964 case(2)
965 ! 2 - COSMIC version
966 call nml_neutrons2%read(file_namelist_param)
967 cosmic_n0 = nml_neutrons2%COSMIC_N0
968 cosmic_n1 = nml_neutrons2%COSMIC_N1
969 cosmic_n2 = nml_neutrons2%COSMIC_N2
970 cosmic_alpha0 = nml_neutrons2%COSMIC_alpha0
971 cosmic_alpha1 = nml_neutrons2%COSMIC_alpha1
972 cosmic_l30 = nml_neutrons2%COSMIC_L30
973 cosmic_l31 = nml_neutrons2%COSMIC_L31
974 cosmic_lw0 = nml_neutrons2%COSMIC_LW0
975 cosmic_lw1 = nml_neutrons2%COSMIC_LW1
976
977 processmatrix(10,2) = 9_i4
978 processmatrix(10,3) = sum(processmatrix(1:10, 2))
979 call append(global_parameters, reshape(cosmic_n0, (/1, ncolpars/)))
980 call append(global_parameters, reshape(cosmic_n1, (/1, ncolpars/)))
981 call append(global_parameters, reshape(cosmic_n2, (/1, ncolpars/)))
982 call append(global_parameters, reshape(cosmic_alpha0, (/1, ncolpars/)))
983 call append(global_parameters, reshape(cosmic_alpha1, (/1, ncolpars/)))
984 call append(global_parameters, reshape(cosmic_l30, (/1, ncolpars/)))
985 call append(global_parameters, reshape(cosmic_l31, (/1, ncolpars/)))
986 call append(global_parameters, reshape(cosmic_lw0, (/1, ncolpars/)))
987 call append(global_parameters, reshape(cosmic_lw1, (/1, ncolpars/)))
988
989 call append(global_parameters_name, (/ &
990 'COSMIC_N0 ', &
991 'COSMIC_N1 ', &
992 'COSMIC_N2 ', &
993 'COSMIC_alpha0 ', &
994 'COSMIC_alpha1 ', &
995 'COSMIC_L30 ', &
996 'COSMIC_L31 ', &
997 'COSMIC_LW0 ', &
998 'COSMIC_LW1 '/))
999 ! check if parameter are in range
1000 if (.not. in_bound(global_parameters)) &
1001 call error_message('***ERROR: parameter in namelist "neutrons2" out of bound in ', trim(adjustl(file_namelist_param)))
1002
1003 case DEFAULT
1004 call error_message('***ERROR: Process description for process "NEUTRON count" does not exist!')
1005 end select
1006
1007 end subroutine mpr_read_config
1008
1009end module mo_mpr_read_config
Provides constants commonly used by mHM, mRM and MPR.
integer(i4), parameter, public ncolpars
real(dp), parameter, public eps_dp
epsilon(1.0) in double precision
integer(i4), parameter, public maxnodomains
real(dp), parameter, public nodata_dp
Provides small utility functions used by multiple parts of the code (mHM, mRM, MPR)
logical function, public in_bound(params)
TODO: add description.
Provides structures needed by mHM, mRM and/or mpr.
real(dp), dimension(:, :), allocatable, target, public global_parameters
character(256), dimension(:), allocatable, public global_parameters_name
type(domain_meta), public domainmeta
integer(i4), dimension(nprocesses, 3), public processmatrix
Provides MPR specific constants.
integer(i4), parameter, public maxgeounit
integer(i4), parameter, public maxnosoilhorizons
Global variables for mpr only.
character(256), dimension(:), allocatable, public dirgridded_lai
integer(i4), public nsoilhorizons_mhm
character(256), public inputformat_gridded_lai
real(dp), dimension(:), allocatable, public horizondepth_mhm
integer(i4), public timestep_lai_input
subroutine, public mpr_read_config(file_namelist, file_namelist_param)
Read the general config of mpr.
Module containing all namelists representations.
type(nml_soilmoisture3_t), public nml_soilmoisture3
'soilmoisture3' namelist content
type(nml_percolation1_t), public nml_percolation1
'percolation1' namelist content
type(nml_pet1_t), public nml_pet1
'pet1' namelist content
type(nml_lcover_mpr_t), public nml_lcover_mpr
'lcover_mpr' namelist content
type(nml_interflow1_t), public nml_interflow1
'interflow1' namelist content
type(nml_snow1_t), public nml_snow1
'snow1' namelist content
type(nml_petminus1_t), public nml_petminus1
'petminus1' namelist content
type(nml_pet0_t), public nml_pet0
'pet0' namelist content
type(nml_pet2_t), public nml_pet2
'pet2' namelist content
type(nml_neutrons1_t), public nml_neutrons1
'neutrons1' namelist content
type(nml_soilmoisture1_t), public nml_soilmoisture1
'soilmoisture1' namelist content
type(nml_lai_data_information_t), public nml_lai_data_information
'lai_data_information' namelist content
type(nml_neutrons2_t), public nml_neutrons2
'neutrons2' namelist content
type(nml_directories_mpr_t), public nml_directories_mpr
'directories_mpr' namelist content
type(nml_soilmoisture2_t), public nml_soilmoisture2
'soilmoisture2' namelist content
type(nml_pet3_t), public nml_pet3
'pet3' namelist content
type(nml_soildata_t), public nml_soildata
'soildata' namelist content
type(nml_geoparameter_t), public nml_geoparameter
'geoparameter' namelist content
type(nml_directrunoff1_t), public nml_directrunoff1
'directrunoff1' namelist content
type(nml_interception1_t), public nml_interception1
'interception1' namelist content
type(nml_soilmoisture4_t), public nml_soilmoisture4
'soilmoisture4' namelist content