Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
scalar_ic Module Reference

Scalar initial condition.

Data Types

interface  set_scalar_ic
 

Functions/Subroutines

subroutine set_scalar_ic_int (s, coef, gs, type, params, i)
 Set scalar initial condition (builtin)
 
subroutine set_scalar_ic_usr (scheme_name, s, coef, gs, user_proc)
 Set scalar intial condition (user defined)
 
subroutine set_scalar_ic_common (s, coef, gs)
 Set scalar initial condition (common)
 
subroutine set_scalar_ic_uniform (s, ic_value)
 Uniform initial condition.
 
subroutine set_scalar_ic_point_zone (s, base_value, zone_name, zone_value)
 Point zone initial condition.
 
subroutine set_scalar_ic_fld (s, file_name, interpolate, tolerance, mesh_file_name, i)
 Set the initial condition of the scalar based on a field. @detail The field is read from an fld file. If enabled, interpolation is also possible. In that case, the mesh coordinates can be read from another file in the fld field series.
 

Function/Subroutine Documentation

◆ set_scalar_ic_common()

subroutine scalar_ic::set_scalar_ic_common ( type(field_t), intent(inout s,
type(coef_t), intent(in coef,
type(gs_t), intent(inout gs 
)
private

Finalize scalar initial condition by distributing the initial condition across elements and multiplying by the coefficient (if any).

Parameters
sScalar field.
coefCoefficient.
gsGather-Scatter object.

Definition at line 163 of file scalar_ic.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_scalar_ic_fld()

subroutine scalar_ic::set_scalar_ic_fld ( type(field_t), intent(inout s,
character(len=*), intent(in file_name,
logical, intent(in interpolate,
real(kind=rp), intent(in tolerance,
character(len=*), intent(inout mesh_file_name,
integer, intent(in i 
)
private
Parameters
sThe scalar field.
file_nameThe name of the "fld" file series.
sample_idxindex of the field file .f000* to read, default is -1.
interpolateFlag to indicate wether or not to interpolate the values onto the current mesh.
toleranceIf interpolation is enabled, tolerance for finding the points in the mesh.
mesh_file_nameIf interpolation is enabled, name of the field file series where the mesh coordinates are located.
iIndex of the scalar field.

Definition at line 255 of file scalar_ic.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_scalar_ic_int()

subroutine scalar_ic::set_scalar_ic_int ( type(field_t), intent(inout s,
type(coef_t), intent(in coef,
type(gs_t), intent(inout gs,
character(len=*)  type,
type(json_file), intent(inout params,
integer, intent(in i 
)
private

Set scalar initial condition using one of the builtin types currently supported:

  • uniform
  • point zone
  • field
    Parameters
    sScalar field.
    coefCoefficient.
    gsGather-Scatter object.
    typeType of initial condition.
    paramsJSON parameters.
    iIndex of the scalar field.

Definition at line 83 of file scalar_ic.f90.

◆ set_scalar_ic_point_zone()

subroutine scalar_ic::set_scalar_ic_point_zone ( type(field_t), intent(inout s,
real(kind=rp), intent(in base_value,
character(len=*), intent(in zone_name,
real(kind=rp), intent(in zone_value 
)
private

Set scalar initial condition to a uniform value across a point zone.

Parameters
sScalar field.
base_valueBase value of the scalar field.
zone_nameName of the point zone.
zone_valueDesired value of the scalar field in the point zone.

Definition at line 214 of file scalar_ic.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_scalar_ic_uniform()

subroutine scalar_ic::set_scalar_ic_uniform ( type(field_t), intent(inout s,
real(kind=rp), intent(in ic_value 
)
private

Set scalar initial condition to a uniform value across the domain.

Parameters
sScalar field.
ic_valueDesired value of the scalar field.

Definition at line 189 of file scalar_ic.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_scalar_ic_usr()

subroutine scalar_ic::set_scalar_ic_usr ( character(len=*), intent(in scheme_name,
type(field_t), intent(inout), target  s,
type(coef_t), intent(in coef,
type(gs_t), intent(inout gs,
procedure(user_initial_conditions_intf user_proc 
)
private

Set scalar initial condition using a user defined function.

Parameters
scheme_nameName of the scheme calling the user routine.
sScalar field.
coefCoefficient.
gsGather-Scatter object.
user_procUser defined initial condition function.

Definition at line 139 of file scalar_ic.f90.