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