Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
fluid_scheme_compressible_euler.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 comm, only : neko_comm
35 use advection, only : advection_t
37 use operators, only : div, rotate_cyc
41 use math, only : col2
42 use device_math, only : device_col2
43 use field, only : field_t
47 use gather_scatter, only : gs_t
48 use num_types, only : rp
49 use mesh, only : mesh_t
50 use checkpoint, only : chkp_t
51 use json_module, only : json_file, json_core, json_value
54 use user_intf, only : user_t
56 use ax_product, only : ax_t, ax_helm_factory
57 use coefs, only : coef_t
58 use euler_residual, only : euler_rhs_t, euler_rhs_factory
61 use bc_list, only : bc_list_t
62 use bc, only : bc_t
64 use logger, only : log_size
65 use time_state, only : time_state_t
73 use mpi_f08, only : mpi_allreduce, mpi_integer, mpi_max
74 use regularization, only : regularization_t, regularization_factory
75 implicit none
76 private
77
78 type, public, extends(fluid_scheme_compressible_t) &
80 type(field_t) :: rho_res, m_x_res, m_y_res, m_z_res, m_e_res
81 type(field_t) :: drho, dm_x, dm_y, dm_z, de
82 type(field_t) :: h
83 real(kind=rp) :: c_avisc_low
84 class(advection_t), allocatable :: adv
85 class(ax_t), allocatable :: ax
86 class(euler_rhs_t), allocatable :: euler_rhs
87 type(runge_kutta_time_scheme_t) :: rk_scheme
88
89 class(regularization_t), allocatable :: regularization
90
91 ! List of boundary conditions for velocity
92 type(bc_list_t) :: bcs_density
93 contains
94 procedure, pass(this) :: init => fluid_scheme_compressible_euler_init
95 procedure, pass(this) :: free => fluid_scheme_compressible_euler_free
96 procedure, pass(this) :: step => fluid_scheme_compressible_euler_step
97 procedure, pass(this) :: restart => fluid_scheme_compressible_euler_restart
99 procedure, pass(this) :: setup_bcs &
101 procedure, pass(this) :: compute_h
102 procedure, pass(this), private :: setup_regularization
104
105 interface
106
113 module subroutine density_bc_factory(object, scheme, json, coef, user)
114 class(bc_t), pointer, intent(inout) :: object
115 type(fluid_scheme_compressible_euler_t), intent(in) :: scheme
116 type(json_file), intent(inout) :: json
117 type(coef_t), intent(in) :: coef
118 type(user_t), intent(in) :: user
119 end subroutine density_bc_factory
120 end interface
121
122 interface
123
130 module subroutine pressure_bc_factory(object, scheme, json, coef, user)
131 class(bc_t), pointer, intent(inout) :: object
132 type(fluid_scheme_compressible_euler_t), intent(inout) :: scheme
133 type(json_file), intent(inout) :: json
134 type(coef_t), intent(in) :: coef
135 type(user_t), intent(in) :: user
136 end subroutine pressure_bc_factory
137 end interface
138
139 interface
140
147 module subroutine velocity_bc_factory(object, scheme, json, coef, user)
148 class(bc_t), pointer, intent(inout) :: object
149 type(fluid_scheme_compressible_euler_t), intent(in) :: scheme
150 type(json_file), intent(inout) :: json
151 type(coef_t), intent(in) :: coef
152 type(user_t), intent(in) :: user
153 end subroutine velocity_bc_factory
154 end interface
155
156contains
164 subroutine fluid_scheme_compressible_euler_init(this, msh, lx, params, user, &
165 chkp)
166 class(fluid_scheme_compressible_euler_t), target, intent(inout) :: this
167 type(mesh_t), target, intent(inout) :: msh
168 integer, intent(in) :: lx
169 type(json_file), target, intent(inout) :: params
170 type(user_t), target, intent(in) :: user
171 type(chkp_t), target, intent(inout) :: chkp
172 character(len=12), parameter :: scheme = 'compressible'
173 integer :: rk_order
174
175 call this%free()
176
177 ! Initialize base class
178 call this%scheme_init(msh, lx, params, scheme, user)
179
180 call euler_rhs_factory(this%euler_rhs)
181
182 associate(xh_lx => this%Xh%lx, xh_ly => this%Xh%ly, xh_lz => this%Xh%lz, &
183 dm_xh => this%dm_Xh, nelv => this%msh%nelv)
184
185 call this%drho%init(dm_xh, 'drho')
186 call this%dm_x%init(dm_xh, 'dm_x')
187 call this%dm_y%init(dm_xh, 'dm_y')
188 call this%dm_z%init(dm_xh, 'dm_z')
189 call this%dE%init(dm_xh, 'dE')
190 call this%h%init(dm_xh, 'h')
191
192 end associate
193
194 if (neko_bcknd_device .eq. 1) then
195 associate(p => this%p, rho => this%rho, &
196 u => this%u, v => this%v, w => this%w, &
197 m_x => this%m_x, m_y => this%m_y, m_z => this%m_z, &
198 effective_visc => this%effective_visc)
199 call device_memcpy(p%x, p%x_d, p%dof%size(), &
200 host_to_device, sync = .false.)
201 call device_memcpy(rho%x, rho%x_d, rho%dof%size(), &
202 host_to_device, sync = .false.)
203 call device_memcpy(u%x, u%x_d, u%dof%size(), &
204 host_to_device, sync = .false.)
205 call device_memcpy(v%x, v%x_d, v%dof%size(), &
206 host_to_device, sync = .false.)
207 call device_memcpy(w%x, w%x_d, w%dof%size(), &
208 host_to_device, sync = .false.)
209 call device_memcpy(m_x%x, m_x%x_d, m_x%dof%size(), &
210 host_to_device, sync = .false.)
211 call device_memcpy(m_y%x, m_y%x_d, m_y%dof%size(), &
212 host_to_device, sync = .false.)
213 call device_memcpy(m_z%x, m_z%x_d, m_z%dof%size(), &
214 host_to_device, sync = .false.)
215 call device_memcpy(effective_visc%x, effective_visc%x_d, &
216 effective_visc%dof%size(), host_to_device, sync = .false.)
217 end associate
218 end if
219
220 ! Initialize the diffusion operator
221 call ax_helm_factory(this%Ax, full_formulation = .false.)
222
223 ! Compute h
224 call this%compute_h()
225
226 ! Initialize regularization
227 call this%setup_regularization(params)
228
229 ! Initialize Runge-Kutta scheme
230 call json_get_or_default(params, 'case.numerics.time_order', rk_order, 4)
231 call this%rk_scheme%init(rk_order)
232
233 ! Set up boundary conditions
234 call this%setup_bcs(user, params)
235
236 end subroutine fluid_scheme_compressible_euler_init
237
240 subroutine fluid_scheme_compressible_euler_free(this)
241 class(fluid_scheme_compressible_euler_t), intent(inout) :: this
242
243 call this%scheme_free()
244
245 if (allocated(this%Ax)) then
246 deallocate(this%Ax)
247 end if
248
249 if (allocated(this%euler_rhs)) then
250 deallocate(this%euler_rhs)
251 end if
252
253 call this%drho%free()
254 call this%dm_x%free()
255 call this%dm_y%free()
256 call this%dm_z%free()
257 call this%dE%free()
258 call this%h%free()
259
260 if (allocated(this%regularization)) then
261 call this%regularization%free()
262 deallocate(this%regularization)
263 end if
264
265 call this%bcs_density%free()
266
267 end subroutine fluid_scheme_compressible_euler_free
268
274 subroutine fluid_scheme_compressible_euler_step(this, time, dt_controller)
276 class(fluid_scheme_compressible_euler_t), target, intent(inout) :: this
277 type(time_state_t), intent(in) :: time
278 type(time_step_controller_t), intent(in) :: dt_controller
279 type(field_t), pointer :: temp
280 integer :: temp_indices(1)
281 ! number of degrees of freedom
282 integer :: n
283 integer :: i
284 class(bc_t), pointer :: b
285
286 n = this%dm_Xh%size()
287 call neko_scratch_registry%request_field(temp, temp_indices(1), .false.)
288 b => null()
289
290 call profiler_start_region('Fluid compressible', 1)
291 associate(u => this%u, v => this%v, w => this%w, p => this%p, &
292 m_x=> this%m_x, m_y => this%m_y, m_z => this%m_z, &
293 xh => this%Xh, msh => this%msh, ax => this%Ax, &
294 c_xh => this%c_Xh, dm_xh => this%dm_Xh, gs_xh => this%gs_Xh, &
295 e => this%E, rho => this%rho, mu => this%mu, &
296 f_x => this%f_x, f_y => this%f_y, f_z => this%f_z, &
297 drho => this%drho, dm_x => this%dm_x, dm_y => this%dm_y, &
298 dm_z => this%dm_z, de => this%dE, &
299 euler_rhs => this%euler_rhs, h => this%h, &
300 t => time%t, tstep => time%tstep, dt => time%dt, &
301 c_avisc_low => this%c_avisc_low, rk_scheme => this%rk_scheme)
302
303 ! Compute artificial viscosity
304 call this%regularization%compute(time, time%tstep, time%dt)
305
306 ! Execute RHS step with effective viscosity field
307 call euler_rhs%step(rho, m_x, m_y, m_z, e, &
308 p, u, v, w, ax, &
309 c_xh, gs_xh, h, this%effective_visc, &
310 rk_scheme, dt)
311
313 call this%bcs_density%apply(rho, time)
314
316 ! Update u, v, w
317 if (neko_bcknd_device .eq. 1) then
318 call compressible_ops_device_update_uvw(u%x_d, v%x_d, w%x_d, &
319 m_x%x_d, m_y%x_d, m_z%x_d, rho%x_d, n)
320 else
321 call compressible_ops_cpu_update_uvw(u%x, v%x, w%x, &
322 m_x%x, m_y%x, m_z%x, rho%x, n)
323 end if
324
326 call this%bcs_vel%apply_vector(u%x, v%x, w%x, &
327 dm_xh%size(), time, strong = .true.)
328
330 if (neko_bcknd_device .eq. 1) then
331 call compressible_ops_device_update_mxyz_p_ruvw(m_x%x_d, m_y%x_d, &
332 m_z%x_d, p%x_d, temp%x_d, u%x_d, v%x_d, w%x_d, e%x_d, &
333 rho%x_d, this%gamma, n)
334 else
335 call compressible_ops_cpu_update_mxyz_p_ruvw(m_x%x, m_y%x, m_z%x, &
336 p%x, temp%x, u%x, v%x, w%x, e%x, rho%x, this%gamma, n)
337 end if
338
340 call this%bcs_prs%apply(p, time)
341
342
344 if (neko_bcknd_device .eq. 1) then
345 call compressible_ops_device_update_e(e%x_d, p%x_d, &
346 temp%x_d, this%gamma, n)
347 else
348 call compressible_ops_cpu_update_e(e%x, p%x, temp%x, this%gamma, n)
349 end if
350
351
353 call this%compute_entropy()
354
356 if (allocated(this%regularization)) then
357 select type (reg => this%regularization)
358 type is (entropy_viscosity_t)
359 call reg%update_lag()
360 end select
361 end if
362
364 call this%compute_max_wave_speed()
365
366 do i = 1, this%bcs_vel%size()
367 b => this%bcs_vel%get(i)
368 b%updated = .false.
369 end do
370
371 do i = 1, this%bcs_prs%size()
372 b => this%bcs_prs%get(i)
373 b%updated = .false.
374 end do
375
376 do i = 1, this%bcs_density%size()
377 b => this%bcs_density%get(i)
378 b%updated = .false.
379 end do
380 nullify(b)
381
382 end associate
383 call profiler_end_region('Fluid compressible', 1)
384
385 call neko_scratch_registry%relinquish_field(temp_indices)
386
387 end subroutine fluid_scheme_compressible_euler_step
388
393 subroutine fluid_scheme_compressible_euler_setup_bcs(this, user, params)
394 class(fluid_scheme_compressible_euler_t), target, intent(inout) :: this
395 type(user_t), target, intent(in) :: user
396 type(json_file), intent(inout) :: params
397 integer :: i, n_bcs, zone_index, j, zone_size, global_zone_size, ierr
398 class(bc_t), pointer :: bc_i
399 type(json_core) :: core
400 type(json_value), pointer :: bc_object
401 type(json_file) :: bc_subdict
402 logical :: found
403 integer, allocatable :: zone_indices(:)
404 character(len=LOG_SIZE) :: log_buf
405
406 ! Process boundary conditions
407 if (params%valid_path('case.fluid.boundary_conditions')) then
408 call params%info('case.fluid.boundary_conditions', n_children = n_bcs)
409 call params%get_core(core)
410 call params%get('case.fluid.boundary_conditions', bc_object, found)
411
412 !
413 ! Velocity bcs
414 !
415 call this%bcs_vel%init(n_bcs)
416
417 do i = 1, n_bcs
418 ! Extract BC configuration
419 call json_extract_item(core, bc_object, i, bc_subdict)
420 call json_get(bc_subdict, "zone_indices", zone_indices)
421
422 ! Validate zones
423 do j = 1, size(zone_indices)
424 zone_size = this%msh%labeled_zones(zone_indices(j))%size
425 call mpi_allreduce(zone_size, global_zone_size, 1, &
426 mpi_integer, mpi_max, neko_comm, ierr)
427
428 if (global_zone_size .eq. 0) then
429 write(log_buf, '(A,I0,A)') "Error: Zone ", zone_indices(j), &
430 " has zero size"
431 call neko_error(log_buf)
432 end if
433 end do
434
435 ! Create BC
436 bc_i => null()
437 call velocity_bc_factory(bc_i, this, bc_subdict, this%c_Xh, user)
438
439 ! Add to appropriate lists
440 if (associated(bc_i)) then
441 call this%bcs_vel%append(bc_i)
442 end if
443 end do
444
445 !
446 ! Pressure bcs
447 !
448 call this%bcs_prs%init(n_bcs)
449
450 do i = 1, n_bcs
451 ! Create a new json containing just the subdict for this bc
452 call json_extract_item(core, bc_object, i, bc_subdict)
453 bc_i => null()
454 call pressure_bc_factory(bc_i, this, bc_subdict, this%c_Xh, user)
455
456 ! Not all bcs require an allocation for pressure in particular,
457 ! so we check.
458 if (associated(bc_i)) then
459 call this%bcs_prs%append(bc_i)
460 end if
461 end do
462
463 !
464 ! Density bcs
465 !
466 call this%bcs_density%init(n_bcs)
467
468 do i = 1, n_bcs
469 ! Create a new json containing just the subdict for this bc
470 call json_extract_item(core, bc_object, i, bc_subdict)
471 bc_i => null()
472 call density_bc_factory(bc_i, this, bc_subdict, this%c_Xh, user)
473
474 ! Not all bcs require an allocation for pressure in particular,
475 ! so we check.
476 if (associated(bc_i)) then
477 call this%bcs_density%append(bc_i)
478 end if
479 end do
480 else
481 ! Check that there are no labeled zones, i.e. all are periodic.
482 do i = 1, size(this%msh%labeled_zones)
483 if (this%msh%labeled_zones(i)%size .gt. 0) then
484 call neko_error("No boundary_conditions entry in the case file!")
485 end if
486 end do
487
488 ! For a pure periodic case, we still need to initilise the bc lists
489 ! to a zero size to avoid issues with apply() in step()
490 call this%bcs_prs%init()
491 call this%bcs_vel%init()
492 call this%bcs_density%init()
493
494 end if
495 end subroutine fluid_scheme_compressible_euler_setup_bcs
496
501 subroutine compute_h(this)
502 class(fluid_scheme_compressible_euler_t), intent(inout) :: this
503 integer :: e, i, j, k
504 integer :: im, ip, jm, jp, km, kp
505 real(kind=rp) :: di, dj, dk, ndim_inv
506 integer :: lx_half, ly_half, lz_half
507
508 lx_half = this%c_Xh%Xh%lx / 2
509 ly_half = this%c_Xh%Xh%ly / 2
510 lz_half = this%c_Xh%Xh%lz / 2
511
512 do concurrent(e = 1:this%c_Xh%msh%nelv)
513 do concurrent(k = 1:this%c_Xh%Xh%lz, &
514 j = 1:this%c_Xh%Xh%ly, i = 1:this%c_Xh%Xh%lx)
515 km = max(1, k-1)
516 kp = min(this%c_Xh%Xh%lz, k+1)
517
518 jm = max(1, j-1)
519 jp = min(this%c_Xh%Xh%ly, j+1)
520
521 im = max(1, i-1)
522 ip = min(this%c_Xh%Xh%lx, i+1)
523
524 di = (this%c_Xh%dof%x(ip, j, k, e) - &
525 this%c_Xh%dof%x(im, j, k, e))**2 &
526 + (this%c_Xh%dof%y(ip, j, k, e) - &
527 this%c_Xh%dof%y(im, j, k, e))**2 &
528 + (this%c_Xh%dof%z(ip, j, k, e) - &
529 this%c_Xh%dof%z(im, j, k, e))**2
530
531 dj = (this%c_Xh%dof%x(i, jp, k, e) - &
532 this%c_Xh%dof%x(i, jm, k, e))**2 &
533 + (this%c_Xh%dof%y(i, jp, k, e) - &
534 this%c_Xh%dof%y(i, jm, k, e))**2 &
535 + (this%c_Xh%dof%z(i, jp, k, e) - &
536 this%c_Xh%dof%z(i, jm, k, e))**2
537
538 dk = (this%c_Xh%dof%x(i, j, kp, e) - &
539 this%c_Xh%dof%x(i, j, km, e))**2 &
540 + (this%c_Xh%dof%y(i, j, kp, e) - &
541 this%c_Xh%dof%y(i, j, km, e))**2 &
542 + (this%c_Xh%dof%z(i, j, kp, e) - &
543 this%c_Xh%dof%z(i, j, km, e))**2
544
545 di = sqrt(di) / (ip - im)
546 dj = sqrt(dj) / (jp - jm)
547 dk = sqrt(dk) / (kp - km)
548 this%h%x(i,j,k,e) = (di * dj * dk)**(1.0_rp / 3.0_rp)
549
550 end do
551 end do
552
553 if (neko_bcknd_device .eq. 1) then
554 call device_memcpy(this%h%x, this%h%x_d, this%h%dof%size(),&
555 host_to_device, sync = .false.)
556 call this%gs_Xh%op(this%h, gs_op_add)
557 call device_col2(this%h%x_d, this%c_Xh%mult_d, this%h%dof%size())
558 else
559 call this%gs_Xh%op(this%h, gs_op_add)
560 call col2(this%h%x, this%c_Xh%mult, this%h%dof%size())
561 end if
562
563 end subroutine compute_h
564
569 subroutine fluid_scheme_compressible_euler_restart(this, chkp)
570 class(fluid_scheme_compressible_euler_t), target, intent(inout) :: this
571 type(chkp_t), intent(inout) :: chkp
572 end subroutine fluid_scheme_compressible_euler_restart
573
574 subroutine setup_regularization(this, params)
577 class(fluid_scheme_compressible_euler_t), target, intent(inout) :: this
578 type(json_file), intent(inout) :: params
579 type(json_file) :: reg_json
580 type(json_core) :: json_core_inst
581 type(json_value), pointer :: reg_params
582 character(len=:), allocatable :: buffer
583 real(kind=rp) :: c_avisc_entropy_val
584 character(len=:), allocatable :: regularization_type
585
586 call json_get_or_default(params, 'case.numerics.c_avisc_low', &
587 this%c_avisc_low, 0.5_rp)
588 call json_get_or_default(params, 'case.numerics.c_avisc_entropy', &
589 c_avisc_entropy_val, 1.0_rp)
590
591 call json_core_inst%initialize()
592 call json_core_inst%create_object(reg_params, '')
593 call json_core_inst%add(reg_params, 'c_avisc_entropy', c_avisc_entropy_val)
594 call json_core_inst%add(reg_params, 'c_avisc_low', this%c_avisc_low)
595 call json_core_inst%print_to_string(reg_params, buffer)
596 call json_core_inst%destroy(reg_params)
597
598 call reg_json%initialize()
599 call reg_json%load_from_string(buffer)
600
601 regularization_type = 'entropy_viscosity'
602
603 call regularization_factory(this%regularization, regularization_type, &
604 reg_json, this%c_Xh, this%dm_Xh, this%effective_visc)
605
606 select type (reg => this%regularization)
607 type is (entropy_viscosity_t)
608 call entropy_viscosity_set_fields(reg, this%S, this%u, this%v, this%w, &
609 this%h, this%max_wave_speed, this%msh, this%Xh, this%gs_Xh)
610 end select
611
612 call reg_json%destroy()
613
614 end subroutine setup_regularization
615
Copy data between host and device (or device and device)
Definition device.F90:71
Abstract interface to evaluate rhs.
Definition euler_res.f90:51
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Retrieves a parameter by name or throws an error.
Subroutines to add advection terms to the RHS of a transport equation.
Definition advection.f90:34
Defines a Matrix-vector product.
Definition ax.f90:34
Defines a list of bc_t.
Definition bc_list.f90:34
Defines a boundary condition.
Definition bc.f90:34
Backward-differencing scheme for time integration.
Generic buffer that is extended with buffers of varying rank.
Definition buffer.F90:34
Defines a checkpoint.
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
CPU implementation of compressible flow operations.
subroutine, public compressible_ops_cpu_update_uvw(u, v, w, m_x, m_y, m_z, rho, n)
Update u,v,w fields.
subroutine, public compressible_ops_cpu_update_e(e, p, ruvw, gamma, n)
Update E field.
subroutine, public compressible_ops_cpu_update_mxyz_p_ruvw(m_x, m_y, m_z, p, ruvw, u, v, w, e, rho, gamma, n)
Update m_x, m_y, m_z, p, ruvw, fields.
Device implementation of compressible flow operations.
subroutine, public compressible_ops_device_update_uvw(u_d, v_d, w_d, m_x_d, m_y_d, m_z_d, rho_d, n)
Update u,v,w fields.
subroutine, public compressible_ops_device_update_mxyz_p_ruvw(m_x_d, m_y_d, m_z_d, p_d, ruvw_d, u_d, v_d, w_d, e_d, rho_d, gamma, n)
Update m_x, m_y, m_z, p, ruvw, fields.
subroutine, public compressible_ops_device_update_e(e_d, p_d, ruvw_d, gamma, n)
Update E field.
subroutine, public device_col2(a_d, b_d, n, strm)
Vector multiplication .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:47
subroutine, public entropy_viscosity_set_fields(this, s, u, v, w, h, max_wave_speed, msh, xh, gs)
Contains the field_serties_t type.
Defines a field.
Definition field.f90:34
subroutine compute_h(this)
Copied from les_model_compute_delta in les_model.f90 TODO: move to a separate module Compute characte...
subroutine fluid_scheme_compressible_euler_init(this, msh, lx, params, user, chkp)
Boundary condition factory for density.
subroutine fluid_scheme_compressible_euler_step(this, time, dt_controller)
Advance the fluid simulation one timestep.
subroutine fluid_scheme_compressible_euler_free(this)
Free allocated memory and cleanup.
subroutine fluid_scheme_compressible_euler_restart(this, chkp)
Restart the simulation from saved state.
subroutine fluid_scheme_compressible_euler_setup_bcs(this, user, params)
Set up boundary conditions for the fluid scheme.
Gather-scatter.
Defines Gather-scatter operations.
Definition gs_ops.f90:34
integer, parameter, public gs_op_add
Definition gs_ops.f90:36
integer, parameter, public gs_op_max
Definition gs_ops.f90:36
integer, parameter, public gs_op_min
Definition gs_ops.f90:36
Utilities for retrieving parameters from the case files.
Logging routines.
Definition log.f90:34
integer, parameter, public log_size
Definition log.f90:46
Definition math.f90:60
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:854
Defines a mesh.
Definition mesh.f90:34
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators.
Definition operators.f90:34
subroutine, public div(res, ux, uy, uz, coef)
Compute the divergence of a vector field.
Profiling interface.
Definition profiler.F90:34
subroutine, public profiler_start_region(name, region_id)
Started a named (name) profiler region.
Definition profiler.F90:78
subroutine, public profiler_end_region(name, region_id)
End the most recently started profiler region.
Definition profiler.F90:107
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.
Compound scheme for the advection and diffusion operators in a transport equation.
Module with things related to the simulation time.
Implements type time_step_controller.
Interfaces for user interaction with NEKO.
Definition user_intf.f90:34
Utilities.
Definition utils.f90:35
subroutine, public neko_type_error(base_type, wrong_type, known_types)
Reports an error allocating a type for a particular base pointer class.
Definition utils.f90:313
Base abstract type for computing the advection operator.
Definition advection.f90:46
Base type for a matrix-vector product providing .
Definition ax.f90:43
Base type for a boundary condition.
Definition bc.f90:62
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:48
Implicit backward-differencing scheme for time integration.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:56
Abstract type to compute rhs.
Definition euler_res.f90:44
Stores a series (sequence) of fields, logically connected to a base field, and arranged according to ...
Implements the logic to compute the time coefficients for the advection and diffusion operators in a ...
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...
#define max(a, b)
Definition tensor.cu:40