5.13.2-dev0
mHM
The mesoscale Hydrological Model
Loading...
Searching...
No Matches
mo_mrm_riv_temp_class.f90
Go to the documentation of this file.
1!> \file mo_mrm_riv_temp_class.f90
2!> \brief \copybrief mo_mrm_riv_temp_class
3!> \details \copydetails mo_mrm_riv_temp_class
4
5!> \brief Class for the river temperature calculations
6!> \details River temperature routing on top of mRM.
7!> \warning This feature is still experimental! River freezing is still missing.
8!> \version 0.1
9!> \authors Sebastian Mueller
10!> \date Sep 2020
11!> \copyright Copyright 2005-\today, the mHM Developers, Luis Samaniego, Sabine Attinger: All rights reserved.
12!! mHM is released under the LGPLv3+ license \license_note
13!> \ingroup f_mrm
15
16 use mo_kind, only: dp, i4
17
18 implicit none
19
20 private
21
22 !> \class riv_temp_type
23 !> \brief This is a container to define the river temperature routing in the current time step
24 !> \details This class provides all procedures to rout river tmperature through the river network.
25 !> \warning This feature is still experimental!
26 !! This first version doesn't provide ice covering, which means, that the river temperature
27 !! can drop below 0 [deg C] in winter.
28 type, public :: riv_temp_type
29 ! config settings
30 logical :: active = .false. !< state if this process is active
31 integer(i4) :: case = 0_i4 !< the selected process-case option
32 character(256) :: nml_name = 'config_riv_temp' !< namelist name in mhm namelist
33 ! riv geometry
34 character(256), dimension(:), allocatable :: dir_riv_widths !< Directory where river widths are stored
35 character(256) :: riv_widths_file !< file name for river widths
36 character(256) :: riv_widths_name !< variable name for river widths
37 real(dp), public, dimension(:), allocatable :: l11_riv_widths !< river widths in L11
38 real(dp), public, dimension(:), allocatable :: l11_riv_areas !< river area in L11
39 ! PET related vars
40 real(dp) :: albedo_water !< albedo of open water
41 real(dp) :: pt_a_water !< priestley taylor alpha parameter for PET on open water
42 ! sensible heat flux related
43 real(dp) :: emissivity_water !< emissivity of water
44 real(dp) :: turb_heat_ex_coeff !< lateral heat exchange coefficient water <-> air
45 !> cutoff value for temperature
46 real(dp) :: delta_t = 0.1_dp
47 ! controlling variables for iterative solver
48 integer(i4) :: max_iter !< input: maximal number of iterations done
49 real(dp) :: delta_iter !< input: convergence criteria for iterative solver
50 real(dp) :: step_iter !< input: step-size for linear search
51 logical :: first_iter !< whether it is at the first iteration (to determine search direction)
52 logical :: up_iter !< whether the search direction is upwards
53 logical :: bisect_iter !< whether to do the bisection search part (after the interval is found)
54 real(dp) :: up_bnd_iter !< upper bound for the current bisection step
55 real(dp) :: low_bnd_iter !< lower bound for the current bisection step
56 ! accumulated later fluxes (in current time-step)
57 real(dp), dimension(:), allocatable :: l1_runoff_e !< runoff energy at L1 level
58 real(dp), dimension(:), allocatable :: l1_acc_ssrd !< accumulated shortwave radiation at L1 level
59 real(dp), dimension(:), allocatable :: l1_acc_strd !< accumulated longwave radiation at L1 level
60 real(dp), dimension(:), allocatable :: l1_acc_temp !< accumulated air temperature at L1 level
61 real(dp), dimension(:), allocatable :: l1_ssrd_calc !< current shortwave radiation at L1 level
62 real(dp), dimension(:), allocatable :: l1_strd_calc !< current longwave radiation at L1 level
63 real(dp), dimension(:), allocatable :: l1_tann_calc !< current mean air temperature at L1 level
64 integer(i4) :: ts_cnt !< sub time-step counter for accumulation of meteo
65 ! vars for routing
66 integer(i4) :: s11 !< starting index for current L11 domain
67 integer(i4) :: e11 !< ending index for current L11 domain
68 real(dp), dimension(:,:), allocatable :: netnode_e_in !< Total energy inputs at t-1 and t
69 real(dp), dimension(:,:), allocatable :: netnode_e_r !< energy leaving at t-1 and t
70 real(dp), dimension(:), allocatable :: netnode_e_mod !< Simulated routed energy
71 real(dp), dimension(:), allocatable :: netnode_e_out !< total energy source from cell in L11
72 real(dp), dimension(:), allocatable :: l11_srad_net !< net short wave radiation at L11
73 real(dp), dimension(:), allocatable :: l11_lrad_in !< incoming long wave radiation at L11
74 real(dp), dimension(:), allocatable :: l11_air_temp !< air temp at L11
75 ! variable containing river temp for each domain at L11 level
76 real(dp), dimension(:), allocatable :: river_temp !< resulting river temp at L11 in [deg C]
77 contains
78 ! config and inits
79 !> \copydoc mo_mrm_riv_temp_class::config
80 procedure :: config !< \see mo_mrm_riv_temp_class::config
81 !> \copydoc mo_mrm_riv_temp_class::init
82 procedure :: init !< \see mo_mrm_riv_temp_class::init
83 !> \copydoc mo_mrm_riv_temp_class::init_area
84 procedure :: init_area !< \see mo_mrm_riv_temp_class::init_area
85 !> \copydoc mo_mrm_riv_temp_class::init_riv_temp
86 procedure :: init_riv_temp !< \see mo_mrm_riv_temp_class::init_riv_temp
87
88 ! source accumulations
89 !> \copydoc mo_mrm_riv_temp_class::acc_source_e
90 procedure :: acc_source_e !< \see mo_mrm_riv_temp_class::acc_source_e
91 !> \copydoc mo_mrm_riv_temp_class::finalize_source_e
92 procedure :: finalize_source_e !< \see mo_mrm_riv_temp_class::finalize_source_e
93
94 ! temp-energy routing routines
95 !> \copydoc mo_mrm_riv_temp_class::get_lrad_out
96 procedure :: get_lrad_out !< \see mo_mrm_riv_temp_class::get_lrad_out
97 !> \copydoc mo_mrm_riv_temp_class::get_lat_heat
98 procedure :: get_lat_heat !< \see mo_mrm_riv_temp_class::get_lat_heat
99 !> \copydoc mo_mrm_riv_temp_class::get_sens_heat
100 procedure :: get_sens_heat !< \see mo_mrm_riv_temp_class::get_sens_heat
101 !> \copydoc mo_mrm_riv_temp_class::get_e_io
102 procedure :: get_e_io !< \see mo_mrm_riv_temp_class::get_e_io
103 !> \copydoc mo_mrm_riv_temp_class::l11_routing_e
104 procedure :: l11_routing_e !< \see mo_mrm_riv_temp_class::l11_routing_e
105
106 ! helper for iterative solver
107 !> \copydoc mo_mrm_riv_temp_class::init_iter
108 procedure :: init_iter !< \see mo_mrm_riv_temp_class::init_iter
109 !> \copydoc mo_mrm_riv_temp_class::next_iter
110 procedure :: next_iter !< \see mo_mrm_riv_temp_class::next_iter
111
112 ! care taker
113 !> \copydoc mo_mrm_riv_temp_class::reset_timestep
114 procedure :: reset_timestep !< \see mo_mrm_riv_temp_class::reset_timestep
115 !> \copydoc mo_mrm_riv_temp_class::alloc_lateral
116 procedure :: alloc_lateral !< \see mo_mrm_riv_temp_class::alloc_lateral
117 !> \copydoc mo_mrm_riv_temp_class::dealloc_lateral
118 procedure :: dealloc_lateral !< \see mo_mrm_riv_temp_class::dealloc_lateral
119 !> \copydoc mo_mrm_riv_temp_class::clean_up
120 procedure :: clean_up !< \see mo_mrm_riv_temp_class::clean_up
121
122 end type riv_temp_type
123
124contains
125
126
127 !> \brief clean up
128 subroutine clean_up( &
129 self &
130 )
131 implicit none
132
133 class(riv_temp_type), intent(inout) :: self
134
135 if ( allocated(self%L1_runoff_E) ) deallocate(self%L1_runoff_E)
136 if ( allocated(self%L1_acc_strd) ) deallocate(self%L1_acc_strd)
137 if ( allocated(self%L1_acc_ssrd) ) deallocate(self%L1_acc_ssrd)
138 if ( allocated(self%L1_acc_temp) ) deallocate(self%L1_acc_temp)
139 if ( allocated(self%dir_riv_widths) ) deallocate(self%dir_riv_widths)
140 if ( allocated(self%L11_riv_widths) ) deallocate(self%L11_riv_widths)
141 if ( allocated(self%L11_riv_areas) ) deallocate(self%L11_riv_areas)
142 if ( allocated(self%netNode_E_IN) ) deallocate(self%netNode_E_IN)
143 if ( allocated(self%netNode_E_R) ) deallocate(self%netNode_E_R)
144 if ( allocated(self%netNode_E_mod) ) deallocate(self%netNode_E_mod)
145 if ( allocated(self%netNode_E_out) ) deallocate(self%netNode_E_out)
146 if ( allocated(self%L11_srad_net) ) deallocate(self%L11_srad_net)
147 if ( allocated(self%L11_lrad_in) ) deallocate(self%L11_lrad_in)
148 if ( allocated(self%L11_air_temp) ) deallocate(self%L11_air_temp)
149 if ( allocated(self%river_temp) ) deallocate(self%river_temp)
150 ! meteo arrays
151 if ( allocated(self%L1_ssrd_calc) ) deallocate(self%L1_ssrd_calc)
152 if ( allocated(self%L1_strd_calc) ) deallocate(self%L1_strd_calc)
153 if ( allocated(self%L1_tann_calc) ) deallocate(self%L1_tann_calc)
154
155 end subroutine clean_up
156
157 !> \brief configure the \ref riv_temp_type class from the mhm namelist
158 subroutine config( &
159 self, &
160 file_namelist, &
161 unamelist, &
162 file_namelist_param, &
163 unamelist_param &
164 )
165
168 use mo_nml, only : close_nml, open_nml, position_nml
169 use mo_check, only : check_dir
170 USE mo_string_utils, ONLY : num2str
171
172 implicit none
173
174 class(riv_temp_type), intent(inout) :: self
175 character(*), intent(in) :: file_namelist !< mhm namelist file
176 character(*), intent(in) :: file_namelist_param !< mhm parameter namelist file
177 integer, intent(in) :: unamelist
178 integer, intent(in) :: unamelist_param
179
180 ! parameter to read in
181 real(dp) :: albedo_water ! albedo of open water
182 real(dp) :: pt_a_water ! priestley taylor alpha parameter for PET on open water
183 real(dp) :: emissivity_water ! emissivity of water
184 real(dp) :: turb_heat_ex_coeff ! lateral heat exchange coefficient water <-> air
185 ! controlling variables for iterative solver
186 integer(i4) :: max_iter
187 real(dp) :: delta_iter
188 real(dp) :: step_iter
189 ! files for river widths
190 character(256), dimension(maxNoDomains) :: dir_riv_widths
191 character(256) :: riv_widths_file ! file name for river widths
192 character(256) :: riv_widths_name ! variable name for river widths
193
194 integer(i4) :: iDomain, domainID
195
196 ! namelist for river temperature configuration
197 namelist /config_riv_temp/ &
198 albedo_water, &
199 pt_a_water, &
200 emissivity_water, &
201 turb_heat_ex_coeff, &
202 max_iter, &
203 delta_iter, &
204 step_iter, &
205 riv_widths_file, &
206 riv_widths_name, &
207 dir_riv_widths
208
209 ! allocate the directory arrays
210 allocate(self%dir_riv_widths(domainmeta%nDomains))
211
212 ! default values
213 albedo_water = 0.15_dp
214 pt_a_water = 1.26_dp
215 emissivity_water = 0.96_dp
216 turb_heat_ex_coeff = 20.0_dp
217 max_iter = 20_i4
218 delta_iter = 1.0e-02_dp
219 step_iter = 5.0_dp
220
221 ! open the namelist file
222 call open_nml(file_namelist, unamelist, quiet=.true.)
223 ! find the river-temp config namelist
224 call position_nml(self%nml_name, unamelist)
225 ! read the river-temp config namelist
226 read(unamelist, nml=config_riv_temp)
227
228 self%albedo_water = albedo_water
229 self%pt_a_water = pt_a_water
230 self%emissivity_water = emissivity_water
231 self%turb_heat_ex_coeff = turb_heat_ex_coeff
232 self%delta_iter = delta_iter
233 self%max_iter = max_iter
234 self%step_iter = step_iter
235 self%riv_widths_file = riv_widths_file
236 self%riv_widths_name = riv_widths_name
237 self%dir_riv_widths = dir_riv_widths
238
239 ! closing the namelist file
240 call close_nml(unamelist)
241
242 do idomain = 1, domainmeta%nDomains
243 domainid = domainmeta%indices(idomain)
244 call check_dir( &
245 path=self%dir_riv_widths(idomain), &
246 text="(domain "//trim(num2str(domainid,'(I3)'))//") River widths directory:", &
247 raise=.true., &
248 tab=4, &
249 text_length=40 &
250 )
251 end do
252
253 end subroutine config
254
255 !> \brief initalize the \ref riv_temp_type class for the current domain
256 subroutine init( &
257 self, &
258 nCells &
259 )
260 use mo_append, only : append
262 ! use mo_common_variables, only : level0, domainMeta
263
264 implicit none
265
266 class(riv_temp_type), intent(inout) :: self
267 integer(i4), intent(in) :: nCells !< number of level-11 cells for the current domain
268
269 real(dp), dimension(:), allocatable :: dummy_Vector11
270 real(dp), dimension(:, :), allocatable :: dummy_Matrix11_IT
271
272 ! dummy vector and matrix
273 allocate(dummy_vector11(ncells))
274 allocate(dummy_matrix11_it(ncells, nroutingstates))
275 dummy_vector11(:) = 0.0_dp
276 dummy_matrix11_it(:, :) = 0.0_dp
277
278 ! simulated energy flux at each node
279 call append(self%netNode_E_mod, dummy_vector11)
280 ! simulated river temperature at each node
281 call append(self%river_temp, dummy_vector11)
282 ! Total outflow from cells L11 at time tt
283 call append(self%netNode_E_out, dummy_vector11)
284 ! Total discharge inputs at t-1 and t
285 call append(self%netNode_E_IN, dummy_matrix11_it)
286 ! Routed outflow leaving a node
287 call append(self%netNode_E_R, dummy_matrix11_it)
288
289 ! free space
290 deallocate(dummy_vector11)
291 deallocate(dummy_matrix11_it)
292
293 end subroutine init
294
295 !> \brief initialize the river area of \ref riv_temp_type class for the current domain
296 subroutine init_area( &
297 self, &
298 iDomain, &
299 L11_netPerm, &
300 L11_fromN, &
301 L11_length, &
302 nLinks, &
303 nCells, &
304 nrows, &
305 ncols, &
306 L11_mask &
307 )
308 use mo_append, only : append
309 use mo_read_nc, only: read_const_nc
310
311 implicit none
312
313 class(riv_temp_type), intent(inout) :: self
314 integer(i4), intent(in) :: iDomain !< Domain ID
315 integer(i4), dimension(:), intent(in) :: L11_netPerm !< L11 routing order
316 integer(i4), dimension(:), intent(in) :: L11_fromN !< L11 source grid cell order
317 real(dp), dimension(:), intent(in) :: L11_length !< L11 link length
318 integer(i4), intent(in) :: nLinks !< number of L11 links in the current domain
319 integer(i4), intent(in) :: nCells !< number of L11 cells of the current domain
320 integer(i4), intent(in) :: ncols !< Number of columns
321 integer(i4), intent(in) :: nrows !< Number of rows
322 logical, dimension(:, :), intent(in) :: L11_mask !< the mask for valid cells in the original grid (nrows, ncols)
323
324 real(dp), dimension(:,:), allocatable :: L11_data ! read data from file
325 real(dp), dimension(:), allocatable :: L11_riv_widths, L11_riv_areas
326
327 integer(i4) :: i, k, iNode
328
329 call read_const_nc(&
330 trim(self%dir_riv_widths(idomain)), &
331 nrows, &
332 ncols, &
333 self%riv_widths_name, &
334 l11_data, &
335 self%riv_widths_file &
336 )
337
338 allocate(l11_riv_widths(ncells))
339 l11_riv_widths(:) = pack(l11_data(:,:), mask=l11_mask)
340 call append(self%L11_riv_widths, l11_riv_widths)
341
342 allocate(l11_riv_areas(ncells))
343 ! at area at Outlet is 0
344 l11_riv_areas = 0.0_dp
345 do k = 1, nlinks
346 ! get LINK routing order -> i
347 i = l11_netperm(k)
348 inode = l11_fromn(i)
349 l11_riv_areas(inode) = l11_length(i) * l11_riv_widths(inode)
350 end do
351 call append(self%L11_riv_areas, l11_riv_areas)
352
353 deallocate(l11_data)
354 deallocate(l11_riv_widths)
355 deallocate(l11_riv_areas)
356
357 end subroutine init_area
358
359 !> \brief initialize the river temperature of \ref riv_temp_type class for the current domain
360 subroutine init_riv_temp( &
361 self, &
362 temp_air, &
363 efecarea, &
364 L1_L11_Id, &
365 L11_areacell, &
366 L11_L1_Id, &
367 map_flag &
368 )
369
371
372 implicit none
373
374 class(riv_temp_type), intent(inout) :: self
375 real(dp), dimension(:), intent(in) :: temp_air !< air temperature [degC] for current timestep
376 real(dp), intent(in), dimension(:) :: efecarea !< effective area in [km2] at Level 1
377 integer(i4), intent(in), dimension(:) :: L1_L11_Id !< L11 Ids mapped on L1
378 real(dp), intent(in), dimension(:) :: L11_areacell !< effective area in [km2] at Level 11
379 integer(i4), intent(in), dimension(:) :: L11_L1_Id !< L1 Ids mapped on L11
380 logical, intent(in) :: map_flag !< Flag indicating whether routing resolution is higher than hydrologic one
381
382 ! map temperature from L1 to L11
383 call l11_meteo_acc( &
384 temp_air, efecarea, l1_l11_id, l11_areacell, l11_l1_id, map_flag, self%river_temp(self%s11 : self%e11))
385
386 ! assure positive temperature
387 self%river_temp(self%s11 : self%e11) = max(self%delta_T, self%river_temp(self%s11 : self%e11))
388
389 end subroutine init_riv_temp
390
391 !> \brief reset \ref riv_temp_type class for next timestep
392 subroutine reset_timestep(self)
393 implicit none
394
395 class(riv_temp_type), intent(inout) :: self
396
397 ! set sub time-step counter to 0
398 self%ts_cnt = 0_i4
399
400 self%L1_runoff_E = 0.0_dp
401 self%L1_acc_strd = 0.0_dp
402 self%L1_acc_ssrd = 0.0_dp
403 self%L1_acc_temp = 0.0_dp
404
405 end subroutine reset_timestep
406
407 !> \brief allocate lateral temp components of \ref riv_temp_type class for current domain
408 subroutine alloc_lateral( &
409 self, &
410 nCells &
411 )
412 implicit none
413
414 class(riv_temp_type), intent(inout) :: self
415 integer(i4), intent(in) :: nCells !< number of level-1 cells for the current domain
416
417 allocate(self%L1_runoff_E(ncells))
418 allocate(self%L1_acc_strd(ncells))
419 allocate(self%L1_acc_ssrd(ncells))
420 allocate(self%L1_acc_temp(ncells))
421 ! meteo arrays
422 allocate(self%L1_ssrd_calc(ncells))
423 allocate(self%L1_strd_calc(ncells))
424 allocate(self%L1_tann_calc(ncells))
425 ! init these arrays to 0
426 call self%reset_timestep()
427
428 end subroutine alloc_lateral
429
430 !> \brief deallocate lateral temp components of \ref riv_temp_type
431 subroutine dealloc_lateral( &
432 self &
433 )
434 implicit none
435
436 class(riv_temp_type), intent(inout) :: self
437
438 deallocate(self%L1_runoff_E)
439 deallocate(self%L1_acc_strd)
440 deallocate(self%L1_acc_ssrd)
441 deallocate(self%L1_acc_temp)
442 deallocate(self%L1_ssrd_calc)
443 deallocate(self%L1_strd_calc)
444 deallocate(self%L1_tann_calc)
445
446 end subroutine dealloc_lateral
447
448 !> \brief accumulate energy sources of \ref riv_temp_type
449 subroutine acc_source_e( &
450 self, &
451 fSealed_area_fraction, &
452 fast_interflow, &
453 slow_interflow, &
454 baseflow, &
455 direct_runoff, &
456 temp_air &
457 )
458
460
461 implicit none
462
463 class(riv_temp_type), intent(inout) :: self
464 real(dp), dimension(:), intent(in) :: fSealed_area_fraction !< sealed area fraction [1]
465 real(dp), dimension(:), intent(in) :: fast_interflow !< \f$ q_0 \f$ Fast runoff component [mm TS-1]
466 real(dp), dimension(:), intent(in) :: slow_interflow !< \f$ q_1 \f$ Slow runoff component [mm TS-1]
467 real(dp), dimension(:), intent(in) :: baseflow !< \f$ q_2 \f$ Baseflow [mm TS-1]
468 real(dp), dimension(:), intent(in) :: direct_runoff !< \f$ q_D \f$ Direct runoff from impervious areas [mm TS-1]
469 real(dp), dimension(:), intent(in) :: temp_air !< air temperature [K]
470
471 ! increase the sub time-step counter
472 self%ts_cnt = self%ts_cnt + 1_i4
473
474 ! caclucate the temperature energy of the runoffs at L1 in [K mm]
475 ! automatically accumulate them
476 call calc_l1_runoff_e( &
477 fsealed_area_fraction, &
478 fast_interflow, slow_interflow, baseflow, direct_runoff, &
479 temp_air, self%L1_tann_calc, &
480 self%L1_runoff_E & ! will be added here
481 )
482 ! accumulate meteo forcings (will be averaged with sub time-step counter later)
483 self%L1_acc_ssrd = self%L1_acc_ssrd + self%L1_ssrd_calc
484 self%L1_acc_strd = self%L1_acc_strd + self%L1_strd_calc
485 self%L1_acc_temp = self%L1_acc_temp + temp_air
486
487 end subroutine acc_source_e
488
489 !> \brief finalize energy sources of \ref riv_temp_type
490 subroutine finalize_source_e( &
491 self, &
492 efecarea, &
493 L1_L11_Id, &
494 L11_areacell, &
495 L11_L1_Id, &
496 timestep, &
497 map_flag &
498 )
499
501
502 implicit none
503
504 class(riv_temp_type), intent(inout) :: self
505 !> effective area in [km2] at Level 1
506 real(dp), intent(in), dimension(:) :: efecarea
507 !> L11 Ids mapped on L1
508 integer(i4), intent(in), dimension(:) :: L1_L11_Id
509 !> effective area in [km2] at Level 11
510 real(dp), intent(in), dimension(:) :: L11_areacell
511 !> L1 Ids mapped on L11
512 integer(i4), intent(in), dimension(:) :: L11_L1_Id
513 !> simulation timestep in [h]
514 integer(i4), intent(in) :: timestep
515 !> Flag indicating whether routing resolution is higher than hydrologic one
516 logical, intent(in) :: map_flag
517
518 ! prepare temporal variables for temp-routing
519 if ( allocated(self%L11_srad_net) ) deallocate(self%L11_srad_net)
520 if ( allocated(self%L11_lrad_in) ) deallocate(self%L11_lrad_in)
521 if ( allocated(self%L11_air_temp) ) deallocate(self%L11_air_temp)
522 allocate(self%L11_srad_net(size(l11_l1_id, dim = 1)))
523 allocate(self%L11_lrad_in(size(l11_l1_id, dim = 1)))
524 allocate(self%L11_air_temp(size(l11_l1_id, dim = 1)))
525
526 ! convert the L1 runoff temperature energy from [K mm] to [K m3 s-1] on L11
527 call l11_runoff_acc( &
528 self%L1_runoff_E, &
529 efecarea, &
530 l1_l11_id, &
531 l11_areacell, &
532 l11_l1_id, &
533 timestep, &
534 map_flag, &
535 self%netNode_E_out(self%s11 : self%e11) &
536 )
537 ! all meteo forcings are accumulated over sub time-steps and need to be divided by the sub time-step counter
538 ! short wave radiation on L11
539 self%L1_acc_ssrd = self%L1_acc_ssrd / real(self%ts_cnt, dp)
540 call l11_meteo_acc(self%L1_acc_ssrd, efecarea, l1_l11_id, l11_areacell, l11_l1_id, map_flag, self%L11_srad_net)
541 ! calculate the net-radiation from incoming/outging short radiation
542 self%L11_srad_net = self%L11_srad_net * (1._dp - self%albedo_water)
543 ! long wave radiation on L11
544 self%L1_acc_strd = self%L1_acc_strd / real(self%ts_cnt, dp)
545 call l11_meteo_acc(self%L1_acc_strd, efecarea, l1_l11_id, l11_areacell, l11_l1_id, map_flag, self%L11_lrad_in)
546 ! air temperature on L11
547 self%L1_acc_temp = self%L1_acc_temp / real(self%ts_cnt, dp)
548 call l11_meteo_acc(self%L1_acc_temp, efecarea, l1_l11_id, l11_areacell, l11_l1_id, map_flag, self%L11_air_temp)
549
550 end subroutine finalize_source_e
551
552 !> \brief get outgoing longwave radiation of \ref riv_temp_type
553 !> \return outgoing longwave radiation
554 real(dp) function get_lrad_out(self, riv_temp) result(lrad_out)
555
556 use mo_constants, only: sigma_dp
557
558 implicit none
559
560 class(riv_temp_type), intent(in) :: self
561 !> river temperature in K
562 real(dp), intent(in) :: riv_temp
563
564 ! outgoing longwave radiation from Boltzmann equation
565 lrad_out = self%emissivity_water * sigma_dp * riv_temp ** 4_i4
566
567 end function get_lrad_out
568
569 !> \brief latent heat flux of \ref riv_temp_type
570 !> \return latent heat flux
571 real(dp) function get_lat_heat(self, air_temp, netrad) result(lat_heat)
572
573 use mo_constants, only : psychro_dp
574 use mo_pet, only : slope_satpressure
575
576 implicit none
577
578 class(riv_temp_type), intent(in) :: self
579 !> air temperature in deg C
580 real(dp), intent(in) :: air_temp
581 !> net radiation in W * m-2
582 real(dp), intent(in) :: netrad
583
584 ! save slope of saturation vapor pressure curve
585 real(dp) :: delta
586
587 delta = slope_satpressure(air_temp) ! slope of saturation vapor pressure curve
588 ! latent heat for priestley taylor PET formula
589 lat_heat = self%pt_a_water * delta / (psychro_dp + delta) * max(netrad, 0._dp)
590
591 end function get_lat_heat
592
593 !> \brief sensible heat flux of \ref riv_temp_type
594 !> \return sensible heat flux
595 real(dp) function get_sens_heat(self, air_temp, riv_temp) result(sens_heat)
596
597 implicit none
598
599 class(riv_temp_type), intent(in) :: self
600 real(dp), intent(in) :: air_temp !< air temperature in [deg C]
601 real(dp), intent(in) :: riv_temp !< river temperature in [deg C]
602
603 ! sensible heat flux resulting for temp. diff.: river <-> air
604 sens_heat = self%turb_heat_ex_coeff * (riv_temp - air_temp)
605
606 end function get_sens_heat
607
608 !> \brief get complete energy source of \ref riv_temp_type at given cell
609 !> \return energy IO
610 real(dp) function get_e_io(self, riv_temp, cell) result(e_io)
611
612 use mo_constants, only : t0_dp, cp_w_dp
613 use mo_mhm_constants, only : h2odens
614
615 implicit none
616
617 class(riv_temp_type), intent(in) :: self
618 !> given river temperature in K to calculate heat fluxes
619 real(dp), intent(in) :: riv_temp
620 !> cell index in the current domain
621 integer(i4), intent(in) :: cell
622 !> net radiation calc from short and longwave radiation in/out
623 real(dp) :: netrad, sens_heat, lat_heat
624
625 ! net radiation
626 netrad = self%L11_srad_net(cell) + self%L11_lrad_in(cell) - self%get_lrad_out(riv_temp)
627 ! sensible heat flux
628 sens_heat = self%get_sens_heat(self%L11_air_temp(cell), riv_temp - t0_dp)
629 ! latent heat flux
630 lat_heat = self%get_lat_heat(self%L11_air_temp(cell), netrad)
631
632 ! convert energy flux [W m-2] to [K m s-1]
633 ! needs to be multiplied with river-area to result in temp-energy flux [K m3 s-1]
634 e_io = (netrad - sens_heat - lat_heat) / (h2odens * cp_w_dp)
635
636 end function get_e_io
637
638 !> \brief execute the temperature routing of \ref riv_temp_type
639 subroutine l11_routing_e( &
640 self, &
641 nLinks, &
642 netPerm, &
643 netLink_fromN, &
644 netLink_toN, &
645 netLink_C1, &
646 netLink_C2, &
647 nInflowGauges, &
648 InflowHeadwater, &
649 InflowNodeList, &
650 L11_qTR, &
651 L11_Qmod &
652 )
653 use mo_constants, only : t0_dp
654
655 implicit none
656
657 class(riv_temp_type), intent(inout) :: self
658 !> number of stream segment (reaches)
659 integer(i4), intent(in) :: nLinks
660 !> routing order of a given domain (permutation)
661 integer(i4), dimension(:), intent(in) :: netPerm
662 !> from node
663 integer(i4), dimension(:), intent(in) :: netLink_fromN
664 !> to node
665 integer(i4), dimension(:), intent(in) :: netLink_toN
666 !> routing parameter C1 (\cite CMM1988 p. 25-41)
667 real(dp), dimension(:), intent(in) :: netLink_C1
668 !> routing parameters C2 (id)
669 real(dp), dimension(:), intent(in) :: netLink_C2
670 !> [-] number of inflow points
671 integer(i4), intent(in) :: nInflowGauges
672 !> [-] if to consider headwater cells of inflow gauge
673 logical, dimension(:), intent(in) :: InflowHeadwater
674 !> [-] L11 ID of inflow points
675 integer(i4), dimension(:), intent(in) :: InflowNodeList
676 !> [m3 s-1] Transformed outflow leaving node I at current timestep(Muskingum)
677 real(dp), intent(in), dimension(:) :: L11_qTR
678 !> [m3 s-1] Simulated routed discharge
679 real(dp), intent(in), dimension(:) :: L11_Qmod
680
681 integer(i4) :: i, k, m, iNode, tNode, L11in, L11to
682 real(dp) :: E_IO_in, E_IO, T_est, T_rout
683
684 ! initialize total input at point time (t+1) in all nodes
685 self%netNode_E_IN(self%s11 : self%e11, 2) = 0.0_dp
686
687 riv_loop: do k = 1, nlinks
688 ! get LINK routing order -> i
689 i = netperm(k)
690 inode = netlink_fromn(i)
691 tnode = netlink_ton(i)
692 ! indices for concaternated domains
693 l11in = inode + self%s11 - 1_i4
694 l11to = tnode + self%s11 - 1_i4
695
696 ! accumulate all inputs in iNode
697 self%netNode_E_IN(l11in, 2) = self%netNode_E_IN(l11in, 2) + self%netNode_E_out(l11in)
698 ! calculate the river temp at the upstream node
699 self%river_temp(l11in) = self%netNode_E_IN(l11in, 2) / l11_qmod(inode) - t0_dp
700
701 ! first guess for the routed temperature is the upstream temp
702 t_est = self%river_temp(l11in) + t0_dp
703 ! calculate the meteo energy source at the upstream node
704 e_io_in = self%get_E_IO(t_est, inode)
705
706 ! initialize iteration
707 call self%init_iter()
708 ! perform iteration
709 iterloop: do m=1, self%max_iter
710 ! meteo energy flux depending on the estimated routed temperature (avg. IN- and OUT-node)
711 e_io = (e_io_in + self%get_E_IO(t_est, tnode)) * 0.5_dp * self%L11_riv_areas(l11in)
712 ! routing iNode
713 self%netNode_E_R(l11in, 2) = self%netNode_E_R(l11in, 1) &
714 + netlink_c1(i) * (self%netNode_E_IN(l11in, 1) - self%netNode_E_R(l11in, 1)) &
715 + netlink_c2(i) * (self%netNode_E_IN(l11in, 2) + e_io - self%netNode_E_IN(l11in, 1))
716 ! calculate the routed temperature
717 t_rout = self%netNode_E_R(l11in, 2) / l11_qtr(inode)
718 ! check for convergence
719 if ( abs(t_rout - t_est) < self%delta_iter ) exit iterloop
720 ! get new estimate for next iteration
721 call self%next_iter(t_est, t_rout)
722 end do iterloop
723
724 ! add the meteo energy flux to the accumulated incoming energy at the IN-node
725 self%netNode_E_IN(l11in, 2) = self%netNode_E_IN(l11in, 2) + e_io
726
727 ! check if the inflow from upstream cells should be deactivated
728 if (ninflowgauges .GT. 0) then
729 do i = 1, ninflowgauges
730 ! check if downstream Node (tNode) is inflow gauge and headwaters should be ignored
731 if ((tnode == inflownodelist(i)) .AND. (.NOT. inflowheadwater(i))) &
732 self%netNode_E_R(l11in, 2) = 0.0_dp
733 end do
734 end if
735
736 ! add routed temp energy to downstream node
737 self%netNode_E_IN(l11to, 2) = self%netNode_E_IN(l11to, 2) + self%netNode_E_R(l11in, 2)
738
739 end do riv_loop
740
741 ! Accumulate all inputs in tNode (netNode_E_out) ONLY for last link
742 tnode = netlink_ton(netperm(nlinks))
743 l11to = tnode + self%s11 - 1_i4
744 self%netNode_E_IN(l11to, 2) = self%netNode_E_IN(l11to, 2) + self%netNode_E_out(l11to)
745 ! calculate riv temp at last link
746 self%river_temp(l11to) = self%netNode_E_IN(l11to, 2) / l11_qmod(tnode) - t0_dp
747
748 ! backflow t-> t-1
749 self%netNode_E_R(self%s11 : self%e11, 1) = self%netNode_E_R(self%s11 : self%e11, 2)
750 self%netNode_E_IN(self%s11 : self%e11, 1) = self%netNode_E_IN(self%s11 : self%e11, 2)
751
752 end subroutine l11_routing_e
753
754 !> \brief initialize iterative solver of \ref riv_temp_type
755 subroutine init_iter(self)
756
757 implicit none
758
759 class(riv_temp_type), intent(inout) :: self
760
761 ! first iteration is used to determine the search direction for the interval
762 self%first_iter = .true.
763 ! after the interval is found, we will perform a biseciton to nail down the temperature
764 self%bisect_iter = .false.
765
766 end subroutine init_iter
767
768 !> \brief execute next iteration with iterative solver of \ref riv_temp_type
769 subroutine next_iter(self, T_est, T_rout)
770
771 implicit none
772
773 class(riv_temp_type), intent(inout) :: self
774 real(dp), intent(inout) :: T_est !< estimated river temperature
775 real(dp), intent(in) :: T_rout !< calculated (routed) river temperature
776
777 ! before performing a bisection we need to search for the interval (with given step-size)
778 if ( .not. self%bisect_iter ) then
779 ! to determine the direction of interval search, we have to wait for the first iteration
780 if ( self%first_iter ) then
781 self%first_iter = .false.
782 ! determine search direction
783 self%up_iter = ( t_rout > t_est )
784 end if
785 ! check if we need to take another step
786 if ( self%up_iter .eqv. ( t_rout > t_est ) ) then
787 if ( self%up_iter ) then
788 t_est = t_est + self%step_iter
789 else
790 t_est = t_est - self%step_iter
791 end if
792 ! if direction changes, we have found the interval for bisection
793 else
794 ! start bisection in next iteration
795 self%bisect_iter = .true.
796 ! set interval bounds for bisection
797 if ( self%up_iter ) then
798 self%up_bnd_iter = t_est
799 self%low_bnd_iter = t_est - self%step_iter
800 else
801 self%up_bnd_iter = t_est + self%step_iter
802 self%low_bnd_iter = t_est
803 end if
804 ! set estimation to interval center
805 t_est = (self%up_bnd_iter + self%low_bnd_iter) / 2.0_dp
806 end if
807 ! perform bisection
808 else
809 ! if routed temp is lower than estimated, use lower interval
810 if ( t_rout < t_est ) then
811 self%up_bnd_iter = t_est
812 ! if routed temp is higher than estimated, use upper interval
813 else
814 self%low_bnd_iter = t_est
815 end if
816 ! set estimation to interval center
817 t_est = (self%up_bnd_iter + self%low_bnd_iter) / 2.0_dp
818 end if
819
820 end subroutine next_iter
821
822end module mo_mrm_riv_temp_class
Input checking routines.
Definition mo_check.f90:12
subroutine, public check_dir(path, text, raise, tab, text_length)
Check if a given directory exists.
Definition mo_check.f90:31
Provides constants commonly used by mHM, mRM and MPR.
integer(i4), parameter, public maxnodomains
integer(i4), parameter, public nodata_i4
Provides structures needed by mHM, mRM and/or mpr.
type(domain_meta), public domainmeta
Provides mHM specific constants.
real(dp), parameter, public h2odens
Provides mRM specific constants.
integer(i4), parameter, public nroutingstates
Performs pre-processing for routing for mHM at level L11.
subroutine, public l11_meteo_acc(meteo_all, efecarea, l1_l11_id, l11_areacell, l11_l1_id, map_flag, meteo_acc)
meteo forcing accumulation at L11 for temperature routing.
subroutine, public l11_runoff_acc(qall, efecarea, l1_l11_id, l11_areacell, l11_l1_id, ts, map_flag, qacc)
total runoff accumulation at L11.
subroutine, public calc_l1_runoff_e(fsealed_area_fraction, fast_interflow, slow_interflow, baseflow, direct_runoff, temp_air, mean_temp_air, lateral_e)
calculate lateral temperature energy from runoff components.
Class for the river temperature calculations.
real(dp) function get_sens_heat(self, air_temp, riv_temp)
sensible heat flux of riv_temp_type
subroutine config(self, file_namelist, unamelist, file_namelist_param, unamelist_param)
configure the riv_temp_type class from the mhm namelist
real(dp) function get_lat_heat(self, air_temp, netrad)
latent heat flux of riv_temp_type
subroutine init_riv_temp(self, temp_air, efecarea, l1_l11_id, l11_areacell, l11_l1_id, map_flag)
initialize the river temperature of riv_temp_type class for the current domain
real(dp) function get_lrad_out(self, riv_temp)
get outgoing longwave radiation of riv_temp_type
subroutine dealloc_lateral(self)
deallocate lateral temp components of riv_temp_type
subroutine finalize_source_e(self, efecarea, l1_l11_id, l11_areacell, l11_l1_id, timestep, map_flag)
finalize energy sources of riv_temp_type
subroutine alloc_lateral(self, ncells)
allocate lateral temp components of riv_temp_type class for current domain
subroutine clean_up(self)
clean up
subroutine next_iter(self, t_est, t_rout)
execute next iteration with iterative solver of riv_temp_type
subroutine init_iter(self)
initialize iterative solver of riv_temp_type
subroutine init_area(self, idomain, l11_netperm, l11_fromn, l11_length, nlinks, ncells, nrows, ncols, l11_mask)
initialize the river area of riv_temp_type class for the current domain
real(dp) function get_e_io(self, riv_temp, cell)
get complete energy source of riv_temp_type at given cell
subroutine acc_source_e(self, fsealed_area_fraction, fast_interflow, slow_interflow, baseflow, direct_runoff, temp_air)
accumulate energy sources of riv_temp_type
subroutine reset_timestep(self)
reset riv_temp_type class for next timestep
subroutine init(self, ncells)
initalize the riv_temp_type class for the current domain
subroutine l11_routing_e(self, nlinks, netperm, netlink_fromn, netlink_ton, netlink_c1, netlink_c2, ninflowgauges, inflowheadwater, inflownodelist, l11_qtr, l11_qmod)
execute the temperature routing of riv_temp_type
Module for calculating reference/potential evapotranspiration [mm d-1].
Definition mo_pet.f90:15
elemental pure real(dp) function, public slope_satpressure(tavg)
slope of saturation vapour pressure curve
Definition mo_pet.f90:385
Reads forcing input data.
subroutine, public read_const_nc(folder, nrows, ncols, varname, data, filename)
Reads time independent forcing input in NetCDF file format.
This is a container to define the river temperature routing in the current time step.