Line data Source code
1 : !> \file mo_optimization_types.f90
2 : !> \brief \copybrief mo_optimization_types
3 : !> \details \copydetails mo_optimization_types
4 :
5 : !> \brief Type definitions for optimization routines
6 : !> \author Maren Kaluza
7 : !> \date Nov 2019
8 : !> \copyright Copyright 2005-\today, the CHS Developers, Sabine Attinger: All rights reserved.
9 : !! FORCES is released under the LGPLv3+ license \license_note
10 : MODULE mo_optimization_types
11 : use mo_kind, only : i4, dp
12 :
13 : IMPLICIT NONE
14 :
15 : public :: optidata, optidata_sim
16 :
17 : private
18 :
19 : !> \brief optional data, such as sm, neutrons, et, tws
20 : !> \details data type for observed data, providing metadata
21 : !! for simulated data
22 : !! dim1 = number grid cells L1
23 : !! dim2 = number of meteorological time steps
24 : type optidata
25 : real(dp), dimension(:, :), allocatable :: dataObs !< observed data
26 : logical, dimension(:, :), allocatable :: maskObs !< mask of observed data
27 : character(256) :: dir !< directory where to read opti data
28 : integer(i4) :: timeStepInput !< time step of optional data
29 : character(256) :: varname !< variable name
30 : end type optidata
31 :
32 : !> \brief type for simulated optional data
33 : type optidata_sim
34 : real(dp), dimension(:, :), allocatable :: dataSim !< simulation data
35 : integer(i4) :: averageTimestep !< the current timestep
36 : !< the simulated opti data is written to
37 : integer(i4) :: averageCounter !< set to 0 on average, incremented on add
38 :
39 : contains
40 : procedure :: init => optidata_sim_init
41 : procedure :: destroy => optidata_sim_destroy
42 : procedure :: increment_counter => optidata_sim_increment_counter
43 : procedure :: add => optidata_sim_add
44 : procedure :: average => optidata_sim_average
45 : procedure :: average_per_timestep => optidata_sim_average_per_timestep
46 : procedure :: average_add => optidata_sim_average_add
47 : end type optidata_sim
48 :
49 : contains
50 :
51 0 : subroutine optidata_sim_init(this, optidataObs)
52 : class(optidata_sim), intent(inout) :: this
53 : type(optidata), intent(in) :: optidataObs
54 :
55 0 : allocate(this%dataSim(size(optidataObs%dataObs, dim = 1), size(optidataObs%dataObs, dim = 2)))
56 0 : this%dataSim(:, :) = 0.0_dp ! has to be intialized with zero because later summation
57 0 : this%averageTimestep = 1
58 0 : this%averageCounter = 0
59 0 : end subroutine optidata_sim_init
60 :
61 0 : subroutine optidata_sim_destroy(this)
62 : class(optidata_sim), intent(inout) :: this
63 :
64 0 : deallocate(this%dataSim)
65 0 : end subroutine optidata_sim_destroy
66 :
67 0 : subroutine optidata_sim_increment_counter(this, timeStepInput, is_new_day, is_new_month, is_new_year)
68 : class(optidata_sim), intent(inout) :: this
69 : integer(i4), intent(in) :: timeStepInput
70 : logical, intent(in) :: is_new_day
71 : logical, intent(in) :: is_new_month
72 : logical, intent(in) :: is_new_year
73 :
74 0 : select case(timeStepInput)
75 : case(-1) ! daily
76 0 : if (is_new_day) then
77 0 : this%averageTimestep = this%averageTimestep + 1
78 : end if
79 : case(-2) ! monthly
80 0 : if (is_new_month) then
81 0 : this%averageTimestep = this%averageTimestep + 1
82 : end if
83 : case(-3) ! yearly
84 0 : if (is_new_year) then
85 0 : this%averageTimestep = this%averageTimestep + 1
86 : end if
87 : end select
88 :
89 0 : end subroutine optidata_sim_increment_counter
90 :
91 0 : subroutine optidata_sim_add(this, data_sim)
92 : class(optidata_sim), intent(inout) :: this
93 : real(dp), dimension(:), intent(in) :: data_sim
94 :
95 : this%dataSim(:, this%averageTimestep) = &
96 0 : this%dataSim(:, this%averageTimestep) + data_sim(:)
97 0 : end subroutine optidata_sim_add
98 :
99 0 : subroutine optidata_sim_average(this)
100 : class(optidata_sim), intent(inout) :: this
101 :
102 : this%dataSim(:, this%averageTimestep) = &
103 0 : this%dataSim(:, this%averageTimestep) / real(this%averageCounter, dp)
104 0 : this%averageTimestep = this%averageTimestep + 1
105 0 : this%averageCounter = 0
106 0 : end subroutine optidata_sim_average
107 :
108 0 : subroutine optidata_sim_average_per_timestep(this, timeStepInput, is_new_day, is_new_month, is_new_year)
109 : class(optidata_sim), intent(inout) :: this
110 : integer(i4), intent(in) :: timeStepInput
111 : logical, intent(in) :: is_new_day
112 : logical, intent(in) :: is_new_month
113 : logical, intent(in) :: is_new_year
114 :
115 0 : select case(timeStepInput)
116 : case(-1) ! daily
117 0 : if (is_new_day) then
118 0 : call this%average()
119 : end if
120 : case(-2) ! monthly
121 0 : if (is_new_month) then
122 0 : call this%average()
123 : end if
124 : case(-3) ! yearly
125 0 : if (is_new_year) then
126 0 : call this%average()
127 : end if
128 : end select
129 0 : end subroutine optidata_sim_average_per_timestep
130 :
131 0 : subroutine optidata_sim_average_add(this, data_sim)
132 : class(optidata_sim), intent(inout) :: this
133 : real(dp), dimension(:), intent(in) :: data_sim
134 :
135 0 : call this%add(data_sim(:))
136 0 : this%averageCounter = this%averageCounter + 1
137 0 : end subroutine optidata_sim_average_add
138 :
139 0 : END MODULE mo_optimization_types
|