5.13.2-dev0
mHM
The mesoscale Hydrological Model
Loading...
Searching...
No Matches
mo_read_optional_data.f90
Go to the documentation of this file.
1!> \file mo_read_optional_data.f90
2!> \brief \copybrief mo_read_optional_data
3!> \details \copydetails mo_read_optional_data
4
5!> \brief Read optional data for mHM calibration.
6!> \details Data have to be provided in resolution of the hydrology.
7!> \authors Matthias Zink
8!> \date Mar 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_mhm
13
14 USE mo_kind, ONLY : i4, dp
15
16 IMPLICIT NONE
17
18 PRIVATE
19
20 PUBLIC :: readoptidataobs
21
22 ! ------------------------------------------------------------------
23
24CONTAINS
25
26 ! ------------------------------------------------------------------
27
28 ! NAME
29 ! readOptidataObs
30
31 ! PURPOSE
32 !> \brief Read evapotranspiration data from NetCDF file for calibration
33
34 !> \details This routine reads oberved evapotranspiration fields which are used for model
35 !> calibration. The evapotranspiration file is expected to be called "et.nc" with
36 !> a variable "et" inside. The data are read only for the evaluation period
37 !> they are intended to be used for calibration. Evapotranspiration data are only
38 !> read if one of the corresponding objective functions is chosen.
39
40 ! INTENT(IN)
41 !> \param[in] "integer(i4) :: iDomain" domain Id
42
43 ! HISTORY
44 !> \authors Johannes Brenner
45
46 !> \date Feb 2017
47
48 ! Modifications:
49 ! Robert Schweppe Jun 2018 - refactoring and reformatting
50 ! Maren Kaluza Oct 2019 - copied from evapotranspiration and adopted for tws
51 ! Maren Kaluza Nov 2019 - removed redundant code, reading works for any gridded optidata
52
53 subroutine readoptidataobs(iDomain, domainID, L1_optiObs)
54
56 use mo_common_variables, only : level1
57 use mo_optimization_types, only: optidata
58 use mo_message, only : message
59 use mo_read_nc, only : read_nc
60 use mo_string_utils, only : num2str
61 use mo_timer, only : timer_get, timer_start, &
62 timer_stop
63
64 implicit none
65
66 ! domain Id
67 integer(i4), intent(in) :: idomain
68
69 ! domain Id
70 integer(i4), intent(in) :: domainid
71
72 ! opti data
73 type(optidata), intent(inout) :: l1_optiobs
74
75 ! loop vars packing L1_data to L1_data_packed
76 integer(i4) :: t
77
78 ! level 1 number of culomns and rows
79 integer(i4) :: nrows1, ncols1
80
81 ! mask of level 1 for packing
82 logical, dimension(:, :), allocatable :: mask1
83
84 ! ncells1 of level 1
85 integer(i4) :: ncells1
86
87 ! data at level-1
88 real(dp), dimension(:, :, :), allocatable :: l1_data
89
90 ! mask at level-1
91 logical, dimension(:, :, :), allocatable :: l1_mask
92
93 integer(i4) :: ntimesteps_l1_opti ! [-] number of time steps in L1_opti_mask
94
95 ! get basic domain information at level-1
96 nrows1 = level1(idomain)%nrows
97 ncols1 = level1(idomain)%ncols
98 ncells1 = level1(idomain)%ncells
99 mask1 = level1(idomain)%mask
100
101 ! domain characteristics and read meteo header
102 call message(' Reading', trim(l1_optiobs%varname) ,'for domain: ', trim(adjustl(num2str(domainid))), ' ...')
103 call timer_start(1)
104 call read_nc(l1_optiobs%dir, nrows1, ncols1, trim(l1_optiobs%varname), mask1, l1_data, &
105 target_period = evalper(idomain), nctimestep = l1_optiobs%timeStepInput, nocheck = .true., maskout = l1_mask)
106
107 ! pack variables
108 ntimesteps_l1_opti = size(l1_data, 3)
109 allocate(l1_optiobs%dataObs(ncells1, ntimesteps_l1_opti))
110 allocate(l1_optiobs%maskObs(ncells1, ntimesteps_l1_opti))
111 do t = 1, ntimesteps_l1_opti
112 l1_optiobs%dataObs(:, t) = pack(l1_data(:, :, t), mask = mask1(:, :))
113 l1_optiobs%maskObs(:, t) = pack(l1_mask(:, :, t), mask = mask1(:, :))
114 end do
115
116 !free space
117 deallocate(l1_data)
118
119 call timer_stop(1)
120 call message(' in ', trim(num2str(timer_get(1), '(F9.3)')), ' seconds.')
121
122 end subroutine readoptidataobs
123
124END MODULE mo_read_optional_data
Provides structures needed by mHM, mRM and/or mpr.
type(period), dimension(:), allocatable, public evalper
Provides structures needed by mHM, mRM and/or mpr.
type(grid), dimension(:), allocatable, target, public level1
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.
Read optional data for mHM calibration.
subroutine, public readoptidataobs(idomain, domainid, l1_optiobs)
Read evapotranspiration data from NetCDF file for calibration.