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
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
60 type(field_t) :: u_old, v_old, w_old, p_old, s_old
63 logical :: have_scalar = .false.
66 type(csv_file_t) :: logger
68 integer :: log_frequency
70 type(vector_t) :: log_data
74 procedure,
public, pass(this) :: init => steady_simcomp_init_from_json
76 procedure,
public, pass(this) :: init_from_attributes => &
77 steady_simcomp_init_from_attributes
79 procedure,
public, pass(this) :: free => steady_simcomp_free
81 procedure,
public, pass(this) :: compute_ => steady_simcomp_compute
87 subroutine steady_simcomp_init_from_json(this, json, case)
89 type(json_file),
intent(inout) :: json
90 class(case_t),
intent(inout),
target :: case
92 integer :: log_frequency
94 call this%init_base(json, case)
97 call json_get_or_default(json,
"tol", tol, 1.0e-6_dp)
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.)
102 call this%init_from_attributes(tol, log_frequency)
104 end subroutine steady_simcomp_init_from_json
107 subroutine steady_simcomp_init_from_attributes(this, tol, log_frequency)
109 real(kind=dp),
intent(in) :: tol
110 integer,
intent(in) :: log_frequency
113 this%log_frequency = log_frequency
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)
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)
127 if (this%have_scalar)
then
128 call this%s_old%init(this%case%scalar%s%dof)
131 end subroutine steady_simcomp_init_from_attributes
134 subroutine steady_simcomp_free(this)
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()
144 call this%free_base()
146 end subroutine steady_simcomp_free
149 subroutine steady_simcomp_compute(this, time)
151 type(time_state_t),
intent(in) :: time
156 real(kind=rp),
dimension(5) :: normed_diff
157 type(field_t),
pointer :: u, v, w, p, s
160 if (this%case%fluid%freeze)
return
165 dt = this%case%time%dt
173 u => this%case%fluid%u
174 v => this%case%fluid%v
175 w => this%case%fluid%w
176 p => this%case%fluid%p
178 if (this%have_scalar)
then
179 s => this%case%scalar%s
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)
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)
202 normed_diff(5) = 0.0_rp
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)
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)
231 this%case%fluid%freeze = .true.
232 if (this%have_scalar)
then
240 end subroutine steady_simcomp_compute
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
254 if (neko_bcknd_device .eq. 1)
then
255 tmp = device_glsc3(delta_fld%x_d, delta_fld%x_d, coef%B_d, n)
257 tmp = glsc3(delta_fld%x, delta_fld%x, coef%B, n)
259 energy_norm = sqrt(tmp) / dt / coef%volume
261 end function energy_norm
Implements the steady_simcomp_t type.
The steady_simcomp_t type is a simulation component that terminates a simulation when the normed diff...