51 use precon,
only :
pc_t, precon_factory, precon_destroy
58 use json_module,
only : json_file
76 character(len=:),
allocatable :: name
92 type(
gs_t),
pointer :: gs_xh
102 integer :: ksp_maxiter
104 integer :: projection_dim
106 integer :: projection_activ_step
112 type(json_file),
pointer :: params
118 character(len=:),
allocatable :: nut_field_name
120 character(len=:),
allocatable :: alphat_field_name
128 type(
field_t),
pointer :: lambda_tot => null()
130 real(kind=
rp) :: pr_turb
134 user_material_properties => null()
136 logical :: freeze = .false.
145 procedure, pass(this) :: set_material_properties => &
148 procedure, pass(this) :: update_material_properties => &
163 numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
175 type(
mesh_t),
target,
intent(in) :: msh
176 type(
coef_t),
target,
intent(in) :: coef
177 type(
gs_t),
target,
intent(inout) :: gs
178 type(json_file),
target,
intent(inout) :: params
179 type(json_file),
target,
intent(inout) :: numerics_params
180 type(
user_t),
target,
intent(in) :: user
181 type(
chkp_t),
target,
intent(inout) :: chkp
184 type(
field_t),
target,
intent(in) :: rho
195 type(
chkp_t),
intent(inout) :: chkp
238 procedure, pass(this) :: move_from => &
241 procedure, pass(this) :: is_allocated => &
260 module subroutine scalar_scheme_factory(object, msh, coef, gs, params, &
262 class(scalar_scheme_t),
allocatable,
intent(inout) :: object
263 type(
mesh_t),
target,
intent(in) :: msh
264 type(
coef_t),
target,
intent(in) :: coef
265 type(
gs_t),
target,
intent(inout) :: gs
266 type(json_file),
target,
intent(inout) :: params
267 type(json_file),
target,
intent(inout) :: numerics_params
269 type(
chkp_t),
target,
intent(inout) :: chkp
272 type(
field_t),
target,
intent(in) :: rho
273 end subroutine scalar_scheme_factory
280 module subroutine scalar_scheme_allocator(object, type_name)
281 class(scalar_scheme_t),
allocatable,
intent(inout) :: object
282 character(len=*),
intent(in):: type_name
283 end subroutine scalar_scheme_allocator
294 subroutine scalar_scheme_allocate(obj)
295 import scalar_scheme_t
296 class(scalar_scheme_t),
allocatable,
intent(inout) :: obj
297 end subroutine scalar_scheme_allocate
302 module subroutine register_scalar_scheme(type_name, allocator)
303 character(len=*),
intent(in) :: type_name
304 procedure(scalar_scheme_allocate),
pointer,
intent(in) :: allocator
305 end subroutine register_scalar_scheme
310 type scalar_scheme_allocator_entry
311 character(len=20) :: type_name
312 procedure(scalar_scheme_allocate),
pointer,
nopass :: allocator
313 end type scalar_scheme_allocator_entry
316 type(scalar_scheme_allocator_entry),
allocatable,
private :: &
317 scalar_scheme_registry(:)
320 integer,
private :: scalar_scheme_registry_size = 0
332 subroutine scalar_scheme_init(this, msh, c_Xh, gs_Xh, params, scheme, user, &
334 class(scalar_scheme_t),
target,
intent(inout) :: this
335 type(
mesh_t),
target,
intent(in) :: msh
336 type(
coef_t),
target,
intent(in) :: c_Xh
337 type(
gs_t),
target,
intent(inout) :: gs_Xh
338 type(json_file),
target,
intent(inout) :: params
339 character(len=*),
intent(in) :: scheme
340 type(
user_t),
target,
intent(in) :: user
341 type(
field_t),
target,
intent(in) :: rho
343 character(len=LOG_SIZE) :: log_buf
345 logical :: logical_val
346 real(kind=
rp) :: real_val, solver_abstol
347 integer :: integer_val, ierr
348 character(len=:),
allocatable :: solver_type, solver_precon
349 type(json_file) :: precon_params
350 type(json_file) :: json_subdict
351 logical :: nut_dependency
361 call json_get(params,
'name', this%name)
367 call json_get(params,
'solver.type', solver_type)
368 call json_get(params,
'solver.preconditioner.type', &
370 call json_get(params,
'solver.preconditioner', precon_params)
375 'solver.projection_space_size', &
376 this%projection_dim, 0)
378 'solver.projection_hold_steps', &
379 this%projection_activ_step, 5)
382 write(log_buf,
'(A, A)')
'Type : ', trim(scheme)
384 write(log_buf,
'(A, A)')
'Name : ', trim(this%name)
386 call neko_log%message(
'Ksp scalar : ('// trim(solver_type) // &
387 ', ' // trim(solver_precon) //
')')
388 write(log_buf,
'(A,ES13.6)')
' `-abs tol :', solver_abstol
392 this%dm_Xh => this%u%dof
393 this%params => params
397 ignore_existing = .true.)
401 call this%slag%init(this%s, 2)
409 call this%set_material_properties(params,
user)
415 this%alphat_field_name =
""
416 this%nut_field_name =
""
417 if (params%valid_path(
'alphat'))
then
418 call json_get(this%params,
'alphat', json_subdict)
419 call json_get(json_subdict,
'nut_dependency', nut_dependency)
420 if (nut_dependency)
then
422 call json_get(json_subdict,
'nut_field', this%nut_field_name)
424 call json_get(json_subdict,
'alphat_field', this%alphat_field_name)
432 call this%f_Xh%init(this%dm_Xh, fld_name =
"scalar_rhs")
435 call this%source_term%init(this%f_Xh, this%c_Xh,
user, this%name)
436 call this%source_term%add(params,
'source_terms')
440 'solver.max_iterations', &
444 logical_val, .false.)
445 call scalar_scheme_solver_factory(this%ksp, this%dm_Xh%size(), &
446 solver_type, integer_val, solver_abstol, logical_val)
447 call scalar_scheme_precon_factory(this%pc, this%ksp, &
448 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs, &
449 solver_precon, precon_params)
453 end subroutine scalar_scheme_init
457 subroutine scalar_scheme_free(this)
458 class(scalar_scheme_t),
intent(inout) :: this
459 class(
bc_t),
pointer :: bc
470 if (
allocated(this%ksp))
then
475 if (
allocated(this%pc))
then
476 call precon_destroy(this%pc)
480 if (
allocated(this%name))
then
481 deallocate(this%name)
484 call this%source_term%free()
486 if (
associated(this%f_Xh))
then
487 call this%f_Xh%free()
488 deallocate(this%f_Xh)
491 do i = 1, this%bcs%size()
492 bc => this%bcs%get(i)
493 if (
associated(
bc))
then
500 call this%slag%free()
501 call this%material_properties%free()
503 if (
allocated(this%nut_field_name))
then
504 deallocate(this%nut_field_name)
507 if (
allocated(this%alphat_field_name))
then
508 deallocate(this%alphat_field_name)
513 nullify(this%lambda_tot)
516 end subroutine scalar_scheme_free
520 subroutine scalar_scheme_validate(this)
521 class(scalar_scheme_t),
target,
intent(inout) :: this
523 if ( (.not.
allocated(this%u%x)) .or. &
524 (.not.
allocated(this%v%x)) .or. &
525 (.not.
allocated(this%w%x)) .or. &
526 (.not.
allocated(this%s%x)))
then
530 if (.not.
allocated(this%ksp))
then
531 call neko_error(
'No Krylov solver for velocity defined')
534 if (.not.
associated(this%Xh))
then
538 if (.not.
associated(this%dm_Xh))
then
542 if (.not.
associated(this%c_Xh))
then
546 if (.not.
associated(this%f_Xh))
then
550 if (.not.
associated(this%params))
then
554 if (.not.
associated(this%rho))
then
558 end subroutine scalar_scheme_validate
562 subroutine scalar_scheme_solver_factory(ksp, n, solver, max_iter, &
564 class(
ksp_t),
allocatable,
target,
intent(inout) :: ksp
565 integer,
intent(in),
value :: n
566 integer,
intent(in) :: max_iter
567 character(len=*),
intent(in) :: solver
568 real(kind=
rp) :: abstol
569 logical,
intent(in) :: monitor
571 call krylov_solver_factory(ksp, n, solver, max_iter, &
572 abstol, monitor = monitor)
574 end subroutine scalar_scheme_solver_factory
577 subroutine scalar_scheme_precon_factory(pc, ksp, coef, dof, gs, bclst, &
579 class(
pc_t),
allocatable,
target,
intent(inout) :: pc
580 class(
ksp_t),
target,
intent(inout) :: ksp
581 type(
coef_t),
target,
intent(in) :: coef
582 type(
dofmap_t),
target,
intent(in) :: dof
583 type(
gs_t),
target,
intent(inout) :: gs
584 type(
bc_list_t),
target,
intent(inout) :: bclst
585 character(len=*) :: pctype
586 type(json_file),
intent(inout) :: pcparams
588 call precon_factory(pc, pctype)
590 select type (pcp => pc)
592 call pcp%init(coef, dof, gs)
594 call pcp%init(coef, dof, gs)
596 call pcp%init(coef, dof, gs)
598 call pcp%init(coef, bclst, pcparams)
603 end subroutine scalar_scheme_precon_factory
609 subroutine scalar_scheme_update_material_properties(this, time)
610 class(scalar_scheme_t),
intent(inout) :: this
612 type(
field_t),
pointer :: nut, alphat
615 type(
field_t),
pointer :: lambda_factor
617 call this%user_material_properties(this%name, this%material_properties, &
621 if (len_trim(this%nut_field_name) .gt. 0 &
622 .and. len_trim(this%alphat_field_name) .eq. 0 )
then
627 call field_cmult2(lambda_factor, nut, 1.0_rp / this%pr_turb)
630 call field_add3(this%lambda_tot, this%lambda, lambda_factor)
633 else if (len_trim(this%alphat_field_name) .gt. 0 &
634 .and. len_trim(this%nut_field_name) .eq. 0 )
then
639 call field_col3(lambda_factor, this%cp, alphat)
641 call field_add3(this%lambda_tot, this%lambda, lambda_factor)
644 else if (len_trim(this%alphat_field_name) .gt. 0 &
645 .and. len_trim(this%nut_field_name) .gt. 0 )
then
646 call neko_error(
"Conflicting definition of eddy diffusivity " // &
647 "for the scalar equation")
659 end subroutine scalar_scheme_update_material_properties
664 subroutine scalar_scheme_set_material_properties(this, params, user)
665 class(scalar_scheme_t),
intent(inout) :: this
666 type(json_file),
intent(inout) :: params
667 type(
user_t),
target,
intent(in) :: user
668 character(len=LOG_SIZE) :: log_buf
671 real(kind=
rp) :: const_cp, const_lambda
679 call neko_registry%add_field(this%dm_Xh, this%name //
"_lambda")
680 call neko_registry%add_field(this%dm_Xh, this%name //
"_lambda_tot")
682 this%lambda =>
neko_registry%get_field(this%name //
"_lambda")
683 this%lambda_tot =>
neko_registry%get_field(this%name //
"_lambda_tot")
686 call this%material_properties%init(2)
687 call this%material_properties%assign(1, this%cp)
688 call this%material_properties%assign(2, this%lambda)
690 if (.not.
associated(
user%material_properties, dummy_mp_ptr))
then
692 write(log_buf,
'(A)')
"Material properties must be set in the user " // &
695 this%user_material_properties =>
user%material_properties
697 call user%material_properties(this%name, this%material_properties, time)
700 if (params%valid_path(
'Pe') .and. &
701 (params%valid_path(
'lambda') .or. &
702 params%valid_path(
'cp')))
then
703 call neko_error(
"To set the material properties for the scalar, " // &
704 "either provide Pe OR lambda and cp in the case file.")
706 else if (params%valid_path(
'Pe'))
then
707 write(log_buf,
'(A)')
'Non-dimensional scalar material properties' //&
710 write(log_buf,
'(A)')
'Specific heat capacity will be set to 1, '
712 write(log_buf,
'(A)')
'conductivity to 1/Pe. Assumes density is 1.'
717 write(log_buf,
'(A,ES13.6)')
'Pe :', const_lambda
723 const_lambda = 1.0_rp/const_lambda
732 if (
associated(
user%material_properties, dummy_mp_ptr))
then
737 write(log_buf,
'(A,ES13.6)')
'lambda :', const_lambda
739 write(log_buf,
'(A,ES13.6)')
'cp :', const_cp
754 end subroutine scalar_scheme_set_material_properties
760 subroutine scalar_scheme_wrapper_init(this, msh, coef, gs, params, &
761 numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
762 class(scalar_scheme_wrapper_t),
intent(inout) :: this
763 type(
mesh_t),
target,
intent(in) :: msh
764 type(
coef_t),
target,
intent(in) :: coef
765 type(
gs_t),
target,
intent(inout) :: gs
766 type(json_file),
target,
intent(inout) :: params
767 type(json_file),
target,
intent(inout) :: numerics_params
768 type(
user_t),
target,
intent(in) :: user
769 type(
chkp_t),
target,
intent(inout) :: chkp
772 type(
field_t),
target,
intent(in) :: rho
775 call scalar_scheme_factory(this%scalar, msh, coef, gs, params, &
778 end subroutine scalar_scheme_wrapper_init
781 subroutine scalar_scheme_wrapper_free(this)
782 class(scalar_scheme_wrapper_t),
intent(inout) :: this
784 if (
allocated(this%scalar))
then
785 call this%scalar%free()
786 deallocate(this%scalar)
789 end subroutine scalar_scheme_wrapper_free
795 subroutine scalar_scheme_wrapper_move_from(this, other)
796 class(scalar_scheme_wrapper_t),
intent(inout) :: this
797 class(scalar_scheme_wrapper_t),
intent(inout) :: other
800 call move_alloc(other%scalar, this%scalar)
802 end subroutine scalar_scheme_wrapper_move_from
806 function scalar_scheme_wrapper_is_allocated(this)
result(is_alloc)
807 class(scalar_scheme_wrapper_t),
intent(in) :: this
809 is_alloc =
allocated(this%scalar)
810 end function scalar_scheme_wrapper_is_allocated
Copy data between host and device (or device and device)
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 to dealocate a scalar formulation.
Abstract interface to initialize a scalar formulation.
Abstract interface to restart a scalar formulation.
Abstract interface to compute a time-step.
Abstract interface for setting material properties.
Defines a boundary condition.
Jacobi preconditioner accelerator backend.
Device abstraction, common interface for various accelerators.
integer, parameter, public device_to_host
Defines a mapping of the degrees of freedom.
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_col3(a, b, c, n)
Vector multiplication with 3 vectors .
subroutine, public field_copy(a, b, n)
Copy a vector .
subroutine, public field_add3(a, b, c, n)
Vector addition .
Contains the field_serties_t type.
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
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Defines a registry for storing solution fields.
type(registry_t), target, public neko_registry
Global field registry.
Contains the scalar_scheme_t type.
subroutine scalar_scheme_init(this, msh, c_xh, gs_xh, params, scheme, user, rho)
Initialize all related components of the current scheme.
logical function scalar_scheme_wrapper_is_allocated(this)
Return allocation status.
subroutine scalar_scheme_wrapper_free(this)
Destructor. Just deallocates the pointer.
subroutine scalar_scheme_free(this)
Deallocate a scalar formulation.
subroutine scalar_scheme_validate(this)
Validate that all fields, solvers etc necessary for performing time-stepping are defined.
subroutine scalar_scheme_update_material_properties(this, time)
Call user material properties routine and update the values of lambda if necessary.
subroutine scalar_scheme_set_material_properties(this, params, user)
Set lamdba and cp.
subroutine scalar_scheme_wrapper_init(this, msh, coef, gs, params, numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
Constructor. Initializes the object.
subroutine scalar_scheme_wrapper_move_from(this, other)
Move assignment operator for the wrapper, needed for storing schemes in lists and arrays.
Implements the scalar_source_term_t type.
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.
Jacobi preconditioner SX-Aurora backend.
Compound scheme for the advection and diffusion operators in a transport equation.
Base class for time integration schemes.
Module with things related to the simulation time.
Implements type time_step_controller.
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.
field_list_t, To be able to group fields together
Stores a series (sequence) of fields, logically connected to a base field, and arranged according to ...
Defines a jacobi preconditioner.
Type for storing initial and final residuals in a Krylov solver.
Base abstract type for a canonical Krylov method, solving .
Defines a canonical Krylov preconditioner.
Base type for a scalar advection-diffusion solver.
A helper type that is needed to have an array of polymorphic objects.
Wrapper contaning and executing the scalar source terms.
The function space for the SEM solution fields.
Defines a jacobi preconditioner for SX-Aurora.
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.
Provides a tool to set time step dt.
A type collecting all the overridable user routines and flag to suppress type injection from custom m...