52 use precon,
only :
pc_t, precon_factory, precon_destroy
63 use json_module,
only : json_file
79 class(
ksp_t),
allocatable :: ksp_vel
80 class(
ksp_t),
allocatable :: ksp_prs
81 class(
pc_t),
allocatable :: pc_vel
82 class(
pc_t),
allocatable :: pc_prs
83 integer :: vel_projection_dim
84 integer :: pr_projection_dim
85 integer :: vel_projection_activ_step
86 integer :: pr_projection_activ_step
87 logical :: pr_projection_reorthogonalize_basis
88 logical :: strict_convergence
89 logical :: allow_stabilization
96 logical :: forced_flow_rate = .false.
99 character(len=:),
allocatable :: nut_field_name
105 integer(kind=i8) :: glb_n_points
107 integer(kind=i8) :: glb_unique_points
121 procedure, pass(this) :: set_material_properties => &
135 module subroutine fluid_scheme_factory(object, type_name)
137 character(len=*) :: type_name
138 end subroutine fluid_scheme_factory
141 public :: fluid_scheme_incompressible_t, fluid_scheme_factory
146 subroutine fluid_scheme_init_base(this, msh, lx, params, scheme, user, &
149 class(fluid_scheme_incompressible_t),
target,
intent(inout) :: this
150 type(
mesh_t),
target,
intent(inout) :: msh
151 integer,
intent(in) :: lx
152 character(len=*),
intent(in) :: scheme
153 type(json_file),
target,
intent(inout) :: params
154 type(
user_t),
target,
intent(in) :: user
155 logical,
intent(in) :: kspv_init
157 character(len=LOG_SIZE) :: log_buf
158 real(kind=
rp),
allocatable :: real_vec(:)
159 real(kind=
rp) :: real_val, kappa, b, z0
160 logical :: logical_val
161 integer :: integer_val, ierr
162 type(json_file) :: wm_json
163 character(len=:),
allocatable :: string_val1, string_val2
164 type(json_file) :: json_subdict
172 if (msh%gdim .eq. 2)
then
173 call this%Xh%init(
gll, lx, lx)
175 call this%Xh%init(
gll, lx, lx, lx)
178 call this%dm_Xh%init(msh, this%Xh)
180 call this%gs_Xh%init(this%dm_Xh)
182 call this%c_Xh%init(this%gs_Xh)
195 write(log_buf,
'(A, A)')
'Type : ', trim(scheme)
197 write(log_buf,
'(A, A)')
'Name : ', trim(this%name)
211 call this%set_material_properties(params,
user)
215 'case.fluid.velocity_solver.projection_space_size', &
216 this%vel_projection_dim, 0)
218 'case.fluid.pressure_solver.projection_space_size', &
219 this%pr_projection_dim, 0)
221 'case.fluid.velocity_solver.projection_hold_steps', &
222 this%vel_projection_activ_step, 5)
224 'case.fluid.pressure_solver.projection_hold_steps', &
225 this%pr_projection_activ_step, 5)
227 'case.fluid.pressure_solver.projection_reorthogonalize_basis', &
228 this%pr_projection_reorthogonalize_basis, .false.)
232 if (params%valid_path(
"case.fluid.flow_rate_force"))
then
233 this%forced_flow_rate = .true.
238 write(log_buf,
'(A, I1)')
'Poly order : ', lx-1
239 else if (lx .ge. 10)
then
240 write(log_buf,
'(A, I2)')
'Poly order : ', lx-1
242 write(log_buf,
'(A, I3)')
'Poly order : ', lx-1
245 this%glb_n_points = int(this%msh%glb_nelv,
i8)*int(this%Xh%lxyz,
i8)
246 this%glb_unique_points = int(
glsum(this%c_Xh%mult, this%dm_Xh%size()),
i8)
248 write(log_buf,
'(A, I0)')
'GLL points : ', this%glb_n_points
250 write(log_buf,
'(A, I0)')
'Unique pts.: ', this%glb_unique_points
254 call json_get(params,
'case.numerics.dealias', logical_val)
255 write(log_buf,
'(A, L1)')
'Dealias : ', logical_val
261 write(log_buf,
'(A, L1)')
'Save bdry : ', logical_val
265 logical_val, .false.)
266 write(log_buf,
'(A, L1)')
'Full stress: ', logical_val
276 call this%f_x%init(this%dm_Xh, fld_name =
"fluid_rhs_x")
277 call this%f_y%init(this%dm_Xh, fld_name =
"fluid_rhs_y")
278 call this%f_z%init(this%dm_Xh, fld_name =
"fluid_rhs_z")
282 call neko_log%section(
"Velocity solver")
284 'case.fluid.velocity_solver.max_iterations', &
286 call json_get(params,
'case.fluid.velocity_solver.type', string_val1)
287 call json_get(params,
'case.fluid.velocity_solver.preconditioner.type', &
290 'case.fluid.velocity_solver.preconditioner', json_subdict)
292 'case.fluid.velocity_solver.absolute_tolerance', &
295 'case.fluid.velocity_solver.monitor', &
296 logical_val, .false.)
298 call neko_log%message(
'Type : ('// trim(string_val1) // &
299 ', ' // trim(string_val2) //
')')
301 write(log_buf,
'(A,ES13.6)')
'Abs tol :', real_val
303 call this%solver_factory(this%ksp_vel, this%dm_Xh%size(), &
304 string_val1, integer_val, real_val, logical_val)
305 call this%precon_factory_(this%pc_vel, this%ksp_vel, &
306 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs_vel, &
307 string_val2, json_subdict)
313 this%strict_convergence, .false.)
316 this%allow_stabilization, .false.)
320 call this%ulag%init(this%u, 2)
321 call this%vlag%init(this%v, 2)
322 call this%wlag%init(this%w, 2)
332 call neko_log%section(
'Fluid Source term')
333 call this%source_term%init(this%f_x, this%f_y, this%f_z, this%c_Xh,
user, &
335 call this%source_term%add(params,
'case.fluid.source_terms')
338 end subroutine fluid_scheme_init_base
340 subroutine fluid_scheme_free(this)
341 class(fluid_scheme_incompressible_t),
intent(inout) :: this
342 class(
bc_t),
pointer :: bc
348 if (
allocated(this%ksp_vel))
then
349 call this%ksp_vel%free()
350 deallocate(this%ksp_vel)
353 if (
allocated(this%ksp_prs))
then
354 call this%ksp_prs%free()
355 deallocate(this%ksp_prs)
358 if (
allocated(this%pc_vel))
then
359 call precon_destroy(this%pc_vel)
360 deallocate(this%pc_vel)
363 if (
allocated(this%pc_prs))
then
364 call precon_destroy(this%pc_prs)
365 deallocate(this%pc_prs)
368 do i = 1, this%bcs_vel%size()
369 bc => this%bcs_vel%get(i)
372 call this%bcs_vel%free()
374 do i = 1, this%bcs_prs%size()
375 bc => this%bcs_prs%get(i)
378 call this%bcs_prs%free()
380 call this%source_term%free()
382 call this%gs_Xh%free()
384 call this%c_Xh%free()
395 call this%ulag%free()
396 call this%vlag%free()
397 call this%wlag%free()
400 if (
associated(this%f_x))
then
405 if (
associated(this%f_y))
then
410 if (
associated(this%f_z))
then
422 call this%dm_Xh%free()
426 end subroutine fluid_scheme_free
430 subroutine fluid_scheme_validate(this)
431 class(fluid_scheme_incompressible_t),
target,
intent(inout) :: this
433 logical :: logical_val
435 if ( (.not.
associated(this%u)) .or. &
436 (.not.
associated(this%v)) .or. &
437 (.not.
associated(this%w)) .or. &
438 (.not.
associated(this%p)))
then
442 if ( (.not.
allocated(this%u%x)) .or. &
443 (.not.
allocated(this%v%x)) .or. &
444 (.not.
allocated(this%w%x)) .or. &
445 (.not.
allocated(this%p%x)))
then
449 if (.not.
allocated(this%ksp_vel))
then
450 call neko_error(
'No Krylov solver for velocity defined')
453 if (.not.
allocated(this%ksp_prs))
then
454 call neko_error(
'No Krylov solver for pressure defined')
457 end subroutine fluid_scheme_validate
463 subroutine fluid_scheme_bc_apply_vel(this, time, strong)
464 class(fluid_scheme_incompressible_t),
intent(inout) :: this
466 logical,
intent(in) :: strong
468 class(
bc_t),
pointer :: b
471 call this%bcs_vel%apply_vector(&
472 this%u%x, this%v%x, this%w%x, this%dm_Xh%size(), time, strong)
474 call rotate_cyc(this%u%x, this%v%x, this%w%x, 1, this%c_Xh)
481 call rotate_cyc(this%u%x, this%v%x, this%w%x, 0, this%c_Xh)
484 call this%bcs_vel%apply_vector(&
485 this%u%x, this%v%x, this%w%x, this%dm_Xh%size(), time, strong)
487 call rotate_cyc(this%u%x, this%v%x, this%w%x, 1, this%c_Xh)
494 call rotate_cyc(this%u%x, this%v%x, this%w%x, 0, this%c_Xh)
496 do i = 1, this%bcs_vel%size()
497 b => this%bcs_vel%get(i)
502 end subroutine fluid_scheme_bc_apply_vel
506 subroutine fluid_scheme_bc_apply_prs(this, time)
507 class(fluid_scheme_incompressible_t),
intent(inout) :: this
511 class(
bc_t),
pointer :: b
514 call this%bcs_prs%apply(this%p, time)
518 call this%bcs_prs%apply(this%p, time)
522 do i = 1, this%bcs_prs%size()
523 b => this%bcs_prs%get(i)
528 end subroutine fluid_scheme_bc_apply_prs
532 subroutine fluid_scheme_solver_factory(ksp, n, solver, &
533 max_iter, abstol, monitor)
534 class(
ksp_t),
allocatable,
target,
intent(inout) :: ksp
535 integer,
intent(in),
value :: n
536 character(len=*),
intent(in) :: solver
537 integer,
intent(in) :: max_iter
538 real(kind=
rp),
intent(in) :: abstol
539 logical,
intent(in) :: monitor
541 call krylov_solver_factory(ksp, n, solver, max_iter, abstol, &
544 end subroutine fluid_scheme_solver_factory
547 subroutine fluid_scheme_precon_factory(this, pc, ksp, coef, dof, gs, bclst, &
549 class(fluid_scheme_incompressible_t),
intent(inout) :: this
550 class(
pc_t),
allocatable,
target,
intent(inout) :: pc
551 class(
ksp_t),
target,
intent(inout) :: ksp
552 type(
coef_t),
target,
intent(in) :: coef
553 type(
dofmap_t),
target,
intent(in) :: dof
554 type(
gs_t),
target,
intent(inout) :: gs
555 type(
bc_list_t),
target,
intent(inout) :: bclst
556 character(len=*) :: pctype
557 type(json_file),
intent(inout) :: pcparams
559 call precon_factory(pc, pctype)
561 select type (pcp => pc)
563 call pcp%init(coef, dof, gs)
565 call pcp%init(coef, dof, gs)
567 call pcp%init(coef, dof, gs)
569 call pcp%init(coef, bclst, pcparams)
571 call pcp%init(coef, bclst, pcparams)
576 end subroutine fluid_scheme_precon_factory
579 function fluid_compute_cfl(this, dt)
result(c)
580 class(fluid_scheme_incompressible_t),
intent(in) :: this
581 real(kind=
rp),
intent(in) :: dt
584 c =
cfl(dt, this%u%x, this%v%x, this%w%x, &
585 this%Xh, this%c_Xh, this%msh%nelv, this%msh%gdim)
587 end function fluid_compute_cfl
594 subroutine fluid_scheme_update_material_properties(this, time)
595 class(fluid_scheme_incompressible_t),
intent(inout) :: this
599 call this%user_material_properties(this%name, this%material_properties, &
602 if (len(trim(this%nut_field_name)) > 0)
then
615 call device_memcpy(this%rho%x, this%rho%x_d, this%rho%size(), &
618 end subroutine fluid_scheme_update_material_properties
623 subroutine fluid_scheme_set_material_properties(this, params, user)
624 class(fluid_scheme_incompressible_t),
target,
intent(inout) :: this
625 type(json_file),
intent(inout) :: params
626 type(
user_t),
target,
intent(in) :: user
627 character(len=LOG_SIZE) :: log_buf
630 logical :: nondimensional
631 real(kind=
rp) :: dummy_lambda, dummy_cp
632 real(kind=
rp) :: const_mu, const_rho
639 call neko_registry%add_field(this%dm_Xh, this%name //
"_mu_tot")
640 call neko_registry%add_field(this%dm_Xh, this%name //
"_rho")
642 this%mu_tot =>
neko_registry%get_field(this%name //
"_mu_tot")
645 call this%material_properties%init(2)
646 call this%material_properties%assign(1, this%rho)
647 call this%material_properties%assign(2, this%mu)
649 if (.not.
associated(
user%material_properties, dummy_mp_ptr))
then
651 write(log_buf,
'(A)')
'Material properties must be set in the user' // &
654 this%user_material_properties =>
user%material_properties
656 call user%material_properties(this%name, this%material_properties, &
662 if (params%valid_path(
'case.fluid.Re') .and. &
663 (params%valid_path(
'case.fluid.mu') .or. &
664 params%valid_path(
'case.fluid.rho')))
then
665 call neko_error(
"To set the material properties for the fluid, " // &
666 "either provide Re OR mu and rho in the case file.")
668 else if (params%valid_path(
'case.fluid.Re'))
then
670 write(log_buf,
'(A)')
'Non-dimensional fluid material properties &
673 write(log_buf,
'(A)')
'Density will be set to 1, dynamic viscosity to&
679 write(log_buf,
'(A)')
'Read non-dimensional material properties'
681 write(log_buf,
'(A,ES13.6)')
'Re :', const_mu
687 const_mu = 1.0_rp/const_mu
697 if (
associated(
user%material_properties, dummy_mp_ptr))
then
704 write(log_buf,
'(A,ES13.6)')
'rho :', const_rho
706 write(log_buf,
'(A,ES13.6)')
'mu :', const_mu
718 call device_memcpy(this%rho%x, this%rho%x_d, this%rho%size(), &
722 call device_memcpy(this%mu_tot%x, this%mu_tot%x_d, this%mu%size(), &
725 end subroutine fluid_scheme_set_material_properties
Copy data between host and device (or device and device)
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...
Retrieves a parameter by name or throws an error.
Abstract interface for setting material properties.
Defines a boundary condition.
Jacobi preconditioner accelerator backend.
Device abstraction, common interface for various accelerators.
subroutine, public device_event_sync(event)
Synchronize an event.
integer, parameter, public device_to_host
type(c_ptr), bind(C), public glb_cmd_event
Event for the global command queue.
Defines a dirichlet boundary condition.
Defines a mapping of the degrees of freedom.
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_copy(a, b, n)
Copy a vector .
subroutine fluid_scheme_set_material_properties(this, params, user)
Sets rho and mu.
subroutine fluid_scheme_update_material_properties(this, time)
Call user material properties routine and update the values of mu if necessary.
subroutine fluid_scheme_precon_factory(this, pc, ksp, coef, dof, gs, bclst, pctype, pcparams)
Initialize a Krylov preconditioner.
subroutine fluid_scheme_free(this)
subroutine fluid_scheme_validate(this)
Validate that all fields, solvers etc necessary for performing time-stepping are defined.
real(kind=rp) function fluid_compute_cfl(this, dt)
Compute CFL.
subroutine fluid_scheme_bc_apply_vel(this, time, strong)
Apply all boundary conditions defined for velocity Here we perform additional gs operations to take c...
subroutine fluid_scheme_solver_factory(ksp, n, solver, max_iter, abstol, monitor)
Initialize a linear solver.
subroutine fluid_scheme_bc_apply_prs(this, time)
Apply all boundary conditions defined for pressure.
subroutine fluid_scheme_init_base(this, msh, lx, params, scheme, user, kspv_init)
Initialise a fluid scheme.
Implements the fluid_source_term_t type.
Computes various statistics for the fluid fields. We use the Reynolds decomposition for a field u = ...
Utilities for retrieving parameters from the case files.
Implements the base abstract type for Krylov solvers plus helper types.
integer, parameter, public ksp_max_iter
Maximum number of iters.
integer, parameter, public neko_log_verbose
Verbose log level.
type(log_t), public neko_log
Global log stream.
integer, parameter, public log_size
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
integer, parameter neko_bcknd_device
integer, parameter, public i8
integer, parameter, public rp
Global precision used in computations.
real(kind=rp) function, public cfl(dt, u, v, w, xh, coef, nelv, gdim)
Hybrid ph-multigrid preconditioner.
Defines a registry for storing solution fields.
type(registry_t), target, public neko_registry
Global field registry.
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.
Implements the source_term_t type and a wrapper source_term_wrapper_t.
Defines a function space.
integer, parameter, public gll
Defines a container for all statistics.
Jacobi preconditioner SX-Aurora backend.
Module with things related to the simulation time.
Interfaces for user interaction with NEKO.
subroutine, public dummy_user_material_properties(scheme_name, properties, time)
Base type for a boundary condition.
A list of allocatable `bc_t`. Follows the standard interface of lists.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Defines a jacobi preconditioner.
Generic Dirichlet boundary condition on .
Base type of all fluid formulations.
Base type of all fluid formulations.
Wrapper contaning and executing the fluid source terms.
Defines a jacobi preconditioner.
Base abstract type for a canonical Krylov method, solving .
Defines a canonical Krylov preconditioner.
Defines a jacobi preconditioner for SX-Aurora.
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...