5.13.2-dev0
mHM
The mesoscale Hydrological Model
Loading...
Searching...
No Matches
mo_common_read_data.f90
Go to the documentation of this file.
1!> \file mo_common_read_data.f90
2!> \brief \copybrief mo_common_read_data
3!> \details \copydetails mo_common_read_data
4
5!> \brief Common reading routines
6!> \details Routines to read the DEM and landcover files.
7!> \authors Robert Schweppe
8!> \date Jun 2018
9!> \copyright Copyright 2005-\today, the mHM Developers, Luis Samaniego, Sabine Attinger: All rights reserved.
10!! mHM is released under the LGPLv3+ license \license_note
11!> \ingroup f_common
13 USE mo_kind, ONLY : i4, dp
14 use mo_message, only: message, error_message
15
16 IMPLICIT NONE
17
18 PRIVATE
19
20 PUBLIC :: read_dem, read_lcover
21
22 ! ------------------------------------------------------------------
23
24CONTAINS
25
26 ! NAME
27 ! read_dem
28
29 ! PURPOSE
30 !> \brief TODO: add description
31
32 !> \details TODO: add description
33
34 ! HISTORY
35 !> \authors Robert Schweppe
36
37 !> \date Jun 2018
38
39 ! Modifications:
40
41 subroutine read_dem
42
43 use mo_append, only : append
45 use mo_common_file, only : file_dem, udem
46 use mo_common_types, only: grid
49 use mo_grid, only : set_domain_indices
51 use mo_string_utils, only : num2str
52
53 implicit none
54
55 ! loop variables
56 integer(i4) :: domainid, idomain
57
58 ! file name of file to read
59 character(256) :: fname
60
61 real(dp), dimension(:, :), allocatable :: data_dp_2d
62
63 type(grid), pointer :: level0_idomain
64
65
66 ! ************************************************
67 ! READ SPATIAL DATA FOR EACH DOMAIN
68 ! ************************************************
69 ! allocate necessary variables at Level0
70 allocate(level0(domainmeta%nDomains))
71
72 do idomain = 1, domainmeta%nDomains
73 domainid = domainmeta%indices(idomain)
74
75 level0_idomain => level0(domainmeta%L0DataFrom(idomain))
76
77 ! check whether L0 data is shared
78 if (idomain .gt. 1) then
79 if (domainmeta%L0DataFrom(idomain) < idomain) then
80 !
81 call message(' Using dem of domain ', &
82 trim(adjustl(num2str(domainmeta%indices(domainmeta%L0DataFrom(idomain))))), ' for domain: ',&
83 trim(adjustl(num2str(idomain))), '...')
84
85 ! DO NOT read L0 data
86 cycle
87
88 end if
89 end if
90
91 call message(' Reading dem for domain: ', trim(adjustl(num2str(domainid))), ' ...')
92
93 ! Header (to check consistency)
94 fname = trim(adjustl(dirmorpho(idomain))) // trim(adjustl(file_dem))
95 call read_header_ascii(trim(fname), udem, &
96 level0_idomain%nrows, level0_idomain%ncols, level0_idomain%xllcorner, &
97 level0_idomain%yllcorner, level0_idomain%cellsize, level0_idomain%nodata_value)
98
99 ! check for L0 and L1 scale consistency
100 if(resolutionhydrology(idomain) .LT. level0_idomain%cellsize) then
101 call error_message('***ERROR: resolutionHydrology (L1) should be smaller than the input data resolution (L0)', raise=.false.)
102 call error_message(' check set-up (in mhm.nml) for domain: ', trim(adjustl(num2str(domainid))), ' ...')
103 end if
104
105 ! DEM + overall mask creation
106 fname = trim(adjustl(dirmorpho(idomain))) // trim(adjustl(file_dem))
107 call read_spatial_data_ascii(trim(fname), udem, &
108 level0_idomain%nrows, level0_idomain%ncols, level0_idomain%xllcorner, &
109 level0_idomain%yllcorner, level0_idomain%cellsize, data_dp_2d, level0_idomain%mask)
110
111 ! put global nodata value into array (probably not all grid cells have values)
112 data_dp_2d = merge(data_dp_2d, nodata_dp, level0_idomain%mask)
113 ! put data in variable
114 call append(l0_elev, pack(data_dp_2d, level0_idomain%mask))
115 ! deallocate arrays
116 deallocate(data_dp_2d)
117
118 level0_idomain%nCells = count(level0_idomain%mask)
119
120 end do
121
122 call set_domain_indices(level0, indices=domainmeta%L0DataFrom)
123
124 end subroutine read_dem
125
126 ! NAME
127 ! read_lcover
128
129 ! PURPOSE
130 !> \brief TODO: add description
131
132 !> \details TODO: add description
133
134 ! HISTORY
135 !> \authors Robert Schweppe
136
137 !> \date Jun 2018
138
139 ! Modifications:
140
141 subroutine read_lcover
142
143 use mo_append, only : append, paste
145 use mo_common_file, only : ulcoverclass
146 use mo_common_types, only: grid
149 use mo_string_utils, only : num2str
150
151 implicit none
152
153 ! loop variables
154 integer(i4) :: domainid, idomain, ivar
155
156 ! file name of file to read
157 character(256) :: fname
158
159 integer(i4), dimension(:, :), allocatable :: data_i4_2d
160
161 integer(i4), dimension(:, :), allocatable :: datamatrix_i4
162
163 logical, dimension(:, :), allocatable :: mask_2d
164
165 type(grid), pointer :: level0_idomain
166
167
168 do idomain = 1, domainmeta%nDomains
169 domainid = domainmeta%indices(idomain)
170
171 level0_idomain => level0(domainmeta%L0DataFrom(idomain))
172
173 ! check whether L0 data is shared
174 ! ToDo: check change
175 if (domainmeta%L0DataFrom(idomain) < idomain) then
176 call message(' Using lcover of domain ', &
177 trim(adjustl(num2str(domainmeta%L0DataFrom(idomain)))), ' for domain: ',&
178 trim(adjustl(num2str(domainid))), '...')
179 ! DO NOT read L0 data
180 cycle
181
182 end if
183
184 call message(' Reading lcover for domain: ', trim(adjustl(num2str(domainid))), ' ...')
185
186 ! LCover read in is realized seperated because of unknown number of scenes
187 do ivar = 1, nlcoverscene
188 fname = trim(adjustl(dirlcover(idomain))) // trim(adjustl(lcfilename(ivar)))
189 call read_spatial_data_ascii(trim(fname), ulcoverclass, &
190 level0_idomain%nrows, level0_idomain%ncols, level0_idomain%xllcorner, &
191 level0_idomain%yllcorner, level0_idomain%cellsize, data_i4_2d, mask_2d)
192 ! put global nodata value into array (probably not all grid cells have values)
193 data_i4_2d = merge(data_i4_2d, nodata_i4, mask_2d)
194 call paste(datamatrix_i4, pack(data_i4_2d, level0_idomain%mask), nodata_i4)
195 deallocate(data_i4_2d)
196 end do
197 call append(l0_lcover, datamatrix_i4)
198 deallocate(datamatrix_i4)
199
200 end do
201
202 end subroutine read_lcover
203
204end module mo_common_read_data
Reads spatial data files of ASCII format.
Provides constants commonly used by mHM, mRM and MPR.
real(dp), parameter, public nodata_dp
integer(i4), parameter, public nodata_i4
Provides file names and units for mRM.
integer, parameter ulcoverclass
Unit for LCover input data file.
character(len= *), parameter file_dem
DEM input data file.
integer, parameter udem
Unit for DEM input data file.
Common reading routines.
subroutine, public read_lcover
TODO: add description.
subroutine, public read_dem
TODO: add description.
Provides common types needed by mHM, mRM and/or mpr.
Provides structures needed by mHM, mRM and/or mpr.
real(dp), dimension(:), allocatable, public resolutionhydrology
character(256), dimension(:), allocatable, public lcfilename
type(domain_meta), public domainmeta
integer(i4), public nlcoverscene
character(256), dimension(:), allocatable, public dirlcover
real(dp), dimension(:), allocatable, public l0_elev
character(256), dimension(:), allocatable, public dirmorpho
integer(i4), dimension(:, :), allocatable, public l0_lcover
type(grid), dimension(:), allocatable, target, public level0
gridding tools
Definition mo_grid.f90:12
subroutine, public set_domain_indices(grids, indices)
TODO: add description.
Definition mo_grid.f90:205
Reads spatial input data.
subroutine, public read_header_ascii(filename, fileunit, header_ncols, header_nrows, header_xllcorner, header_yllcorner, header_cellsize, header_nodata)
Reads header lines of ASCII files.