Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
scalar_scheme::scalar_scheme_t Type Referenceabstract

Base type for a scalar advection-diffusion solver. More...

Inheritance diagram for scalar_scheme::scalar_scheme_t:
Collaboration diagram for scalar_scheme::scalar_scheme_t:

Public Member Functions

procedure, pass(thisscheme_init (this, msh, c_xh, gs_xh, params, scheme, user, rho)
 Constructor for the base type.
 
procedure, pass(thisscheme_free (this)
 Destructor for the base type.
 
procedure, pass(thisvalidate (this)
 Validate successful initialization.
 
procedure, pass(thisset_material_properties (this, params, user)
 Set lambda and cp.
 
procedure, pass(thisupdate_material_properties (this, time)
 Update variable material properties.
 
procedure(scalar_scheme_init_intrf), deferred, pass init (this, msh, coef, gs, params, numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
 Constructor.
 
procedure(scalar_scheme_free_intrf), deferred, pass free (this)
 Destructor.
 
procedure(scalar_scheme_step_intrf), deferred, pass step (this, time, ext_bdf, dt_controller, ksp_results)
 Solve for the current timestep.
 
procedure(scalar_scheme_restart_intrf), deferred, pass restart (this, chkp)
 Restart from a checkpoint.
 

Public Attributes

character(len=:), allocatable name
 A name that can be used to distinguish this solver in e.g. user routines.
 
type(field_t), pointer u
 x-component of Velocity
 
type(field_t), pointer v
 y-component of Velocity
 
type(field_t), pointer w
 z-component of Velocity
 
type(field_t), pointer s
 The scalar.
 
type(field_series_tslag
 Lag arrays, i.e. solutions at previous timesteps.
 
type(space_t), pointer xh
 Function space \( X_h \).
 
type(dofmap_t), pointer dm_xh
 Dofmap associated with \( X_h \).
 
type(gs_t), pointer gs_xh
 Gather-scatter associated with \( X_h \).
 
type(coef_t), pointer c_xh
 Coefficients associated with \( X_h \).
 
type(field_t), pointer f_xh => null()
 Right-hand side.
 
type(scalar_source_term_tsource_term
 The source term for equation.
 
class(ksp_t), allocatable ksp
 Krylov solver.
 
integer ksp_maxiter
 Max iterations in the Krylov solver.
 
integer projection_dim
 Projection space size.
 
integer projection_activ_step
 
class(pc_t), allocatable pc
 Preconditioner.
 
type(bc_list_tbcs
 List of boundary conditions, including the user one.
 
type(json_file), pointer params
 Case paramters.
 
type(mesh_t), pointer msh => null()
 Mesh.
 
type(chkp_t), pointer chkp => null()
 Checkpoint for restarts.
 
character(len=:), allocatable nut_field_name
 The turbulent kinematic viscosity field name.
 
type(field_t), pointer rho => null()
 Density.
 
type(field_t), pointer lambda => null()
 Thermal diffusivity.
 
type(field_t), pointer cp => null()
 Specific heat capacity.
 
type(field_t), pointer lambda_tot => null()
 Total diffusivity.
 
real(kind=rp) pr_turb
 Turbulent Prandtl number.
 
type(field_list_tmaterial_properties
 Field list with cp and lambda.
 

Static Public Attributes

procedure(user_material_properties_intf), pointer, nopass user_material_properties => null()
 

Detailed Description

Definition at line 83 of file scalar_scheme.f90.

Member Function/Subroutine Documentation

◆ free()

procedure(scalar_scheme_free_intrf), deferred, pass scalar_scheme::scalar_scheme_t::free ( class(scalar_scheme_t), intent(inout this)
pure virtual

Definition at line 158 of file scalar_scheme.f90.

◆ init()

procedure(scalar_scheme_init_intrf), deferred, pass scalar_scheme::scalar_scheme_t::init ( class(scalar_scheme_t), intent(inout), target  this,
type(mesh_t), intent(in), target  msh,
type(coef_t), intent(in), target  coef,
type(gs_t), intent(inout), target  gs,
type(json_file), intent(inout), target  params,
type(json_file), intent(inout), target  numerics_params,
type(user_t), intent(in), target  user,
type(chkp_t), intent(inout), target  chkp,
type(field_series_t), intent(in), target  ulag,
type(field_series_t), intent(in), target  vlag,
type(field_series_t), intent(in), target  wlag,
type(time_scheme_controller_t), intent(in), target  time_scheme,
type(field_t), intent(in), target  rho 
)
pure virtual

Definition at line 156 of file scalar_scheme.f90.

◆ restart()

procedure(scalar_scheme_restart_intrf), deferred, pass scalar_scheme::scalar_scheme_t::restart ( class(scalar_scheme_t), intent(inout), target  this,
type(chkp_t), intent(inout chkp 
)
pure virtual

Definition at line 162 of file scalar_scheme.f90.

◆ scheme_free()

procedure, pass(this) scalar_scheme::scalar_scheme_t::scheme_free ( class(scalar_scheme_t), intent(inout this)

Definition at line 146 of file scalar_scheme.f90.

◆ scheme_init()

procedure, pass(this) scalar_scheme::scalar_scheme_t::scheme_init ( class(scalar_scheme_t), intent(inout), target  this,
type(mesh_t), intent(in), target  msh,
type(coef_t), intent(in), 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(field_t), intent(in), target  rho 
)
Parameters
mshThe mesh.
c_XhThe coefficients.
gs_XhThe gather-scatter.
paramsThe parameter dictionary in json.
schemeThe name of the scalar scheme.
userType with user-defined procedures.
rhoThe density of the fluid.

Definition at line 144 of file scalar_scheme.f90.

◆ set_material_properties()

procedure, pass(this) scalar_scheme::scalar_scheme_t::set_material_properties ( class(scalar_scheme_t), intent(inout this,
type(json_file), intent(inout params,
type(user_t), intent(in), target  user 
)
Parameters
paramsThe case file configuration dictionary.
userThe user interface.

Definition at line 150 of file scalar_scheme.f90.

◆ step()

procedure(scalar_scheme_step_intrf), deferred, pass scalar_scheme::scalar_scheme_t::step ( class(scalar_scheme_t), intent(inout this,
type(time_state_t), intent(in time,
type(time_scheme_controller_t), intent(in ext_bdf,
type(time_step_controller_t), intent(in dt_controller,
type(ksp_monitor_t), intent(inout ksp_results 
)
pure virtual

Definition at line 160 of file scalar_scheme.f90.

◆ update_material_properties()

procedure, pass(this) scalar_scheme::scalar_scheme_t::update_material_properties ( class(scalar_scheme_t), intent(inout this,
type(time_state_t), intent(in time 
)
Parameters
tTime value.
tstepCurrent time step.

Definition at line 153 of file scalar_scheme.f90.

◆ validate()

procedure, pass(this) scalar_scheme::scalar_scheme_t::validate ( class(scalar_scheme_t), intent(inout), target  this)

Definition at line 148 of file scalar_scheme.f90.

Member Data Documentation

◆ bcs

type(bc_list_t) scalar_scheme::scalar_scheme_t::bcs

Definition at line 119 of file scalar_scheme.f90.

◆ c_xh

type(coef_t), pointer scalar_scheme::scalar_scheme_t::c_xh

Definition at line 103 of file scalar_scheme.f90.

◆ chkp

type(chkp_t), pointer scalar_scheme::scalar_scheme_t::chkp => null()

Definition at line 125 of file scalar_scheme.f90.

◆ cp

type(field_t), pointer scalar_scheme::scalar_scheme_t::cp => null()

Definition at line 133 of file scalar_scheme.f90.

◆ dm_xh

type(dofmap_t), pointer scalar_scheme::scalar_scheme_t::dm_xh

Definition at line 99 of file scalar_scheme.f90.

◆ f_xh

type(field_t), pointer scalar_scheme::scalar_scheme_t::f_xh => null()

Definition at line 105 of file scalar_scheme.f90.

◆ gs_xh

type(gs_t), pointer scalar_scheme::scalar_scheme_t::gs_xh

Definition at line 101 of file scalar_scheme.f90.

◆ ksp

class(ksp_t), allocatable scalar_scheme::scalar_scheme_t::ksp

Definition at line 109 of file scalar_scheme.f90.

◆ ksp_maxiter

integer scalar_scheme::scalar_scheme_t::ksp_maxiter

Definition at line 111 of file scalar_scheme.f90.

◆ lambda

type(field_t), pointer scalar_scheme::scalar_scheme_t::lambda => null()

Definition at line 131 of file scalar_scheme.f90.

◆ lambda_tot

type(field_t), pointer scalar_scheme::scalar_scheme_t::lambda_tot => null()

Definition at line 135 of file scalar_scheme.f90.

◆ material_properties

type(field_list_t) scalar_scheme::scalar_scheme_t::material_properties

Definition at line 139 of file scalar_scheme.f90.

◆ msh

type(mesh_t), pointer scalar_scheme::scalar_scheme_t::msh => null()

Definition at line 123 of file scalar_scheme.f90.

◆ name

character(len=:), allocatable scalar_scheme::scalar_scheme_t::name

Definition at line 85 of file scalar_scheme.f90.

◆ nut_field_name

character(len=:), allocatable scalar_scheme::scalar_scheme_t::nut_field_name

Definition at line 127 of file scalar_scheme.f90.

◆ params

type(json_file), pointer scalar_scheme::scalar_scheme_t::params

Definition at line 121 of file scalar_scheme.f90.

◆ pc

class(pc_t), allocatable scalar_scheme::scalar_scheme_t::pc

Definition at line 117 of file scalar_scheme.f90.

◆ pr_turb

real(kind=rp) scalar_scheme::scalar_scheme_t::pr_turb

Definition at line 137 of file scalar_scheme.f90.

◆ projection_activ_step

integer scalar_scheme::scalar_scheme_t::projection_activ_step

Definition at line 115 of file scalar_scheme.f90.

◆ projection_dim

integer scalar_scheme::scalar_scheme_t::projection_dim

Steps to activate projection for ksp

Definition at line 113 of file scalar_scheme.f90.

◆ rho

type(field_t), pointer scalar_scheme::scalar_scheme_t::rho => null()

Definition at line 129 of file scalar_scheme.f90.

◆ s

type(field_t), pointer scalar_scheme::scalar_scheme_t::s

Definition at line 93 of file scalar_scheme.f90.

◆ slag

type(field_series_t) scalar_scheme::scalar_scheme_t::slag

Definition at line 95 of file scalar_scheme.f90.

◆ source_term

type(scalar_source_term_t) scalar_scheme::scalar_scheme_t::source_term

Definition at line 107 of file scalar_scheme.f90.

◆ u

type(field_t), pointer scalar_scheme::scalar_scheme_t::u

Definition at line 87 of file scalar_scheme.f90.

◆ user_material_properties

procedure(user_material_properties_intf), pointer, nopass scalar_scheme::scalar_scheme_t::user_material_properties => null()
static

Definition at line 140 of file scalar_scheme.f90.

◆ v

type(field_t), pointer scalar_scheme::scalar_scheme_t::v

Definition at line 89 of file scalar_scheme.f90.

◆ w

type(field_t), pointer scalar_scheme::scalar_scheme_t::w

Definition at line 91 of file scalar_scheme.f90.

◆ xh

type(space_t), pointer scalar_scheme::scalar_scheme_t::xh

Definition at line 97 of file scalar_scheme.f90.


The documentation for this type was generated from the following file: