Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
steady_simcomp.f90
1! Copyright (c) 2024, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32
35 use simulation_component, only: simulation_component_t
36 use num_types, only: rp, dp
37 use field, only: field_t
38 use json_module, only: json_file
39 use json_utils, only: json_get_or_default
40 use case, only: case_t
41 use field_math, only: field_sub2, field_copy
42 use math, only: glsc3
43 use device_math, only: device_glsc3
44 use coefs, only: coef_t
45 use neko_config, only : neko_bcknd_device
46 use csv_file, only : csv_file_t
47 use vector, only: vector_t
48 use time_state, only: time_state_t
49 implicit none
50 private
51
56 type, public, extends(simulation_component_t) :: steady_simcomp_t
57 private
58
59 ! Old fields
60 type(field_t) :: u_old, v_old, w_old, p_old, s_old
61 real(kind=dp) :: tol
62
63 logical :: have_scalar = .false.
64
66 type(csv_file_t) :: logger
68 integer :: log_frequency
69 ! vector to store the data being logged
70 type(vector_t) :: log_data
71
72 contains
73 ! Constructor from json, wrapping the actual constructor.
74 procedure, public, pass(this) :: init => steady_simcomp_init_from_json
75 ! Actual constructor.
76 procedure, public, pass(this) :: init_from_attributes => &
77 steady_simcomp_init_from_attributes
78 ! Destructor.
79 procedure, public, pass(this) :: free => steady_simcomp_free
80 ! Compute the steady_simcomp field.
81 procedure, public, pass(this) :: compute_ => steady_simcomp_compute
82 end type steady_simcomp_t
83
84contains
85
86 ! Constructor from json.
87 subroutine steady_simcomp_init_from_json(this, json, case)
88 class(steady_simcomp_t), intent(inout) :: this
89 type(json_file), intent(inout) :: json
90 class(case_t), intent(inout), target :: case
91 real(kind=dp) :: tol
92 integer :: log_frequency
93
94 call this%init_base(json, case)
95
96 ! Read the tolerance
97 call json_get_or_default(json, "tol", tol, 1.0e-6_dp)
98 ! Read the log frequency
99 call json_get_or_default(json, "log_frequency", log_frequency, 50)
100 call json_get_or_default(json, "scalar_coupled", this%have_scalar, .false.)
101
102 call this%init_from_attributes(tol, log_frequency)
103
104 end subroutine steady_simcomp_init_from_json
105
106 ! Actual constructor.
107 subroutine steady_simcomp_init_from_attributes(this, tol, log_frequency)
108 class(steady_simcomp_t), intent(inout) :: this
109 real(kind=dp), intent(in) :: tol
110 integer, intent(in) :: log_frequency
111
112 this%tol = tol
113 this%log_frequency = log_frequency
114
115 ! initialize the logger
116 call this%logger%init('steady_state_data.csv')
117 call this%logger%set_header('iter,time,u,v,w,p,t')
118 call this%log_data%init(7)
119
120 ! Point local fields to the scratch fields
121 call this%u_old%init(this%case%fluid%u%dof)
122 call this%v_old%init(this%case%fluid%v%dof)
123 call this%w_old%init(this%case%fluid%w%dof)
124 call this%p_old%init(this%case%fluid%p%dof)
125
126 ! Check if the scalar field is allocated
127 if (this%have_scalar) then
128 call this%s_old%init(this%case%scalar%s%dof)
129 end if
130
131 end subroutine steady_simcomp_init_from_attributes
132
133 ! Destructor.
134 subroutine steady_simcomp_free(this)
135 class(steady_simcomp_t), intent(inout) :: this
136
137 call this%u_old%free()
138 call this%v_old%free()
139 call this%w_old%free()
140 call this%p_old%free()
141 call this%s_old%free()
142 call this%log_data%free()
143
144 call this%free_base()
145
146 end subroutine steady_simcomp_free
147
148 ! Compute the steady_simcomp field.
149 subroutine steady_simcomp_compute(this, time)
150 class(steady_simcomp_t), intent(inout) :: this
151 type(time_state_t), intent(in) :: time
152 integer :: tstep
153 real(kind=rp) :: t
154 real(kind=rp) :: dt
155
156 real(kind=rp), dimension(5) :: normed_diff
157 type(field_t), pointer :: u, v, w, p, s
158
159 ! A frozen field is not interesting to compute differences for.
160 if (this%case%fluid%freeze) return
161
162 ! Get the current time and time step
163 t = time%t
164 tstep = time%tstep
165 dt = this%case%time%dt
166
167 ! ------------------------------------------------------------------------ !
168 ! Computation of the maximal normed difference.
169 !
170 ! Our goal is to freeze the simulation when the normed difference between
171 ! the old and new fields is below a certain tolerance.
172
173 u => this%case%fluid%u
174 v => this%case%fluid%v
175 w => this%case%fluid%w
176 p => this%case%fluid%p
177
178 if (this%have_scalar) then
179 s => this%case%scalar%s
180 else
181 s => null()
182 end if
183
184 ! Compute the difference between the old and new fields
185 call field_sub2(this%u_old, u)
186 call field_sub2(this%v_old, v)
187 call field_sub2(this%w_old, w)
188 call field_sub2(this%p_old, p)
189 if (this%have_scalar) then
190 call field_sub2(this%s_old, s)
191 end if
192
193 ! Here we compute the squared difference between the old and new fields
194 ! and store the result in the `normed_diff` array.
195 normed_diff(1) = energy_norm(this%u_old, this%case%fluid%C_Xh, dt)
196 normed_diff(2) = energy_norm(this%v_old, this%case%fluid%C_Xh, dt)
197 normed_diff(3) = energy_norm(this%w_old, this%case%fluid%C_Xh, dt)
198 normed_diff(4) = energy_norm(this%p_old, this%case%fluid%C_Xh, dt)
199 if (this%have_scalar) then
200 normed_diff(5) = energy_norm(this%s_old, this%case%fluid%C_Xh, dt)
201 else
202 normed_diff(5) = 0.0_rp
203 end if
204
205 ! If the normed difference is below the tolerance, we consider the
206 ! simulation to have converged. Otherwise, we copy the new fields to the
207 ! old fields and continue the simulation.
208 if (maxval(normed_diff) .gt. this%tol) then
209 call field_copy(this%u_old, u)
210 call field_copy(this%v_old, v)
211 call field_copy(this%w_old, w)
212 call field_copy(this%p_old, p)
213 if (this%have_scalar) then
214 call field_copy(this%s_old, s)
215 end if
216
217 ! this could be a csv, but I think it's nicer in the logfile, it can
218 ! always be greped out later if you want to plot.
219 if (mod(tstep, this%log_frequency) .eq. 0) then
220 this%log_data%x(1) = real(tstep, kind=rp)
221 this%log_data%x(2) = t
222 this%log_data%x(3) = normed_diff(1)
223 this%log_data%x(4) = normed_diff(2)
224 this%log_data%x(5) = normed_diff(3)
225 this%log_data%x(6) = normed_diff(4)
226 this%log_data%x(7) = normed_diff(5)
227 call this%logger%write(this%log_data)
228 end if
229
230 else
231 this%case%fluid%freeze = .true.
232 if (this%have_scalar) then
233 ! @todo
234 ! we should implement a freeze for the scalar too
235 !this%case%scalar%freeze = .true.
236 end if
237 end if
238
239
240 end subroutine steady_simcomp_compute
241
246 function energy_norm(delta_fld, coef, dt)
247 type(field_t), intent(in) :: delta_fld
248 type(coef_t), intent(in) :: coef
249 real(kind=rp), intent(in) :: dt
250 real(kind=rp) :: energy_norm, tmp
251 integer :: n
252
253 n = delta_fld%size()
254 if (neko_bcknd_device .eq. 1) then
255 tmp = device_glsc3(delta_fld%x_d, delta_fld%x_d, coef%B_d, n)
256 else
257 tmp = glsc3(delta_fld%x, delta_fld%x, coef%B, n)
258 end if
259 energy_norm = sqrt(tmp) / dt / coef%volume
260
261 end function energy_norm
262
263
264end module steady_simcomp
265
Implements the steady_simcomp_t type.
The steady_simcomp_t type is a simulation component that terminates a simulation when the normed diff...