5.13.3-dev0
mHM
The mesoscale Hydrological Model
Loading...
Searching...
No Matches
mo_common_read_config.F90
Go to the documentation of this file.
1!> \file mo_common_read_config.f90
2!> \brief \copybrief mo_common_read_config
3!> \details \copydetails mo_common_read_config
4
5!> \brief Reading of main model configurations.
6!> \details This routine reads the configurations of namelists commonly used by mHM, mRM and MPR
7!> \authors Matthias Zink
8!> \date Dec 2012
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_common
13
14 USE mo_kind, ONLY : i4, dp
15 use mo_message, only: error_message
16
17 IMPLICIT NONE
18
19 PRIVATE
20
22
23 ! ------------------------------------------------------------------
24
25CONTAINS
26
27
28 !> \brief Read main configurations commonly used by mHM, mRM and MPR
29 !> \details Read the main configurations commonly used by mHM, mRM and MPR, namely:
30 !! project_description, directories_general, mainconfig, processSelection, LCover
31 !> \changelog
32 !! - Robert Schweppe Dec 2018
33 !! - refactoring and restructuring
34 !! - Sebastian Müller Mar 2023
35 !! - added check_L0Domain
36 !> \authors Matthias Zink
37 !> \date Dec 2012
38 subroutine common_read_config(file_namelist)
39
40 use mo_namelists, only : &
53 use mo_message, only : message
54 use mo_string_utils, only : num2str
55
56 implicit none
57
58 !> name of file
59 character(*), intent(in) :: file_namelist
60
61 ! Choosen process description number
62 integer(i4), dimension(nProcesses) :: processcase
63
64 character(256), dimension(maxNoDomains) :: dir_morpho
65
66 character(256), dimension(maxNoDomains) :: mhm_file_restartout
67
68 character(256), dimension(maxNoDomains) :: mrm_file_restartout
69
70 character(256), dimension(maxNoDomains) :: dir_lcover
71
72 character(256), dimension(maxNoDomains) :: dir_out
73
74 character(256), dimension(maxNoDomains) :: file_latlon
75
76 real(dp), dimension(maxNoDomains) :: resolution_hydrology
77
78 integer(i4), dimension(maxNoDomains) :: l0domain
79
80 integer(i4), dimension(maxNoDomains) :: read_opt_domain_data
81
82 ! starting year LCover
83 integer(i4), dimension(maxNLCovers) :: lcoveryearstart
84
85 ! ending year LCover
86 integer(i4), dimension(maxNLCovers) :: lcoveryearend
87
88 ! filename of Lcover file
89 character(256), dimension(maxNLCovers) :: lcoverfname
90
91 integer(i4) :: i, newdomainid, domainid, idomain, ndomains
92
93 ! flag to advance nuniqueL0Domain counter
94 logical :: addcounter
95
96 !===============================================================
97 ! Read namelist specifying the project description
98 !===============================================================
99 call nml_project_description%read(file_namelist)
107
108 !===============================================================
109 ! Read namelist specifying the model configuration
110 !===============================================================
111 call nml_mainconfig%read(file_namelist)
112 iflag_cordinate_sys = nml_mainconfig%iFlag_cordinate_sys
113 resolution_hydrology = nml_mainconfig%resolution_Hydrology
114 ndomains = nml_mainconfig%nDomains
115 l0domain = nml_mainconfig%L0Domain
116 write_restart = nml_mainconfig%write_restart
117 read_opt_domain_data = nml_mainconfig%read_opt_domain_data
118
119 call init_domain_variable(ndomains, read_opt_domain_data(1:ndomains), domainmeta)
120
121 if (ndomains .GT. maxnodomains) then
122 call error_message('***ERROR: Number of domains is resticted to ', trim(num2str(maxnodomains)), '!')
123 end if
124
125 call check_l0domain(l0domain, ndomains)
126
127 ! allocate patharray sizes
128 allocate(resolutionhydrology(domainmeta%nDomains))
129 allocate(dirmorpho(domainmeta%nDomains))
130 allocate(mhmfilerestartout(domainmeta%nDomains))
131 allocate(mrmfilerestartout(domainmeta%nDomains))
132 allocate(dirlcover(domainmeta%nDomains))
133 allocate(dirout(domainmeta%nDomains))
134 allocate(filelatlon(domainmeta%nDomains))
135 allocate(domainmeta%L0DataFrom(domainmeta%nDomains))
136 allocate(domainmeta%optidata(domainmeta%nDomains))
137 allocate(domainmeta%doRouting(domainmeta%nDomains))
138
139 nuniquel0domains = 0_i4
140 do idomain = 1, domainmeta%nDomains
141 domainid = domainmeta%indices(idomain)
142 resolutionhydrology(idomain) = resolution_hydrology(domainid)
143 ! if a domain uses the same L0 data as a previous one, write
144 ! the index into domainMeta%L0DataFrom
145 newdomainid = l0domain(domainid)
146 domainmeta%L0DataFrom(idomain) = idomain
147 !
148 addcounter = .true.
149 do i = 1, idomain - 1
150 if (newdomainid == domainmeta%indices(i)) then
151 domainmeta%L0DataFrom(idomain) = i
152 addcounter = .false.
153 end if
154 end do
155 if (addcounter) nuniquel0domains = nuniquel0domains + 1_i4
156 end do
157
158 ! check for possible options
159 if(.NOT. (iflag_cordinate_sys == 0 .OR. iflag_cordinate_sys == 1)) then
160 call error_message('***ERROR: coordinate system for the model run should be 0 or 1')
161 end if
162
163 !===============================================================
164 ! Read land cover
165 !===============================================================
166 call nml_lcover%read(file_namelist)
167 nlcoverscene = nml_lcover%nLcoverScene
168 lcoveryearstart = nml_lcover%LCoverYearStart
169 lcoveryearend = nml_lcover%LCoverYearEnd
170 lcoverfname = nml_lcover%LCoverfName
171
172 ! put land cover scenes to corresponding file name and LuT
173 ! this is done already here for MPR, which does not check for the time periods
174 allocate(lcfilename(nlcoverscene))
175 allocate(lc_year_start(nlcoverscene))
176 allocate(lc_year_end(nlcoverscene))
177 lcfilename(:) = lcoverfname(1 : nlcoverscene)
178 lc_year_start(:) = lcoveryearstart(1 : nlcoverscene)
179 lc_year_end(:) = lcoveryearend(1 : nlcoverscene)
180
181 !===============================================================
182 ! Read namelist for mainpaths
183 !===============================================================
184 call nml_directories_general%read(file_namelist)
187 dir_morpho = nml_directories_general%dir_Morpho
188 dir_lcover = nml_directories_general%dir_LCover
189 dir_out = nml_directories_general%dir_Out
190 mhm_file_restartout = nml_directories_general%mhm_file_RestartOut
191 mrm_file_restartout = nml_directories_general%mrm_file_RestartOut
192 file_latlon = nml_directories_general%file_LatLon
193
194 do idomain = 1, domainmeta%nDomains
195 domainid = domainmeta%indices(idomain)
196 domainmeta%optidata(idomain) = read_opt_domain_data(domainid)
197 dirmorpho(idomain) = dir_morpho(domainid)
198 mhmfilerestartout(idomain) = mhm_file_restartout(domainid)
199 mrmfilerestartout(idomain) = mrm_file_restartout(domainid)
200 dirlcover(idomain) = dir_lcover(domainid)
201 dirout(idomain) = dir_out(domainid)
202 filelatlon(idomain) = file_latlon(domainid)
203 end do
204
205 ! ToDo: add test if opti_function matches at least one domainMeta%optidata
206 ! as soon as common and common_mRM_mHM are merged, if that is the plan
207
208
209 !===============================================================
210 ! Read process selection list
211 !===============================================================
212 call nml_processselection%read(file_namelist)
213 processcase = nml_processselection%processCase
214
215 processmatrix = 0_i4
216 processmatrix(:, 1) = processcase
217 if (processmatrix(8, 1) == 0) then
218 domainmeta%doRouting(:) = .false.
219 else
220 domainmeta%doRouting(:) = .true.
221 end if
222
223 end subroutine common_read_config
224
225
226 !> \brief Set land cover scenes IDs
227 !> \changelog
228 !! - Robert Schweppe Dec 2018
229 !! - refactoring and restructuring
230 !> \authors Matthias Zink
231 !> \date Dec 2012
232 subroutine set_land_cover_scenes_id(sim_Per, LCyear_Id)
233
235 use mo_common_types, only: period
237 use mo_string_utils, only : num2str
238
239 implicit none
240
241 type(period), dimension(:), intent(in) :: sim_per !< simulation period
242 integer(i4), dimension(:, :), allocatable, intent(inout) :: lcyear_id !< land cover year ID
243
244 integer(i4) :: ii, idomain
245
246
247 ! countercheck if land cover covers simulation period
248 if (lc_year_start(1) .GT. minval(sim_per(1 : domainmeta%nDomains)%yStart)) then
249 call error_message('***ERROR: Land cover for warming period is missing!', raise=.false.)
250 call error_message(' SimStart : ', trim(num2str(minval(sim_per(1 : domainmeta%nDomains)%yStart))), raise=.false.)
251 call error_message(' LCoverStart: ', trim(num2str(lc_year_start(1))))
252 end if
253 if (lc_year_end(nlcoverscene) .LT. maxval(sim_per(1 : domainmeta%nDomains)%yEnd)) then
254 call error_message('***ERROR: Land cover period shorter than modelling period!', raise=.false.)
255 call error_message(' SimEnd : ', trim(num2str(maxval(sim_per(1 : domainmeta%nDomains)%yEnd))), raise=.false.)
256 call error_message(' LCoverEnd: ', trim(num2str(lc_year_end(nlcoverscene))))
257 end if
258 !
259 allocate(lcyear_id(minval(sim_per(1 : domainmeta%nDomains)%yStart) : maxval(sim_per(1 : domainmeta%nDomains)%yEnd), &
260 domainmeta%nDomains))
261 lcyear_id = nodata_i4
262 do idomain = 1, domainmeta%nDomains
263 do ii = 1, nlcoverscene
264 ! land cover before model period or land cover after model period
265 if ((lc_year_end(ii) .LT. sim_per(idomain)%yStart) .OR. &
266 (lc_year_start(ii) .GT. sim_per(idomain)%yEnd)) then
267 cycle
268 ! land cover period fully covers model period
269 else if ((lc_year_start(ii) .LE. sim_per(idomain)%yStart) .AND. &
270 (lc_year_end(ii) .GE. sim_per(idomain)%yEnd)) then
271 lcyear_id(sim_per(idomain)%yStart : sim_per(idomain)%yEnd, idomain) = ii
272 exit
273 ! land cover period covers beginning of model period
274 else if ((lc_year_start(ii) .LE. sim_per(idomain)%yStart) .AND. &
275 (lc_year_end(ii) .LT. sim_per(idomain)%yEnd)) then
276 lcyear_id(sim_per(idomain)%yStart : lc_year_end(ii), idomain) = ii
277 ! land cover period covers end of model period
278 else if ((lc_year_start(ii) .GT. sim_per(idomain)%yStart) .AND. &
279 (lc_year_end(ii) .GE. sim_per(idomain)%yEnd)) then
280 lcyear_id(lc_year_start(ii) : sim_per(idomain)%yEnd, idomain) = ii
281 ! land cover period covers part of model_period
282 else
283 lcyear_id(lc_year_start(ii) : lc_year_end(ii), idomain) = ii
284 end if
285 end do
286 end do
287
288
289 end subroutine set_land_cover_scenes_id
290
291
292 !> \brief Initialization of the domain variables
293 !> \details Initialization of the domain variable for all domain loops and if activated for parallelization
294 !! In case of MPI parallelization domainMeta%overAllNumberOfDomains is a
295 !! variable where the number of domains from the namelist is stored. By this
296 !! every process knows the total number of domains. Then, in a loop the
297 !! domains are distributed onto the processes. There is a master process
298 !! and several subprocesses. The master process only reads the confings in the
299 !! mHM driver.
300 !!
301 !! The subprocesses get a number of domains. domainMeta%nDomain refers
302 !! to the number of domains assigned to a specific process. It is a local
303 !! variable and therefore has a different value for each process.
304 !!
305 !! In case more domains are there than processes, currently the domains
306 !! are distributed round robin, i.e. like cards in a card game.
307 !!
308 !! In case less domains than processes exist, all remaining processes
309 !! are assigned to the routing domains round robin. In that case the
310 !! local communicator is of interest: It is a group of processes assigned
311 !! to a routing domain again with a master process
312 !! (domainMeta%isMasterInComLocal) and subprocesses. This communicator can
313 !! in future be passed to the routing parallelization.
314 !> \author Maren Kaluza
315 !> \date Sep 2019
316 subroutine init_domain_variable(nDomains, optiData, domainMeta)
318#ifdef MPI
319 use mo_common_variables, only: comm
320 use mpi_f08
321#endif
322 integer(i4), intent(in) :: nDomains !< number of domains
323 integer(i4), dimension(:), intent(in) :: optiData !< optimization data
324 type(domain_meta), intent(inout) :: domainMeta !< domain meta info
325
326 integer :: ierror
327 integer(i4) :: nproc
328 integer(i4) :: rank
329 integer(i4) :: iDomain
330 integer(i4) :: colDomain, colMasters
331
332 domainmeta%overallNumberOfDomains = ndomains
333#ifdef MPI
334 ! find number of processes nproc
335 call mpi_comm_size(comm, nproc, ierror)
336 ! find the number the process is referred to, called rank
337 call mpi_comm_rank(comm, rank, ierror)
338 if (nproc < 2) then
339 call error_message('at least 2 processes are required')
340 end if
341 ! if there are more processes than domains
342 if (nproc > domainmeta%overallNumberOfDomains + 1) then
343 domainmeta%nDomains = 0
344 ! master reads only metadata of all domains
345 if (rank == 0) then
346 call init_domain_variable_for_master(domainmeta, colmasters, coldomain)
347 ! all other nodes read metadata but also data of assigned domains
348 else
349 ! currently each domain gets one process except it is a routing domain.
350 ! in that case the remaining processes are distributed round robin to
351 ! the routing domains.
352 call distribute_processes_to_domains_according_to_role(optidata, rank, &
353 domainmeta, colmasters, coldomain)
354 end if
355 ! two communicators are created, i.e. groups of processes that talk about
356 ! a certain topic:
357 ! comMaster is the communicator of all processes that need to read all
358 ! data. These are the processes that are masters in the comLocal plus the
359 ! master over all processes. comLocal is a communicator for a group of
360 ! processes assigned to the same domain.
361 call mpi_comm_split(comm, colmasters, rank, domainmeta%comMaster, ierror)
362 call mpi_comm_split(comm, coldomain, rank, domainmeta%comLocal, ierror)
363 call mpi_comm_size(domainmeta%comMaster, nproc, ierror)
364 else
365 ! in case of more domains than processes, distribute domains round robin
366 ! onto the processes
367 call mpi_comm_dup(comm, domainmeta%comMaster, ierror)
368 domainmeta%isMasterInComLocal = .true.
369 domainmeta%nDomains = 0
370 ! master reads only metadata of all domains
371 if (rank == 0) then
372 domainmeta%nDomains = domainmeta%overallNumberOfDomains
373 allocate(domainmeta%indices(domainmeta%nDomains))
374 do idomain = 1, domainmeta%nDomains
375 domainmeta%indices(idomain) = idomain
376 end do
377 ! all other nodes read metadata but also data of assigned domains
378 else
379 call distributedomainsroundrobin(nproc, rank, domainmeta)
380 end if
381 end if ! round robin
382#else
383 domainmeta%nDomains = ndomains
384 allocate(domainmeta%indices(domainmeta%nDomains))
385 do idomain = 1, domainmeta%nDomains
386 domainmeta%indices(idomain) = idomain
387 end do
388#endif
389
390 end subroutine init_domain_variable
391
392#ifdef MPI
393 subroutine init_domain_variable_for_master(domainMeta, colMasters, colDomain)
395 type(domain_meta), intent(inout) :: domainMeta
396 integer(i4), intent(out) :: colMasters
397 integer(i4), intent(out) :: colDomain
398 !local
399 integer(i4) :: iDomain
400
401 domainmeta%nDomains = domainmeta%overallNumberOfDomains
402 allocate(domainmeta%indices(domainmeta%nDomains))
403 do idomain = 1, domainmeta%nDomains
404 domainmeta%indices(idomain) = idomain
405 end do
406 colmasters = 1
407 coldomain = 0
408 domainmeta%isMasterInComLocal = .true.
409
410 end subroutine init_domain_variable_for_master
411
412 subroutine distributedomainsroundrobin(nproc, rank, domainMeta)
414 integer(i4), intent(in) :: nproc
415 integer(i4), intent(in) :: rank
416 type(domain_meta), intent(inout) :: domainMeta
417
418 integer(i4) :: iDomain, iProcDomain
419
420 do idomain = 1 , domainmeta%overallNumberOfDomains
421 if (rank == (modulo(idomain + nproc - 2, (nproc - 1)) + 1)) then
422 domainmeta%nDomains = domainmeta%nDomains + 1
423 end if
424 end do
425 allocate(domainmeta%indices(domainmeta%nDomains))
426 iprocdomain = 0
427 do idomain = 1 , domainmeta%overallNumberOfDomains
428 if (rank == (modulo(idomain + nproc - 2, (nproc - 1)) + 1)) then
429 iprocdomain = iprocdomain + 1
430 domainmeta%indices(iprocdomain) = idomain
431 end if
432 end do
433 end subroutine distributedomainsroundrobin
434
435 subroutine distribute_processes_to_domains_according_to_role(optiData, rank, &
436 domainMeta, colMasters, colDomain)
438 integer(i4), dimension(:), intent(in) :: optiData
439 integer(i4), intent(in) :: rank
440 type(domain_meta), intent(inout) :: domainMeta
441 integer(i4), intent(out) :: colMasters
442 integer(i4), intent(out) :: colDomain
443
444 ! local
445 integer(i4) :: nDomainsAll, nTreeDomains, i, iDomain
446 integer(i4), dimension(:), allocatable :: treeDomainList
447
448 ndomainsall = domainmeta%overallNumberOfDomains
449 ntreedomains = 0
450 do idomain = 1, ndomainsall
451 !ToDo: should also routing but not opti domains be counted?
452 if (optidata(idomain) == 1) then
453 ntreedomains = ntreedomains + 1
454 end if
455 end do
456 allocate(treedomainlist(ntreedomains))
457 i = 0
458 do idomain = 1, ndomainsall
459 if (optidata(idomain) == 1) then
460 i = i + 1
461 treedomainlist(i) = idomain
462 end if
463 end do
464 if (rank < ndomainsall + 1) then
465 colmasters = 1
466 coldomain = rank
467 domainmeta%isMasterInComLocal = .true.
468 domainmeta%nDomains = 1
469 allocate(domainmeta%indices(domainmeta%nDomains))
470 domainmeta%indices(1) = rank
471 else
472 colmasters = 0
473 if (ntreedomains > 0) then
474 coldomain = treedomainlist(mod(rank, ntreedomains) + 1)
475 else
476 coldomain = 1
477 end if
478 domainmeta%isMasterInComLocal = .false.
479 domainmeta%nDomains = 1
480 allocate(domainmeta%indices(domainmeta%nDomains))
481 ! ToDo : temporary solution, this should either not read data at all
482 ! or data corresponding to the master process
483 domainmeta%indices(1) = 1
484 end if
485 deallocate(treedomainlist)
486 end subroutine
487#endif
488
489 ! \brief check the L0Domain variable from the namelist
490 !> \authors Sebastian Müller
491 !> \date Mar 2023
492 subroutine check_l0domain(L0Domain, nDomains)
494 use mo_string_utils, only : num2str
495
496 integer(i4), dimension(maxNoDomains), intent(in) ::L0Domain !< given L0Domain variable
497 integer(i4), intent(in) :: nDomains !< number of domains
498
499 integer(i4) :: i
500
501 do i = 1, ndomains
502 if (l0domain(i) < 0) call error_message( &
503 "L0Domain values need to be positive: ", &
504 "L0Domain(", trim(adjustl(num2str(i))), ") = ", trim(adjustl(num2str(l0domain(i)))))
505 if (l0domain(i) > i) call error_message( &
506 "L0Domain values need to be less or equal to the domain index: ", &
507 "L0Domain(", trim(adjustl(num2str(i))), ") = ", trim(adjustl(num2str(l0domain(i)))))
508 ! check for increasing values
509 if (i > 1) then
510 if (l0domain(i) < l0domain(i-1)) call error_message( &
511 "L0Domain values need to be increasing: ", &
512 "L0Domain(", trim(adjustl(num2str(i-1))), ") = ", trim(adjustl(num2str(l0domain(i-1)))), &
513 ", L0Domain(", trim(adjustl(num2str(i))), ") = ", trim(adjustl(num2str(l0domain(i)))))
514 end if
515 ! if lower, check that the reference domain uses its own L0Data
516 if (l0domain(i) < i) then
517 if (l0domain(l0domain(i)) /= l0domain(i)) call error_message( &
518 "L0Domain values should be taken from a domain with its own L0 data: ", &
519 "L0Domain(", trim(adjustl(num2str(i))), ") = ", trim(adjustl(num2str(l0domain(i)))), &
520 ", L0Domain(", trim(adjustl(num2str(l0domain(i)))), ") = ", trim(adjustl(num2str(l0domain(l0domain(i))))))
521 end if
522 end do
523
524 end subroutine check_l0domain
525
526END MODULE mo_common_read_config
Provides constants commonly used by mHM, mRM and MPR.
integer(i4), parameter, public maxnodomains
integer(i4), parameter, public maxnlcovers
integer(i4), parameter, public nodata_i4
Reading of main model configurations.
subroutine check_l0domain(l0domain, ndomains)
subroutine, public set_land_cover_scenes_id(sim_per, lcyear_id)
Set land cover scenes IDs.
subroutine init_domain_variable(ndomains, optidata, domainmeta)
Initialization of the domain variables.
subroutine, public common_read_config(file_namelist)
Read main configurations commonly used by mHM, mRM and MPR.
Provides common types needed by mHM, mRM and/or mpr.
Provides structures needed by mHM, mRM and/or mpr.
integer(i4), parameter, public nprocesses
character(256), dimension(:), allocatable, public filelatlon
character(256), dimension(:), allocatable, public mhmfilerestartout
character(1024), public history
details on version/creation date
character(1024), public setup_description
any specific description of simulation
real(dp), dimension(:), allocatable, public resolutionhydrology
logical, public write_restart
character(256), dimension(:), allocatable, public lcfilename
type(domain_meta), public domainmeta
character(1024), public project_details
project including funding instituion., PI, etc.
character(1024), public contact
contact details, incl.
character(256), public dirconfigout
character(256), public conventions
convention used for dataset
character(1024), public simulation_type
e.g.
integer(i4), public nuniquel0domains
integer(i4), public nlcoverscene
character(256), dimension(:), allocatable, public dirlcover
integer(i4), dimension(:), allocatable, public lc_year_end
character(256), dimension(:), allocatable, public dirout
character(256), public dircommonfiles
character(256), dimension(:), allocatable, public dirmorpho
integer(i4), public iflag_cordinate_sys
integer(i4), dimension(nprocesses, 3), public processmatrix
character(256), dimension(:), allocatable, public mrmfilerestartout
character(1024), public mhm_details
developing institution, specific mHM revision
integer(i4), dimension(:), allocatable, public lc_year_start
Module containing all namelists representations.
type(nml_lcover_t), public nml_lcover
'LCover' namelist content
type(nml_processselection_t), public nml_processselection
'processSelection' namelist content
type(nml_mainconfig_t), public nml_mainconfig
'mainconfig' namelist content
type(nml_project_description_t), public nml_project_description
'project_description' namelist content
type(nml_directories_general_t), public nml_directories_general
'directories_general' namelist content
DOMAIN general description.