Line data Source code
1 : !> \file mo_mrm_read_config.f90
2 : !> \brief \copybrief mo_mrm_read_config
3 : !> \details \copydetails mo_mrm_read_config
4 :
5 : !> \brief read mRM config
6 : !> \details This module contains all mRM subroutines related to reading the mRM configuration either from file or copy from mHM.
7 : !> \authors Stephan Thober
8 : !> \date Aug 2015
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_mrm
12 : module mo_mrm_read_config
13 :
14 : use mo_kind, only : i4, dp
15 : use mo_message, only: message, error_message
16 :
17 : implicit none
18 :
19 : public :: mrm_read_config
20 :
21 : contains
22 :
23 : ! ------------------------------------------------------------------
24 :
25 : ! NAME
26 : ! mrm_read_config
27 :
28 : ! PURPOSE
29 : !> \brief Read the general config of mRM
30 :
31 : !> \details Depending on the variable mrm_coupling_config, the
32 : !> mRM config is either read from mrm.nml and parameters from
33 : !> mrm_parameter.nml or copied from mHM.
34 :
35 : ! INTENT(IN)
36 : !> \param[in] "character(*) :: file_namelist, file_namelist_param"
37 : !> \param[in] "character(*) :: file_namelist, file_namelist_param"
38 : !> \param[in] "logical :: do_message" - flag for writing mHM standard messages
39 :
40 : ! HISTORY
41 : !> \authors Stephan Thober
42 :
43 : !> \date Aug 2015
44 :
45 : ! Modifications:
46 : ! Stephan Thober Sep 2015 - removed stop condition when routing resolution is smaller than hydrologic resolution
47 : ! Stephan Thober Oct 2015 - added NLoutputResults namelist, fileLatLon to directories_general namelist, and readLatLon flag
48 : ! Robert Schweppe Jun 2018 - refactoring and reformatting
49 :
50 28 : subroutine mrm_read_config(file_namelist, file_namelist_param, do_message)
51 :
52 : use mo_common_constants, only : maxNoDomains, nodata_i4
53 : use mo_common_mHM_mRM_read_config, only : common_check_resolution
54 : use mo_common_variables, only : ALMA_convention, domainMeta, processMatrix
55 : use mo_mrm_constants, only : maxNoGauges
56 : use mo_mrm_file, only : file_defOutput
57 : use mo_mrm_global_variables, only : InflowGauge, domainInfo_mRM, domain_mrm, &
58 : dirGauges, dirTotalRunoff, filenameTotalRunoff, dirBankfullRunoff, gauge, is_start, &
59 : nGaugesTotal, nGaugesLocal, nInflowGaugesTotal, outputFlxState_mrm, &
60 : timeStep_model_outputs_mrm, &
61 : varnameTotalRunoff, gw_coupling, &
62 : output_deflate_level_mrm, output_double_precision_mrm, output_time_reference_mrm, &
63 : readLatLon
64 : use mo_string_utils, only : num2str
65 : use mo_namelists, only : &
66 : nml_mainconfig_mrm, &
67 : nml_directories_mrm, &
68 : nml_evaluation_gauges, &
69 : nml_inflow_gauges, &
70 : nml_mrm_outputs
71 :
72 : implicit none
73 :
74 : character(*), intent(in) :: file_namelist, file_namelist_param
75 :
76 : ! - flag for writing mHM standard messages
77 : logical, intent(in) :: do_message
78 :
79 : integer(i4), dimension(maxNoDomains) :: NoGauges_domain
80 : integer(i4), dimension(maxNoDomains, maxNoGauges) :: Gauge_id
81 : character(256), dimension(maxNoDomains, maxNoGauges) :: Gauge_filename
82 :
83 : integer(i4), dimension(maxNoDomains) :: NoInflowGauges_domain
84 : integer(i4), dimension(maxNoDomains, maxNoGauges) :: InflowGauge_id
85 : character(256), dimension(maxNoDomains, maxNoGauges) :: InflowGauge_filename
86 : logical, dimension(maxNoDomains, maxNoGauges) :: InflowGauge_Headwater
87 :
88 : integer(i4) :: domainID, iDomain
89 :
90 : integer(i4) :: iGauge
91 :
92 : integer(i4) :: idx
93 :
94 : character(256), dimension(maxNoDomains) :: dir_Gauges
95 : character(256), dimension(maxNoDomains) :: dir_Total_Runoff
96 : character(256), dimension(maxNoDomains) :: dir_Bankfull_Runoff
97 :
98 : logical :: file_exists
99 :
100 : type(domainInfo_mRM), pointer :: domain_mrm_iDomain
101 :
102 : !===============================================================
103 : ! INITIALIZATION
104 : !===============================================================
105 14 : is_start = .True.
106 :
107 : !===============================================================
108 : ! Read namelist for mainconfig for mRM
109 : !===============================================================
110 14 : call nml_mainconfig_mrm%read(file_namelist)
111 14 : ALMA_convention = nml_mainconfig_mrm%ALMA_convention
112 14 : filenameTotalRunoff = nml_mainconfig_mrm%filenameTotalRunoff
113 14 : varnameTotalRunoff = nml_mainconfig_mrm%varnameTotalRunoff
114 14 : gw_coupling = nml_mainconfig_mrm%gw_coupling
115 :
116 : !===============================================================
117 : ! Read namelist for mainpaths
118 : !===============================================================
119 14 : call nml_directories_mrm%read(file_namelist)
120 714 : dir_Gauges = nml_directories_mrm%dir_Gauges
121 714 : dir_Total_Runoff = nml_directories_mrm%dir_Total_Runoff
122 714 : dir_Bankfull_Runoff = nml_directories_mrm%dir_Bankfull_Runoff
123 :
124 70 : allocate(dirGauges(domainMeta%nDomains), dirTotalRunoff(domainMeta%nDomains), dirBankfullRunoff(domainMeta%nDomains))
125 40 : do iDomain = 1, domainMeta%nDomains
126 26 : domainID = domainMeta%indices(iDomain)
127 26 : dirGauges(iDomain) = dir_Gauges(domainID)
128 26 : dirTotalRunoff(iDomain) = dir_Total_Runoff(domainID)
129 40 : dirBankfullRunoff(iDomain) = dir_Bankfull_Runoff(domainID)
130 : end do
131 :
132 : !===============================================================
133 : ! READ EVALUATION GAUGES
134 : !===============================================================
135 14 : call nml_evaluation_gauges%read(file_namelist)
136 14 : nGaugesTotal = nml_evaluation_gauges%nGaugesTotal
137 714 : NoGauges_domain = nml_evaluation_gauges%NoGauges_domain
138 142814 : Gauge_id = nml_evaluation_gauges%Gauge_id
139 142814 : gauge_filename = nml_evaluation_gauges%gauge_filename
140 :
141 14 : if (nGaugesTotal .GT. maxNoGauges) then
142 : call error_message('***ERROR: ', trim(file_namelist), ': Total number of evaluation gauges is restricted to', &
143 0 : num2str(maxNoGauges), raise=.false.)
144 0 : call error_message(' Error occured in namlist: evaluation_gauges')
145 : end if
146 :
147 : ! ToDo: check
148 14 : nGaugesLocal = 0
149 40 : do iDomain = 1, domainMeta%nDomains
150 26 : domainID = domainMeta%indices(iDomain)
151 40 : nGaugesLocal = nGaugesLocal + NoGauges_domain(domainID)
152 : end do
153 : ! End ToDo
154 :
155 65 : allocate(gauge%gaugeId(nGaugesLocal)) ; gauge%gaugeId = nodata_i4
156 65 : allocate(gauge%domainId(nGaugesLocal)) ; gauge%domainId = nodata_i4
157 42 : allocate(gauge%fName (nGaugesLocal))
158 14 : if (nGaugesLocal > 0) then
159 14 : gauge%fName(1) = num2str(nodata_i4)
160 : end if
161 68 : allocate(domain_mrm(domainMeta%nDomains))
162 :
163 14 : idx = 0
164 40 : do iDomain = 1, domainMeta%nDomains
165 26 : domainID = domainMeta%indices(iDomain)
166 26 : domain_mrm_iDomain => domain_mrm(iDomain)
167 : ! initialize
168 26 : domain_mrm_iDomain%nGauges = nodata_i4
169 1378 : allocate(domain_mrm_iDomain%gaugeIdList(maxval(NoGauges_domain(:))))
170 55 : domain_mrm_iDomain%gaugeIdList = nodata_i4
171 1378 : allocate(domain_mrm_iDomain%gaugeIndexList(maxval(NoGauges_domain(:))))
172 55 : domain_mrm_iDomain%gaugeIndexList = nodata_i4
173 1378 : allocate(domain_mrm_iDomain%gaugeNodeList(maxval(NoGauges_domain(:))))
174 55 : domain_mrm_iDomain%gaugeNodeList = nodata_i4
175 : ! check if NoGauges_domain has a valid value
176 26 : if (NoGauges_domain(domainID) .EQ. nodata_i4) then
177 : call error_message('***ERROR: ', trim(file_namelist), ': Number of evaluation gauges for subdomain ', &
178 0 : trim(adjustl(num2str(domainID))), ' is not defined!', raise=.false.)
179 0 : call error_message(' Error occured in namelist: evaluation_gauges')
180 : end if
181 :
182 26 : domain_mrm_iDomain%nGauges = NoGauges_domain(domainID)
183 :
184 63 : do iGauge = 1, NoGauges_domain(domainID)
185 : ! check if NoGauges_domain has a valid value
186 23 : if (Gauge_id(domainID, iGauge) .EQ. nodata_i4) then
187 : call error_message('***ERROR: ', trim(file_namelist), ': ID ', &
188 0 : trim(adjustl(num2str(Gauge_id(domainID, iGauge)))), ' of evaluation gauge ', &
189 : trim(adjustl(num2str(iGauge))), ' for subdomain ', &
190 0 : trim(adjustl(num2str(iDomain))), ' is not defined!', raise=.false.)
191 0 : call error_message(' Error occured in namelist: evaluation_gauges')
192 46 : else if (trim(gauge_filename(domainID, iGauge)) .EQ. trim(num2str(nodata_i4))) then
193 : call error_message('***ERROR: ', trim(file_namelist), ': Filename of evaluation gauge ', &
194 : trim(adjustl(num2str(iGauge))), ' for subdomain ', &
195 0 : trim(adjustl(num2str(iDomain))), ' is not defined!', raise=.false.)
196 23 : call error_message(' Error occured in namelist: evaluation_gauges')
197 : end if
198 : !
199 23 : idx = idx + 1
200 23 : gauge%domainId(idx) = iDomain
201 23 : gauge%gaugeId(idx) = Gauge_id(domainID, iGauge)
202 23 : gauge%fname(idx) = trim(dirGauges(iDomain)) // trim(gauge_filename(domainID, iGauge))
203 23 : domain_mrm_iDomain%gaugeIdList(iGauge) = Gauge_id(domainID, iGauge)
204 49 : domain_mrm_iDomain%gaugeIndexList(iGauge) = idx
205 : end do
206 : end do
207 :
208 14 : if (nGaugesLocal .NE. idx) then
209 : call error_message('***ERROR: ', trim(file_namelist), ': Total number of evaluation gauges (', &
210 : trim(adjustl(num2str(nGaugesLocal))), &
211 0 : ') different from sum of gauges in subdomains (', trim(adjustl(num2str(idx))), ')!', raise=.false.)
212 0 : call error_message(' Error occured in namelist: evaluation_gauges')
213 : end if
214 :
215 : !===============================================================
216 : ! Read inflow gauge information
217 : !===============================================================
218 :
219 14 : call nml_inflow_gauges%read(file_namelist)
220 14 : nInflowGaugesTotal = nml_inflow_gauges%nInflowGaugesTotal
221 714 : NoInflowGauges_domain = nml_inflow_gauges%NoInflowGauges_domain
222 142814 : InflowGauge_id = nml_inflow_gauges%InflowGauge_id
223 142814 : InflowGauge_filename = nml_inflow_gauges%InflowGauge_filename
224 142814 : InflowGauge_Headwater = nml_inflow_gauges%InflowGauge_Headwater
225 :
226 14 : if (nInflowGaugesTotal .GT. maxNoGauges) then
227 : call error_message('***ERROR: ', trim(file_namelist), &
228 0 : ':read_gauge_lut: Total number of inflow gauges is restricted to', num2str(maxNoGauges), raise=.false.)
229 0 : call error_message(' Error occured in namlist: inflow_gauges')
230 : end if
231 :
232 : ! allocation - max() to avoid allocation with zero, needed for mhm call
233 42 : allocate(InflowGauge%gaugeId (max(1, nInflowGaugesTotal)))
234 28 : allocate(InflowGauge%domainId (max(1, nInflowGaugesTotal)))
235 42 : allocate(InflowGauge%fName (max(1, nInflowGaugesTotal)))
236 : ! dummy initialization
237 29 : InflowGauge%gaugeId = nodata_i4
238 29 : InflowGauge%domainId = nodata_i4
239 29 : InflowGauge%fName = num2str(nodata_i4)
240 :
241 14 : idx = 0
242 40 : do iDomain = 1, domainMeta%nDomains
243 26 : domainID = domainMeta%indices(iDomain)
244 26 : domain_mrm_iDomain => domain_mrm(iDomain)
245 :
246 1378 : allocate(domain_mrm_iDomain%InflowGaugeIdList (max(1, maxval(NoInflowGauges_domain(:)))))
247 1378 : allocate(domain_mrm_iDomain%InflowGaugeHeadwater (max(1, maxval(NoInflowGauges_domain(:)))))
248 1378 : allocate(domain_mrm_iDomain%InflowGaugeIndexList (max(1, maxval(NoInflowGauges_domain(:)))))
249 1378 : allocate(domain_mrm_iDomain%InflowGaugeNodeList (max(1, maxval(NoInflowGauges_domain(:)))))
250 : ! dummy initialization
251 26 : domain_mrm_iDomain%nInflowGauges = 0
252 52 : domain_mrm_iDomain%InflowGaugeIdList = nodata_i4
253 52 : domain_mrm_iDomain%InflowGaugeHeadwater = .FALSE.
254 52 : domain_mrm_iDomain%InflowGaugeIndexList = nodata_i4
255 52 : domain_mrm_iDomain%InflowGaugeNodeList = nodata_i4
256 : ! no inflow gauge for subdomain i
257 26 : if (NoInflowGauges_domain(domainID) .EQ. nodata_i4) then
258 0 : NoInflowGauges_domain(domainID) = 0
259 : end if
260 :
261 26 : domain_mrm_iDomain%nInflowGauges = NoInflowGauges_domain(domainID)
262 :
263 42 : do iGauge = 1, NoInflowGauges_domain(domainID)
264 : ! check if NoInflowGauges_domain has a valid value
265 2 : if (InflowGauge_id(domainID, iGauge) .EQ. nodata_i4) then
266 : call error_message('***ERROR: ', trim(file_namelist), ':ID of inflow gauge ', &
267 : trim(adjustl(num2str(iGauge))), ' for subdomain ', &
268 0 : trim(adjustl(num2str(iDomain))), ' is not defined!', raise=.false.)
269 0 : call error_message(' Error occured in namlist: inflow_gauges')
270 4 : else if (trim(InflowGauge_filename(domainID, iGauge)) .EQ. trim(num2str(nodata_i4))) then
271 : call error_message('***ERROR: ', trim(file_namelist), ':Filename of inflow gauge ', &
272 : trim(adjustl(num2str(iGauge))), ' for subdomain ', &
273 0 : trim(adjustl(num2str(iDomain))), ' is not defined!', raise=.false.)
274 2 : call error_message(' Error occured in namlist: inflow_gauges')
275 : end if
276 : !
277 2 : idx = idx + 1
278 2 : InflowGauge%domainId(idx) = iDomain
279 2 : InflowGauge%gaugeId(idx) = InflowGauge_id(domainID, iGauge)
280 2 : InflowGauge%fname(idx) = trim(dirGauges(domainID)) // trim(InflowGauge_filename(domainID, iGauge))
281 2 : domain_mrm_iDomain%InflowGaugeIdList(iGauge) = InflowGauge_id(domainID, iGauge)
282 2 : domain_mrm_iDomain%InflowGaugeHeadwater(iGauge) = InflowGauge_Headwater(domainID, iGauge)
283 28 : domain_mrm_iDomain%InflowGaugeIndexList(iGauge) = idx
284 : end do
285 : end do
286 :
287 14 : if (nInflowGaugesTotal .NE. idx) then
288 : call error_message('***ERROR: ', trim(file_namelist), ': Total number of inflow gauges (', &
289 : trim(adjustl(num2str(nInflowGaugesTotal))), &
290 0 : ') different from sum of inflow gauges in subdomains (', trim(adjustl(num2str(idx))), ')!', raise=.false.)
291 0 : call error_message(' Error occured in namlist: inflow_gauges')
292 : end if
293 :
294 14 : call common_check_resolution(do_message, .true.)
295 :
296 : !===============================================================
297 : ! Read namelist global parameters
298 : !===============================================================
299 14 : call read_mrm_routing_params(processMatrix(8, 1), file_namelist_param)
300 :
301 : !===============================================================
302 : ! Read Output specifications for mRM
303 : !===============================================================
304 14 : call nml_mrm_outputs%read(file_defOutput)
305 14 : output_deflate_level_mrm = nml_mrm_outputs%output_deflate_level_mrm
306 14 : output_double_precision_mrm = nml_mrm_outputs%output_double_precision_mrm
307 14 : output_time_reference_mrm = nml_mrm_outputs%output_time_reference_mrm
308 14 : timeStep_model_outputs_mrm = nml_mrm_outputs%timeStep_model_outputs_mrm
309 42 : outputFlxState_mrm = nml_mrm_outputs%outputFlxState_mrm
310 :
311 34 : readLatLon = any(outputFlxState_mrm)
312 :
313 34 : if (any(outputFlxState_mrm)) then
314 4 : call message('')
315 4 : call message(' Following output will be written:')
316 4 : call message(' NetCDF deflate level: ', adjustl(trim(num2str(output_deflate_level_mrm))))
317 4 : if ( output_double_precision_mrm ) then
318 4 : call message(' NetCDF output precision: double')
319 : else
320 0 : call message(' NetCDF output precision: single')
321 : end if
322 4 : select case(output_time_reference_mrm)
323 : case(0)
324 4 : call message(' NetCDF output time reference point: start of time interval')
325 : case(1)
326 0 : call message(' NetCDF output time reference point: center of time interval')
327 : case(2)
328 4 : call message(' NetCDF output time reference point: end of time interval')
329 : end select
330 4 : call message(' FLUXES:')
331 4 : if (outputFlxState_mrm(1)) then
332 4 : call message(' routed streamflow (L11_qMod) [m3 s-1]')
333 : end if
334 4 : if (outputFlxState_mrm(2)) then
335 1 : call message(' river temperature (RivTemp) [deg C]')
336 : end if
337 4 : if (gw_coupling) then
338 0 : call message(' river head (river_head) [m]')
339 : end if
340 : end if
341 :
342 14 : call message('')
343 14 : call message(' FINISHED reading config')
344 14 : call message('')
345 :
346 14 : end subroutine mrm_read_config
347 :
348 : ! ---------------------------------------------------------------------------
349 : ! SUBROUTINE READ_MRM_ROUTING_PARAMS
350 : ! ---------------------------------------------------------------------------
351 : ! NAME
352 : ! read_mrm_routing_params
353 :
354 : ! PURPOSE
355 : !> \brief TODO: add description
356 :
357 : !> \details TODO: add description
358 :
359 : ! INTENT(IN)
360 : !> \param[in] "integer(i4) :: processCase" it is the default case, should be one
361 : !> \param[in] "character(*) :: file_namelist_param" file name containing parameter namelist
362 :
363 : ! HISTORY
364 : !> \authors Robert Schweppe
365 :
366 : !> \date Jun 2018
367 :
368 : ! Modifications:
369 :
370 14 : subroutine read_mrm_routing_params(processCase, file_namelist_param)
371 :
372 14 : use mo_common_constants, only : nColPars
373 : use mo_common_functions, only : in_bound
374 : use mo_common_variables, only : global_parameters, global_parameters_name, processMatrix
375 : use mo_message, only : message
376 : use mo_namelists, only : &
377 : nml_routing1, &
378 : nml_routing2, &
379 : nml_routing3
380 :
381 : implicit none
382 :
383 : ! it is the default case, should be one
384 : integer(i4), intent(in) :: processCase
385 :
386 : ! file name containing parameter namelist
387 : character(*), intent(in) :: file_namelist_param
388 :
389 : ! equals sum of previous parameters
390 : integer(i4) :: start_index
391 :
392 84 : real(dp), dimension(nColPars) :: muskingumTravelTime_constant
393 84 : real(dp), dimension(nColPars) :: muskingumTravelTime_riverLength
394 84 : real(dp), dimension(nColPars) :: muskingumTravelTime_riverSlope
395 84 : real(dp), dimension(nColPars) :: muskingumTravelTime_impervious
396 84 : real(dp), dimension(nColPars) :: muskingumAttenuation_riverSlope
397 :
398 84 : real(dp), dimension(nColPars) :: streamflow_celerity
399 84 : real(dp), dimension(nColPars) :: slope_factor
400 :
401 14 : if (processCase .eq. 1_i4) then
402 9 : call nml_routing1%read(file_namelist_param)
403 54 : muskingumTravelTime_constant = nml_routing1%muskingumTravelTime_constant
404 54 : muskingumTravelTime_riverLength = nml_routing1%muskingumTravelTime_riverLength
405 54 : muskingumTravelTime_riverSlope = nml_routing1%muskingumTravelTime_riverSlope
406 54 : muskingumTravelTime_impervious = nml_routing1%muskingumTravelTime_impervious
407 54 : muskingumAttenuation_riverSlope = nml_routing1%muskingumAttenuation_riverSlope
408 5 : else if (processCase .eq. 2_i4) then
409 1 : call nml_routing2%read(file_namelist_param)
410 6 : streamflow_celerity = nml_routing2%streamflow_celerity
411 4 : else if (processCase .eq. 3_i4) then
412 3 : call nml_routing3%read(file_namelist_param)
413 18 : slope_factor = nml_routing3%slope_factor
414 : end if
415 :
416 : ! -------------------------------------------------------------------------
417 : ! INCLUDE MRM PARAMETERS IN PARAMETERS OF MHM
418 : ! -------------------------------------------------------------------------
419 : ! Muskingum routing parameters with MPR
420 14 : if (processCase .eq. 1_i4) then
421 : ! insert parameter values and names at position required by mhm
422 9 : processMatrix(8, 1) = processCase
423 9 : processMatrix(8, 2) = 5_i4
424 81 : processMatrix(8, 3) = sum(processMatrix(1 : 8, 2))
425 9 : start_index = processMatrix(8, 3) - processMatrix(8, 2)
426 54 : global_parameters(start_index + 1, :) = muskingumTravelTime_constant
427 54 : global_parameters(start_index + 2, :) = muskingumTravelTime_riverLength
428 54 : global_parameters(start_index + 3, :) = muskingumTravelTime_riverSlope
429 54 : global_parameters(start_index + 4, :) = muskingumTravelTime_impervious
430 54 : global_parameters(start_index + 5, :) = muskingumAttenuation_riverSlope
431 :
432 0 : global_parameters_name(start_index + 1 : start_index + processMatrix(8, 2)) = (/ &
433 : 'muskingumTravelTime_constant ', &
434 : 'muskingumTravelTime_riverLength', &
435 : 'muskingumTravelTime_riverSlope ', &
436 : 'muskingumTravelTime_impervious ', &
437 54 : 'muskingumAttenuation_riverSlope'/)
438 : ! adaptive timestep routing
439 5 : else if (processCase .eq. 2_i4) then
440 1 : processMatrix(8, 1) = processCase
441 1 : processMatrix(8, 2) = 1_i4
442 9 : processMatrix(8, 3) = sum(processMatrix(1 : 8, 2))
443 1 : start_index = processMatrix(8, 3) - processMatrix(8, 2)
444 6 : global_parameters(start_index + 1, :) = streamflow_celerity
445 :
446 0 : global_parameters_name(start_index + 1 : start_index + processMatrix(8, 2)) = (/ &
447 2 : 'streamflow_celerity'/)
448 : ! adaptive timestep routing - varying celerity
449 4 : else if (processCase .eq. 3_i4) then
450 : ! insert parameter values and names at position required by mhm
451 3 : processMatrix(8, 1) = processCase
452 3 : processMatrix(8, 2) = 1_i4
453 27 : processMatrix(8, 3) = sum(processMatrix(1:8, 2))
454 3 : start_index = processMatrix(8, 3) - processMatrix(8, 2)
455 18 : global_parameters(start_index + 1, :) = slope_factor
456 :
457 0 : global_parameters_name(start_index + 1 : start_index + processMatrix(8,2)) = (/ &
458 6 : 'slope_factor'/)
459 : end if
460 :
461 : ! check if parameter are in range
462 14 : if (.not. in_bound(global_parameters)) then
463 : call error_message('***ERROR: parameter in routing namelist out of bound in ', &
464 0 : trim(adjustl(file_namelist_param)))
465 : end if
466 :
467 28 : end subroutine read_mrm_routing_params
468 : end module mo_mrm_read_config
|