Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
steady_simcomp.f90
Go to the documentation of this file.
1
34
37 use simulation_component, only: simulation_component_t
38 use num_types, only: rp, dp
39 use field, only: field_t
40 use json_module, only: json_file
41 use json_utils, only: json_get_or_default
42 use case, only: case_t
43 use field_math, only: field_sub2, field_copy
44 use math, only: glsc3
45 use device_math, only: device_glsc3
46 use coefs, only: coef_t
47 use neko_config, only : neko_bcknd_device
48 use csv_file, only : csv_file_t
49 use vector, only: vector_t
50 use time_state, only: time_state_t
51 use utils, only: neko_error
52 implicit none
53 private
55
60 type, public, extends(simulation_component_t) :: steady_simcomp_t
61 private
62
63 ! Old fields
64 type(field_t) :: u_old, v_old, w_old, p_old, s_old
65 real(kind=dp) :: tol
66
67 logical :: have_scalar = .false.
68
70 type(csv_file_t) :: logger
72 integer :: log_frequency
73 ! vector to store the data being logged
74 type(vector_t) :: log_data
75
76 contains
77 ! Constructor from json, wrapping the actual constructor.
78 procedure, public, pass(this) :: init => steady_simcomp_init_from_json
79 ! Actual constructor.
80 procedure, public, pass(this) :: init_from_attributes => &
81 steady_simcomp_init_from_attributes
82 ! Destructor.
83 procedure, public, pass(this) :: free => steady_simcomp_free
84 ! Compute the steady_simcomp field.
85 procedure, public, pass(this) :: compute_ => steady_simcomp_compute
86 end type steady_simcomp_t
87
88contains
89
92 class(simulation_component_t), allocatable, intent(inout) :: obj
93 allocate(steady_simcomp_t::obj)
94 end subroutine steady_simcomp_allocate
95
96 ! Constructor from json.
97 subroutine steady_simcomp_init_from_json(this, json, case)
98 class(steady_simcomp_t), intent(inout), target :: this
99 type(json_file), intent(inout) :: json
100 class(case_t), intent(inout), target :: case
101 real(kind=dp) :: tol
102 integer :: log_frequency
103
104 call this%init_base(json, case)
105
106 ! Read the tolerance
107 call json_get_or_default(json, "tol", tol, 1.0e-6_dp)
108 ! Read the log frequency
109 call json_get_or_default(json, "log_frequency", log_frequency, 50)
110 call json_get_or_default(json, "scalar_coupled", this%have_scalar, .false.)
111
112 call this%init_from_attributes(tol, log_frequency)
113
114 end subroutine steady_simcomp_init_from_json
115
116 ! Actual constructor.
117 subroutine steady_simcomp_init_from_attributes(this, tol, log_frequency)
118 class(steady_simcomp_t), intent(inout) :: this
119 real(kind=dp), intent(in) :: tol
120 integer, intent(in) :: log_frequency
121
122 this%tol = tol
123 this%log_frequency = log_frequency
124
125 ! initialize the logger
126 call this%logger%init('steady_state_data.csv')
127 call this%logger%set_header('iter,time,u,v,w,p,t')
128 call this%log_data%init(7)
129
130 ! Point local fields to the scratch fields
131 call this%u_old%init(this%case%fluid%u%dof)
132 call this%v_old%init(this%case%fluid%v%dof)
133 call this%w_old%init(this%case%fluid%w%dof)
134 call this%p_old%init(this%case%fluid%p%dof)
135
136 ! Check if the scalar field is allocated
137 if (this%have_scalar) then
138 call this%s_old%init(this%case%scalars%scalar_fields(1)%s%dof)
139 end if
140
141 end subroutine steady_simcomp_init_from_attributes
142
143 ! Destructor.
144 subroutine steady_simcomp_free(this)
145 class(steady_simcomp_t), intent(inout) :: this
146
147 call this%u_old%free()
148 call this%v_old%free()
149 call this%w_old%free()
150 call this%p_old%free()
151 call this%s_old%free()
152 call this%log_data%free()
153
154 call this%free_base()
155
156 end subroutine steady_simcomp_free
157
158 ! Compute the steady_simcomp field.
159 subroutine steady_simcomp_compute(this, time)
160 class(steady_simcomp_t), intent(inout) :: this
161 type(time_state_t), intent(in) :: time
162 integer :: tstep
163 real(kind=rp) :: t
164 real(kind=rp) :: dt
165
166 real(kind=rp), dimension(5) :: normed_diff
167 type(field_t), pointer :: u, v, w, p, s
168
169 ! A frozen field is not interesting to compute differences for.
170 if (this%case%fluid%freeze) return
171
172 ! Get the current time and time step
173 t = time%t
174 tstep = time%tstep
175 dt = this%case%time%dt
176
177 ! ------------------------------------------------------------------------ !
178 ! Computation of the maximal normed difference.
179 !
180 ! Our goal is to freeze the simulation when the normed difference between
181 ! the old and new fields is below a certain tolerance.
182
183 u => this%case%fluid%u
184 v => this%case%fluid%v
185 w => this%case%fluid%w
186 p => this%case%fluid%p
187
188 if (this%have_scalar) then
189 if (size(this%case%scalars%scalar_fields) .gt. 1) then
190 call neko_error('steady simcomp only works for a single scalar')
191 end if
192 s => this%case%scalars%scalar_fields(1)%s
193 else
194 s => null()
195 end if
196
197 ! Compute the difference between the old and new fields
198 call field_sub2(this%u_old, u)
199 call field_sub2(this%v_old, v)
200 call field_sub2(this%w_old, w)
201 call field_sub2(this%p_old, p)
202 if (this%have_scalar) then
203 call field_sub2(this%s_old, s)
204 end if
205
206 ! Here we compute the squared difference between the old and new fields
207 ! and store the result in the `normed_diff` array.
208 normed_diff(1) = energy_norm(this%u_old, this%case%fluid%C_Xh, dt)
209 normed_diff(2) = energy_norm(this%v_old, this%case%fluid%C_Xh, dt)
210 normed_diff(3) = energy_norm(this%w_old, this%case%fluid%C_Xh, dt)
211 normed_diff(4) = energy_norm(this%p_old, this%case%fluid%C_Xh, dt)
212 if (this%have_scalar) then
213 normed_diff(5) = energy_norm(this%s_old, this%case%fluid%C_Xh, dt)
214 else
215 normed_diff(5) = 0.0_rp
216 end if
217
218 ! If the normed difference is below the tolerance, we consider the
219 ! simulation to have converged. Otherwise, we copy the new fields to the
220 ! old fields and continue the simulation.
221 if (maxval(normed_diff) .gt. this%tol) then
222 call field_copy(this%u_old, u)
223 call field_copy(this%v_old, v)
224 call field_copy(this%w_old, w)
225 call field_copy(this%p_old, p)
226 if (this%have_scalar) then
227 call field_copy(this%s_old, s)
228 end if
229
230 ! this could be a csv, but I think it's nicer in the logfile, it can
231 ! always be greped out later if you want to plot.
232 if (mod(tstep, this%log_frequency) .eq. 0) then
233 this%log_data%x(1) = real(tstep, kind=rp)
234 this%log_data%x(2) = t
235 this%log_data%x(3) = normed_diff(1)
236 this%log_data%x(4) = normed_diff(2)
237 this%log_data%x(5) = normed_diff(3)
238 this%log_data%x(6) = normed_diff(4)
239 this%log_data%x(7) = normed_diff(5)
240 call this%logger%write(this%log_data)
241 end if
242
243 else
244 this%case%fluid%freeze = .true.
245 if (this%have_scalar) then
246 ! @todo
247 ! we should implement a freeze for the scalar too
248 !this%case%scalar%freeze = .true.
249 end if
250 end if
251
252
253 end subroutine steady_simcomp_compute
254
259 function energy_norm(delta_fld, coef, dt)
260 type(field_t), intent(in) :: delta_fld
261 type(coef_t), intent(in) :: coef
262 real(kind=rp), intent(in) :: dt
263 real(kind=rp) :: energy_norm, tmp
264 integer :: n
265
266 n = delta_fld%size()
267 if (neko_bcknd_device .eq. 1) then
268 tmp = device_glsc3(delta_fld%x_d, delta_fld%x_d, coef%B_d, n)
269 else
270 tmp = glsc3(delta_fld%x, delta_fld%x, coef%B, n)
271 end if
272 energy_norm = sqrt(tmp) / dt / coef%volume
273
274 end function energy_norm
275
276
277end module steady_simcomp
Implements the steady_simcomp_t type.
subroutine, public steady_simcomp_allocate(obj)
Allocator for the steady simulation component.
The steady_simcomp_t type is a simulation component that terminates a simulation when the normed diff...