Line data Source code
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
12 : MODULE mo_read_optional_data
13 :
14 : USE mo_kind, ONLY : i4, dp
15 :
16 : IMPLICIT NONE
17 :
18 : PRIVATE
19 :
20 : PUBLIC :: readOptidataObs
21 :
22 : ! ------------------------------------------------------------------
23 :
24 : CONTAINS
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 9 : subroutine readOptidataObs(iDomain, domainID, L1_optiObs)
54 :
55 : use mo_common_mhm_mrm_variables, only : evalPer
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 9 : logical, dimension(:, :), allocatable :: mask1
83 :
84 : ! ncells1 of level 1
85 : integer(i4) :: ncells1
86 :
87 : ! data at level-1
88 9 : real(dp), dimension(:, :, :), allocatable :: L1_data
89 :
90 : ! mask at level-1
91 9 : 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 9 : nrows1 = level1(iDomain)%nrows
97 9 : ncols1 = level1(iDomain)%ncols
98 9 : ncells1 = level1(iDomain)%ncells
99 765 : mask1 = level1(iDomain)%mask
100 :
101 : ! domain characteristics and read meteo header
102 9 : call message(' Reading', trim(L1_optiObs%varname) ,'for domain: ', trim(adjustl(num2str(domainID))), ' ...')
103 9 : call timer_start(1)
104 : call read_nc(L1_optiObs%dir, nRows1, nCols1, trim(L1_optiObs%varname), mask1, L1_data, &
105 9 : target_period = evalPer(iDomain), nctimestep = L1_optiObs%timeStepInput, nocheck = .TRUE., maskout = L1_mask)
106 :
107 : ! pack variables
108 9 : nTimeSteps_L1_opti = size(L1_data, 3)
109 36 : allocate(L1_optiObs%dataObs(nCells1, nTimeSteps_L1_opti))
110 36 : allocate(L1_optiObs%maskObs(nCells1, nTimeSteps_L1_opti))
111 698 : do t = 1, nTimeSteps_L1_opti
112 51490 : L1_optiObs%dataObs(:, t) = pack(L1_data(:, :, t), MASK = mask1(:, :))
113 51499 : L1_optiObs%maskObs(:, t) = pack(L1_mask(:, :, t), MASK = mask1(:, :))
114 : end do
115 :
116 : !free space
117 9 : deallocate(L1_data)
118 :
119 9 : call timer_stop(1)
120 9 : call message(' in ', trim(num2str(timer_get(1), '(F9.3)')), ' seconds.')
121 :
122 9 : end subroutine readOptidataObs
123 :
124 : END MODULE mo_read_optional_data
|