Neko
0.8.1
A portable framework for high-order spectral element flow simulations
|
Contains the scalar_scheme_t type. More...
Data Types | |
type | scalar_scheme_t |
Base type for a scalar advection-diffusion solver. More... | |
interface | scalar_scheme_init_intrf |
Abstract interface to initialize a scalar formulation. More... | |
interface | scalar_scheme_restart_intrf |
Abstract interface to restart a scalar formulation. More... | |
interface | scalar_scheme_free_intrf |
Abstract interface to dealocate a scalar formulation. More... | |
interface | scalar_scheme_step_intrf |
Abstract interface to compute a time-step. More... | |
Functions/Subroutines | |
subroutine | scalar_scheme_add_bcs (this, zones, bc_labels) |
Initialize boundary conditions. More... | |
subroutine | scalar_scheme_init (this, msh, c_Xh, gs_Xh, params, scheme, user, material_properties) |
Initialize all related components of the current scheme. More... | |
subroutine | scalar_scheme_free (this) |
Deallocate a scalar formulation. More... | |
subroutine | scalar_scheme_validate (this) |
Validate that all fields, solvers etc necessary for performing time-stepping are defined. More... | |
subroutine | scalar_scheme_solver_factory (ksp, n, solver, max_iter, abstol) |
Initialize a linear solver. More... | |
subroutine | scalar_scheme_precon_factory (pc, ksp, coef, dof, gs, bclst, pctype) |
Initialize a Krylov preconditioner. More... | |
subroutine | scalar_scheme_set_user_bc (this, usr_eval) |
Initialize a user defined scalar bc. More... | |
Contains the scalar_scheme_t type.
subroutine scalar_scheme::scalar_scheme_add_bcs | ( | class(scalar_scheme_t), intent(inout) | this, |
type(facet_zone_t), dimension(neko_msh_max_zlbls), intent(inout) | zones, | ||
character(len=neko_msh_max_zlbl_len), dimension(:), intent(in) | bc_labels | ||
) |
Initialize boundary conditions.
zones | List of zones |
bc_labels | List of user specified bcs from the parameter file currently dirichlet 'd=X' and 'user' supported |
Check if user bc on this zone
Definition at line 227 of file scalar_scheme.f90.
subroutine scalar_scheme::scalar_scheme_free | ( | class(scalar_scheme_t), intent(inout) | this | ) |
Deallocate a scalar formulation.
Definition at line 453 of file scalar_scheme.f90.
subroutine scalar_scheme::scalar_scheme_init | ( | class(scalar_scheme_t), intent(inout), target | this, |
type(mesh_t), intent(inout), target | msh, | ||
type(coef_t), intent(inout), target | c_Xh, | ||
type(gs_t), intent(inout), target | gs_Xh, | ||
type(json_file), intent(inout), target | params, | ||
character(len=*), intent(in) | scheme, | ||
type(user_t), intent(in), target | user, | ||
type(material_properties_t), intent(inout), target | material_properties | ||
) |
Initialize all related components of the current scheme.
msh | The mesh. |
c_Xh | The coefficients. |
gs_Xh | The gather-scatter. |
params | The case parameter file in json. |
scheme | The name of the scalar scheme. |
user | Type with user-defined procedures. |
Definition at line 298 of file scalar_scheme.f90.
subroutine scalar_scheme::scalar_scheme_precon_factory | ( | class(pc_t), intent(inout), allocatable, target | pc, |
class(ksp_t), intent(inout), target | ksp, | ||
type(coef_t), intent(inout), target | coef, | ||
type(dofmap_t), intent(inout), target | dof, | ||
type(gs_t), intent(inout), target | gs, | ||
type(bc_list_t), intent(inout), target | bclst, | ||
character(len=*) | pctype | ||
) |
Initialize a Krylov preconditioner.
Definition at line 552 of file scalar_scheme.f90.
subroutine scalar_scheme::scalar_scheme_set_user_bc | ( | class(scalar_scheme_t), intent(inout) | this, |
procedure(usr_scalar_bc_eval) | usr_eval | ||
) |
Initialize a user defined scalar bc.
usr_eval | User specified boundary condition for scalar field |
Definition at line 589 of file scalar_scheme.f90.
subroutine scalar_scheme::scalar_scheme_solver_factory | ( | class(ksp_t), intent(inout), allocatable, target | ksp, |
integer, intent(in), value | n, | ||
character(len=*), intent(in) | solver, | ||
integer, intent(in) | max_iter, | ||
real(kind=rp) | abstol | ||
) |
Initialize a linear solver.
Definition at line 540 of file scalar_scheme.f90.
subroutine scalar_scheme::scalar_scheme_validate | ( | class(scalar_scheme_t), intent(inout), target | this | ) |
Validate that all fields, solvers etc necessary for performing time-stepping are defined.
Definition at line 496 of file scalar_scheme.f90.