5.13.3-dev0
mHM
The mesoscale Hydrological Model
Loading...
Searching...
No Matches
mo_prepare_gridded_lai.f90
Go to the documentation of this file.
1!> \file mo_prepare_gridded_lai.f90
2!> \brief \copybrief mo_prepare_gridded_lai
3!> \details \copydetails mo_prepare_gridded_lai
4
5!> \brief Prepare daily LAI fields (e.g., MODIS data) for mHM
6!> \details Prepare daily LAI fields(e.g., MODIS data) for mHM
7!> \authors John Craven & Rohini Kumar
8!> \date Aug 2013
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_mpr
13
14 ! This module provides routines to read daily gridded LAI data.
15
16 ! Written John Craven & Rohini Kumar, August 2013
17 ! Modified from mo_meteo_forcings
18
19 USE mo_kind, ONLY : i4, dp
20 use mo_message, only: error_message
21
22 IMPLICIT NONE
23
24 PRIVATE
25
28
29 ! ------------------------------------------------------------------
30
31CONTAINS
32
33 ! ------------------------------------------------------------------
34
35 ! NAME
36 ! prepare_gridded_daily_LAI_data
37
38 ! PURPOSE
39 !> \brief Prepare gridded daily LAI data
40
41 !> \details Prepare gridded daily LAI data at Level-0 (e.g., using MODIS datasets)
42
43 ! INTENT(IN)
44 !> \param[in] "integer(i4) :: iDomain, nrows, ncols" domain Id
45 !> \param[in] "integer(i4) :: iDomain, nrows, ncols" domain Id
46 !> \param[in] "integer(i4) :: iDomain, nrows, ncols" domain Id
47 !> \param[in] "logical, dimension(:, :) :: mask"
48
49 ! INTENT(IN), OPTIONAL
50 !> \param[in] "type(period), optional :: LAIPer_iDomain"
51
52 ! HISTORY
53 !> \authors John Craven & Rohini Kumar
54
55 !> \date Aug 2013
56
57 ! Modifications:
58 ! Matthias Cuntz & Juliane Mai Nov 2014 - use meteo reading routines
59 ! Robert Schweppe Jun 2018 - refactoring and reformatting
60
61 subroutine prepare_gridded_daily_lai_data(iDomain, nrows, ncols, mask, LAIPer_iDomain)
62
63 use mo_append, only : append
64 use mo_common_types, only: period
67 use mo_read_nc, only : read_nc
68
69 implicit none
70
71 ! domain Id
72 integer(i4), intent(in) :: idomain, nrows, ncols
73
74 logical, dimension(:, :), intent(in) :: mask
75
76 type(period), intent(in), optional :: laiper_idomain
77
78 integer(i4) :: ncells, ilai
79
80 ! data at level-0 [nRow X nCols X nTimeSteps]
81 real(dp), dimension(:, :, :), allocatable :: lai0_3d
82
83 ! data at level-0 [nCells X nTimeSteps]
84 real(dp), dimension(:, :), allocatable :: lai0_2d
85
86
87 ! select case depending on input data format
88 SELECT CASE(trim(inputformat_gridded_lai))
89
90 ! netcdf file input option
91 CASE('nc')
92 CALL read_nc(dirgridded_lai(idomain), nrows, ncols, &
93 'lai', mask, lai0_3d, target_period = laiper_idomain, &
94 lower = 1.00e-10_dp, upper = 30.0_dp, nctimestep = timestep_lai_input)
95 CASE DEFAULT
96 call error_message('***ERROR: No recognized input format')
97
98 END SELECT
99
100 ! pack variables
101 ncells = count(mask)
102 ! only set if not yet allocated (e.g. domain 1)
103 if (.not. allocated(laiboundaries)) then
104 nlai = size(lai0_3d, 3)
105 allocate(laiboundaries(nlai+1))
106 laiboundaries = [(ilai, ilai=1, nlai+1)]
107 end if
108 allocate(lai0_2d(ncells, nlai))
109
110 do ilai = 1, nlai
111 lai0_2d(:, ilai) = pack(lai0_3d(:, :, ilai), mask = mask(:, :))
112 end do
113
114 ! append to Global variable
115 call append(l0_gridded_lai, lai0_2d(:, :))
116
117 !free space
118 deallocate(lai0_2d, lai0_3d)
119
120 end subroutine prepare_gridded_daily_lai_data
121
122 ! ------------------------------------------------------------------
123
124 ! NAME
125 ! prepare_gridded_mean_monthly_LAI_data
126
127 ! PURPOSE
128 !> \brief prepare_gridded_mean_monthly_LAI_data
129
130 !> \details Long term mean monthly gridded LAI data at Level-0 (e.g., using MODIS datasets)
131 !> The netcdf file should contain 12 (calender months) gridded fields of climatological
132 !> LAI data at the input L0 data resolution.
133
134 ! INTENT(IN)
135 !> \param[in] "integer(i4) :: iDomain, nrows, ncols" domain Id
136 !> \param[in] "integer(i4) :: iDomain, nrows, ncols" domain Id
137 !> \param[in] "integer(i4) :: iDomain, nrows, ncols" domain Id
138 !> \param[in] "logical, dimension(:, :) :: mask"
139
140 ! HISTORY
141 !> \authors Rohini Kumar
142
143 !> \date Dec 2016
144
145 ! Modifications:
146 ! Robert Schweppe Jun 2018 - refactoring and reformatting
147
148 subroutine prepare_gridded_mean_monthly_lai_data(iDomain, nrows, ncols, mask)
149
150 use mo_append, only : append
152 use mo_ncread, only : get_ncdim, get_ncvar, get_ncvaratt
153 use mo_string_utils, only : num2str
154 use mo_utils, only : eq
155
156 implicit none
157
158 ! domain Id
159 integer(i4), intent(in) :: idomain, nrows, ncols
160
161 logical, dimension(:, :), intent(in) :: mask
162
163 integer(i4) :: ncells, ilai
164
165 ! data at level-0 [nRow X nCols X nTimeSteps]
166 real(dp), dimension(:, :, :), allocatable :: lai0_3d
167
168 ! data at level-0 [nCells X nTimeSteps]
169 real(dp), dimension(:, :), allocatable :: lai0_2d
170
171 integer(i4) :: t
172
173 ! name of NetCDF file
174 character(256) :: fname
175
176 ! netcdf attribute values
177 character(256) :: attvalues
178
179 ! datatype of attribute
180 integer(i4) :: datatype
181
182 ! dimension for NetCDF file
183 integer(i4), dimension(5) :: dimen
184
185 ! data nodata value
186 real(dp) :: nodata_value
187
188
189 fname = trim(dirgridded_lai(idomain)) // trim('lai.nc')
190
191 ! get dimensions
192 dimen = get_ncdim(trim(fname), 'lai')
193 if ((dimen(1) .ne. nrows) .or. (dimen(2) .ne. ncols)) then
194 call error_message('***ERROR: read_nc: mHM generated x and y are not matching NetCDF dimensions')
195 end if
196 if (dimen(3) .ne. 12) then
197 call error_message('***ERROR: read_nc: the time dimenion of LAI NetCDF file under the option-1 is not 12')
198 end if
199
200 ! determine no data value
201 call get_ncvaratt(trim(fname), 'lai', '_FillValue', attvalues, dtype = datatype)
202 ! convert to number
203 read(attvalues, *) nodata_value
204
205 call get_ncvar(trim(fname), 'lai', lai0_3d)
206
207 ! start checking values
208 do t = 1, dimen(3)
209 ! checking for nodata values if optional nocheck is given
210 if (any(eq(lai0_3d(:, :, t), nodata_value) .and. (mask))) then
211 call error_message('***ERROR: read_nc: nodata value within domain ', raise=.false.)
212 call error_message(' boundary in variable: ', 'lai', raise=.false.)
213 call error_message(' at timestep : ', trim(num2str(t)))
214 end if
215 ! optional check
216 if (any((lai0_3d(:, :, t) .lt. 0.0_dp) .AND. mask(:, :))) then
217 call error_message('***ERROR: read_nc: values in variable lai are lower than ', trim(num2str(0, '(F7.2)')), raise=.false.)
218 call error_message(' at timestep : ', trim(num2str(t)), raise=.false.)
219 call error_message('File: ', trim(fname), raise=.false.)
220 call error_message('Minval at timestep: ', trim(num2str(minval(lai0_3d(:, :, t)), '(F7.2)')), raise=.false.)
221 call error_message('Total minval: ', trim(num2str(minval(lai0_3d(:, :, :)), '(F7.2)')))
222 end if
223
224 if (any((lai0_3d(:, :, t) .gt. 30.0_dp) .AND. mask(:, :))) then
225 call error_message('***ERROR: read_nc: values in variable lai are greater than ', trim(num2str(30, '(F7.2)')), raise=.false.)
226 call error_message(' at timestep : ', trim(num2str(t)), raise=.false.)
227 call error_message('File: ', trim(fname), raise=.false.)
228 call error_message('Maxval at timestep: ', trim(num2str(maxval(lai0_3d(:, :, t)), '(F7.2)')), raise=.false.)
229 call error_message('Total maxval: ', trim(num2str(maxval(lai0_3d(:, :, :)), '(F7.2)')))
230 end if
231 end do
232
233 ! pack variables
234 ncells = count(mask)
235 ! only set if not yet allocated (e.g. domain 1)
236 if (.not. allocated(laiboundaries)) then
237 nlai = size(lai0_3d, 3)
238 allocate(laiboundaries(nlai+1))
239 laiboundaries = [(ilai, ilai=1, nlai+1)]
240 end if
241 allocate(lai0_2d(ncells, nlai))
242 do ilai = 1, nlai
243 lai0_2d(:, ilai) = pack(lai0_3d(:, :, ilai), mask = mask(:, :))
244 end do
245
246 ! append to Global variable
247 call append(l0_gridded_lai, lai0_2d(:, :))
248
249 !free space
250 deallocate(lai0_2d, lai0_3d)
251
253
254
255END MODULE mo_prepare_gridded_lai
Provides common types needed by mHM, mRM and/or mpr.
Global variables for mpr only.
character(256), dimension(:), allocatable, public dirgridded_lai
real(dp), dimension(:, :), allocatable, public l0_gridded_lai
real(dp), dimension(:), allocatable, public laiboundaries
character(256), public inputformat_gridded_lai
integer(i4), public timestep_lai_input
Prepare daily LAI fields (e.g., MODIS data) for mHM.
subroutine, public prepare_gridded_mean_monthly_lai_data(idomain, nrows, ncols, mask)
prepare_gridded_mean_monthly_LAI_data
subroutine, public prepare_gridded_daily_lai_data(idomain, nrows, ncols, mask, laiper_idomain)
Prepare gridded daily LAI data.
Reads forcing input data.
subroutine, public read_nc(folder, nrows, ncols, varname, mask, data, target_period, lower, upper, nctimestep, filename, nocheck, maskout, is_meteo, bound_error, ntstepforcingday)
Reads forcing input in NetCDF file format.