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