Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
fluid_scheme_compressible.f90
Go to the documentation of this file.
1! Copyright (c) 2025, 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!
34 use field, only : field_t
37
38 use registry, only : neko_registry
40 use json_module, only : json_file
41 use num_types, only : rp
42 use mesh, only : mesh_t
44 use space, only : gll
45 use user_intf, only : user_t
47 use mpi_f08
48 use operators, only : cfl_compressible
54 use time_state, only : time_state_t
55 use logger, only : neko_log, log_size
56 use math, only : glsum
57 use num_types, only : i8
58 implicit none
59 private
60
62 type, public, abstract, extends(fluid_scheme_base_t) :: fluid_scheme_compressible_t
64 type(field_t), pointer :: m_x => null()
65 type(field_t), pointer :: m_y => null()
66 type(field_t), pointer :: m_z => null()
67 type(field_t), pointer :: e => null()
68 type(field_t), pointer :: max_wave_speed => null()
69 type(field_t), pointer :: s => null()
70
71 real(kind=rp) :: gamma
72
74 integer(kind=i8) :: glb_n_points
76 integer(kind=i8) :: glb_unique_points
77
78 contains
80 procedure, pass(this) :: scheme_init => fluid_scheme_compressible_init
82 procedure, pass(this) :: scheme_free => fluid_scheme_compressible_free
83
85 procedure, pass(this) :: validate => fluid_scheme_compressible_validate
87 procedure, pass(this) :: compute_cfl &
90 procedure, pass(this) :: update_material_properties => &
93 procedure, pass(this) :: compute_entropy => &
96 procedure, pass(this) :: compute_max_wave_speed => &
99 procedure, pass(this) :: log_solver_info => &
101
103
104contains
112 subroutine fluid_scheme_compressible_init(this, msh, lx, params, scheme, user)
113 class(fluid_scheme_compressible_t), target, intent(inout) :: this
114 type(mesh_t), target, intent(inout) :: msh
115 integer, intent(in) :: lx
116 character(len=*), intent(in) :: scheme
117 type(json_file), target, intent(inout) :: params
118 type(user_t), target, intent(in) :: user
119
120 !
121 ! SEM simulation fundamentals
122 !
123
124 this%msh => msh
125
126 if (msh%gdim .eq. 2) then
127 call this%Xh%init(gll, lx, lx)
128 else
129 call this%Xh%init(gll, lx, lx, lx)
130 end if
131
132 call this%dm_Xh%init(msh, this%Xh)
133
134 call this%gs_Xh%init(this%dm_Xh)
135
136 call this%c_Xh%init(this%gs_Xh)
137
138 ! Assign Dofmap to scratch registry
139 call neko_scratch_registry%set_dofmap(this%dm_Xh)
140
141 ! Case parameters
142 this%params => params
143
144 ! Assign a name
145 call json_get_or_default(params, 'case.fluid.name', this%name, "fluid")
146
147 ! Fill mu and rho field with the physical value
148 call neko_registry%add_field(this%dm_Xh, this%name // "_mu")
149 call neko_registry%add_field(this%dm_Xh, this%name // "_rho")
150 this%mu => neko_registry%get_field(this%name // "_mu")
151 this%rho => neko_registry%get_field(this%name // "_rho")
152 call field_cfill(this%mu, 0.0_rp, this%mu%size())
153
154 ! Assign momentum fields
155 call neko_registry%add_field(this%dm_Xh, "m_x")
156 call neko_registry%add_field(this%dm_Xh, "m_y")
157 call neko_registry%add_field(this%dm_Xh, "m_z")
158 this%m_x => neko_registry%get_field("m_x")
159 this%m_y => neko_registry%get_field("m_y")
160 this%m_z => neko_registry%get_field("m_z")
161 call this%m_x%init(this%dm_Xh, "m_x")
162 call this%m_y%init(this%dm_Xh, "m_y")
163 call this%m_z%init(this%dm_Xh, "m_z")
164
165 ! Assign energy field
166 call neko_registry%add_field(this%dm_Xh, "E")
167 this%E => neko_registry%get_field("E")
168 call this%E%init(this%dm_Xh, "E")
169
170 ! Assign maximum wave speed field
171 call neko_registry%add_field(this%dm_Xh, "max_wave_speed")
172 this%max_wave_speed => neko_registry%get_field("max_wave_speed")
173 call this%max_wave_speed%init(this%dm_Xh, "max_wave_speed")
174
175 ! Assign entropy field
176 call neko_registry%add_field(this%dm_Xh, "S")
177 this%S => neko_registry%get_field("S")
178 call this%S%init(this%dm_Xh, "S")
179
180 ! ! Assign velocity fields
181 call neko_registry%add_field(this%dm_Xh, "u")
182 call neko_registry%add_field(this%dm_Xh, "v")
183 call neko_registry%add_field(this%dm_Xh, "w")
184 this%u => neko_registry%get_field("u")
185 this%v => neko_registry%get_field("v")
186 this%w => neko_registry%get_field("w")
187 call this%u%init(this%dm_Xh, "u")
188 call this%v%init(this%dm_Xh, "v")
189 call this%w%init(this%dm_Xh, "w")
190 call neko_registry%add_field(this%dm_Xh, 'p')
191 this%p => neko_registry%get_field('p')
192 call this%p%init(this%dm_Xh, "p")
193
194 !
195 ! Setup right-hand side fields.
196 !
197 allocate(this%f_x)
198 allocate(this%f_y)
199 allocate(this%f_z)
200 call this%f_x%init(this%dm_Xh, fld_name = "fluid_rhs_x")
201 call this%f_y%init(this%dm_Xh, fld_name = "fluid_rhs_y")
202 call this%f_z%init(this%dm_Xh, fld_name = "fluid_rhs_z")
203
204 ! Compressible parameters
205 call json_get_or_default(params, 'case.fluid.gamma', this%gamma, 1.4_rp)
206
207 ! Calculate global points for logging
208 this%glb_n_points = int(this%msh%glb_nelv, i8)*int(this%Xh%lxyz, i8)
209 this%glb_unique_points = int(glsum(this%c_Xh%mult, this%dm_Xh%size()), i8)
210
211 !
212 ! Log solver information
213 !
214 call this%log_solver_info(params, scheme, lx)
215 end subroutine fluid_scheme_compressible_init
216
220 class(fluid_scheme_compressible_t), intent(inout) :: this
221 call this%dm_Xh%free()
222 call this%gs_Xh%free()
223 call this%c_Xh%free()
224 call this%Xh%free()
225
226 if (associated(this%m_x)) then
227 call this%m_x%free()
228 end if
229
230 if (associated(this%m_y)) then
231 call this%m_y%free()
232 end if
233
234 if (associated(this%m_z)) then
235 call this%m_z%free()
236 end if
237
238 if (associated(this%E)) then
239 call this%E%free()
240 end if
241
242 if (associated(this%max_wave_speed)) then
243 call this%max_wave_speed%free()
244 end if
245
246 if (associated(this%S)) then
247 call this%S%free()
248 end if
249
250 if (associated(this%f_x)) then
251 call this%f_x%free()
252 deallocate(this%f_x)
253 end if
254
255 if (associated(this%f_y)) then
256 call this%f_y%free()
257 deallocate(this%f_y)
258 end if
259
260 if (associated(this%f_z)) then
261 call this%f_z%free()
262 deallocate(this%f_z)
263 end if
264
265 nullify(this%f_x)
266 nullify(this%f_y)
267 nullify(this%f_z)
268
269 nullify(this%m_x)
270 nullify(this%m_y)
271 nullify(this%m_z)
272 nullify(this%E)
273 nullify(this%max_wave_speed)
274 nullify(this%S)
275
276 nullify(this%u)
277 nullify(this%v)
278 nullify(this%w)
279 nullify(this%p)
280 nullify(this%rho)
281 nullify(this%mu)
282
283 end subroutine fluid_scheme_compressible_free
284
288 class(fluid_scheme_compressible_t), target, intent(inout) :: this
289 integer :: n
290 type(field_t), pointer :: temp
291 integer :: temp_indices(1)
292
293 n = this%dm_Xh%size()
294 call neko_scratch_registry%request_field(temp, temp_indices(1), .false.)
295
297 call field_col3(this%m_x, this%u, this%rho)
298 call field_col3(this%m_y, this%v, this%rho)
299 call field_col3(this%m_z, this%w, this%rho)
300
303 call field_cmult2(this%E, this%p, 1.0_rp/(this%gamma - 1.0_rp), n)
304 call field_col3(temp, this%u, this%u, n)
305 call field_addcol3(temp, this%v, this%v, n)
306 call field_addcol3(temp, this%w, this%w, n)
307 call field_col2(temp, this%rho, n)
308 call field_cmult(temp, 0.5_rp, n)
309 call field_add2(this%E, temp, n)
310
311 call neko_scratch_registry%relinquish_field(temp_indices)
312
314 call this%compute_max_wave_speed()
315
317
322 function fluid_scheme_compressible_compute_cfl(this, dt) result(c)
323 class(fluid_scheme_compressible_t), intent(in) :: this
324 real(kind=rp), intent(in) :: dt
325 real(kind=rp) :: c
326 integer :: n
327
328 associate(u => this%u, v => this%v, w => this%w, p => this%p, &
329 rho => this%rho, xh => this%Xh, c_xh => this%c_Xh, &
330 msh => this%msh, gamma => this%gamma, &
331 max_wave_speed => this%max_wave_speed)
332
333 n = xh%lx * xh%ly * xh%lz * msh%nelv
334
335 ! Use the compressible CFL function with precomputed maximum wave speed
336 c = cfl_compressible(dt, max_wave_speed%x, xh, c_xh, msh%nelv, msh%gdim)
337 end associate
338
340
344 class(fluid_scheme_compressible_t), intent(inout) :: this
345 type(time_state_t), intent(in) :: time
347
351 class(fluid_scheme_compressible_t), intent(inout) :: this
352 integer :: n
353
354 n = this%S%dof%size()
355
357 if (neko_bcknd_device .eq. 1) then
358 call compressible_ops_device_compute_entropy(this%S, this%p, this%rho, this%gamma, n)
359 else
360 call compressible_ops_cpu_compute_entropy(this%S%x, this%p%x, this%rho%x, this%gamma, n)
361 end if
362
364
368 class(fluid_scheme_compressible_t), intent(inout) :: this
369 integer :: n
370
371 n = this%u%dof%size()
372
374 if (neko_bcknd_device .eq. 1) then
375 call compressible_ops_device_compute_max_wave_speed(this%max_wave_speed, &
376 this%u, this%v, this%w, this%gamma, this%p, this%rho, n)
377 else
378 call compressible_ops_cpu_compute_max_wave_speed(this%max_wave_speed%x, &
379 this%u%x, this%v%x, this%w%x, this%gamma, this%p%x, this%rho%x, n)
380 end if
381
383
389 subroutine fluid_scheme_compressible_log_solver_info(this, params, scheme, lx)
390 class(fluid_scheme_compressible_t), intent(inout) :: this
391 type(json_file), intent(inout) :: params
392 character(len=*), intent(in) :: scheme
393 integer, intent(in) :: lx
394 character(len=LOG_SIZE) :: log_buf
395 logical :: logical_val
396 real(kind=rp) :: real_val
397 integer :: integer_val
398
399 call neko_log%section('Fluid')
400 write(log_buf, '(A, A)') 'Type : ', trim(scheme)
401 call neko_log%message(log_buf)
402 write(log_buf, '(A, A)') 'Name : ', trim(this%name)
403 call neko_log%message(log_buf)
404
405 ! Polynomial order
406 if (lx .lt. 10) then
407 write(log_buf, '(A, I1)') 'Poly order : ', lx-1
408 else if (lx .ge. 10) then
409 write(log_buf, '(A, I2)') 'Poly order : ', lx-1
410 else
411 write(log_buf, '(A, I3)') 'Poly order : ', lx-1
412 end if
413 call neko_log%message(log_buf)
414
415 ! Global points information
416 write(log_buf, '(A, I0)') 'GLL points : ', this%glb_n_points
417 call neko_log%message(log_buf)
418 write(log_buf, '(A, I0)') 'Unique pts.: ', this%glb_unique_points
419 call neko_log%message(log_buf)
420
421 ! Material properties
422 write(log_buf, '(A,ES13.6)') 'gamma :', this%gamma
423 call neko_log%message(log_buf)
424
425 ! Compressible-specific parameters
426 call json_get_or_default(params, 'case.numerics.c_avisc_low', real_val, 0.5_rp)
427 write(log_buf, '(A,ES13.6)') 'c_avisc_low:', real_val
428 call neko_log%message(log_buf)
429
430 call json_get_or_default(params, 'case.numerics.time_order', integer_val, 4)
431 write(log_buf, '(A, I0)') 'RK order : ', integer_val
432 call neko_log%message(log_buf)
433
435
Abstract interface to sets rho and mu.
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
CPU implementation of compressible flow operations.
subroutine, public compressible_ops_cpu_compute_entropy(s, p, rho, gamma, n)
Compute entropy field S = 1/(gamma-1) * rho * (log(p) - gamma * log(rho)) on CPU.
subroutine, public compressible_ops_cpu_compute_max_wave_speed(max_wave_speed, u, v, w, gamma, p, rho, n)
Compute maximum wave speed for compressible flows on CPU.
Device implementation of compressible flow operations.
subroutine, public compressible_ops_device_compute_entropy(s, p, rho, gamma, n)
Compute entropy field S = 1/(gamma-1) * rho * (log(p) - gamma * log(rho)) on device.
subroutine, public compressible_ops_device_compute_max_wave_speed(max_wave_speed, u, v, w, gamma, p, rho, n)
Compute maximum wave speed for compressible flows on device.
subroutine, public field_col2(a, b, n)
Vector multiplication .
subroutine, public field_cmult2(a, b, c, n)
Multiplication by constant c .
subroutine, public field_cfill(a, c, n)
Set all elements to a constant c .
subroutine, public field_addcol3(a, b, c, n)
Returns .
subroutine, public field_add2(a, b, n)
Vector addition .
subroutine, public field_col3(a, b, c, n)
Vector multiplication with 3 vectors .
subroutine, public field_cmult(a, c, n)
Multiplication by constant c .
Defines a field.
Definition field.f90:34
subroutine fluid_scheme_compressible_log_solver_info(this, params, scheme, lx)
Log comprehensive solver information.
subroutine fluid_scheme_compressible_free(this)
Free allocated memory and cleanup resources.
subroutine fluid_scheme_compressible_update_material_properties(this, time)
Set rho and mu.
subroutine fluid_scheme_compressible_init(this, msh, lx, params, scheme, user)
Initialize common data for compressible fluid scheme.
real(kind=rp) function fluid_scheme_compressible_compute_cfl(this, dt)
Compute CFL number.
subroutine fluid_scheme_compressible_compute_max_wave_speed(this)
Compute maximum wave speed for compressible flows.
subroutine fluid_scheme_compressible_validate(this)
Validate field initialization and compute derived quantities.
subroutine fluid_scheme_compressible_compute_entropy(this)
Compute entropy field S = 1/(gamma-1) * rho * (log(p) - gamma * log(rho))
Utilities for retrieving parameters from the case files.
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:76
integer, parameter, public log_size
Definition log.f90:46
Definition math.f90:60
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
Definition math.f90:499
Defines a mesh.
Definition mesh.f90:34
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public i8
Definition num_types.f90:7
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators.
Definition operators.f90:34
real(kind=rp) function, public cfl_compressible(dt, max_wave_speed, xh, coef, nelv, gdim)
Defines a registry for storing solution fields.
Definition registry.f90:34
type(registry_t), target, public neko_registry
Global field registry.
Definition registry.f90:128
Defines a registry for storing and requesting temporary objects This can be used when you have a func...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Defines a function space.
Definition space.f90:34
integer, parameter, public gll
Definition space.f90:49
Module with things related to the simulation time.
Interfaces for user interaction with NEKO.
Definition user_intf.f90:34
Base type of all fluid formulations.
A struct that contains all info about the time, expand as needed.
A type collecting all the overridable user routines and flag to suppress type injection from custom m...