Linear advection (Rotation)¶
All the necessary routines that implement a
simple 2D translation, linear advection are defined in
GoverningEquations/LinearAdvectionRotation.py
Governing equations¶
-
class
3_LESARD.SourceCode.ProblemDefinition.GoverningEquations.LinearAdvectionRotation.
Equations
(FluidProp, EOs)¶ Class furnishing all the methods that are required to define numerical schemes and evolve the solution according to the Linear Advection (translation)
Fields
It contains the following initialisation parameters (later available as well as attributes)
FluidProp (list of FluidModel) – the list of fluid property instances (NbFluids)
EoS (list of EOS) – the list of Equation of state instances (NbFluids)
It further contains the following attributes, filled upon initialisation
NbVariables (integer) – number of variables of the model
VariablesNames (list of strings) – ordered name of the variables
VariablesUnits (list of strings) – symbols of the units of the variables
VariablesLatexNames (list of strings) – ordered name of the variables, latex encoding
VariablesLatexUnits (list of strings) – symbols of the units of the variables, latex encoding
XYLatexUnits (list of strings) – units of the x-y coordinates in latex encoding
- Parameters
FluidProp (list of FluidModel) – the list of fluid property instances (NbFluids)
EOs (list of EOS) – the list of Equation of state instances (NbFluids)
Methods
-
ConservativeToPrimary
(Var, FluidIndex, *args)¶ Converts the conservative variables to the primary ones
- Parameters
Var (2D numpy array) – the variables of the problem (NbVars x NbGivenPoints)
FluidIndex (integer numpy array) – the fluid present at each given point (NbGivenPoints)
- Returns
(NbVars x NbGivenPoints) the corresponding primary variables
- Return type
Var (numpy array)
Note
- *args is there only for compatibility reason at call time
-
PrimaryToConservative
(PrimVar, FluidIndex, *args)¶ Converts the primary variables to the conservative ones
- Parameters
Var (2D numpy array) – the primary variables of the problem (NbVars x NbGivenPoints)
FluidIndex (integer numpy array) – the fluid present at each given point (NbGivenPoints)
- Returns
(NbVars x NbGivenPoints) the corresponding conservative variables
- Return type
Var (numpy array)
Note
- *args is there only for compatibility reason at call time
-
Flux
(Var, x, *args)¶ Flux function of the LinearAdvection (rotation case)
- Parameters
Var (float numpy array) – the value of the variables (NbVariables x NbGivenPoints)
x (float numpy array) – (generally optional, required for this problem), x-y coordinates of the given points (NbGivenPoints x 2)
Optional argument
– FluidIndex (integer numpy array) – the fluid present at each given point (NbGivenPoints)
- Returns
the obtained flux (spatialdim x NbVars x NbPoints)
- Return type
Flux (3D numpy array)
Note
- This function is Vectorised
- Fluid index is the fluid index in the initially given list, not the fluid type
- FluidProp is the list of all the fluid properties for each fluid (given the list, not the type)
- *args is there only for compatibility reason at call time
-
GetUV
(Var, x, *argv)¶ Function giving back the x-y velocity at the considered point, given the variational values given for the advection.
- Parameters
Var (float numpy array) – the value of the variables (NbVariables x NbGivenPoints)
x (float numpy array) – (generally optional, required for this problem), x-y coordinates of the given points (NbGivenPoints x 2)
Optional argument
– FluidIndex (integer numpy array) – the fluid present at each given point (NbGivenPoints)
- Returns
UV (float numpy array) – the velocities values at the points (2 x NbGivenPoints)
Note
- This function is Vectorised
- *args is there only for compatibility reason at call time
-
Jacobian
(Var, x, *argv)¶ Computes the Jacobian of the flux for the LinearAdvection (rotation case)
- Parameters
Var (float numpy array) – the value of the variables (NbVariables x NbGivenPoints)
x (float numpy array) – (generally optional, required for this problem), x-y coordinates of the given points (NbGivenPoints x 2)
Optional argument
– FluidIndex (integer numpy array) – the fluid present at each given point (NbGivenPoints)
- Returns
J[:,:,i] gives the jacobian of the flux taking care of the dynamic of the ith spatial coordinate.
- Return type
J (3D numpy array)
Note
- For each flux fi = (fi1,…, fin), the returned Jacobian reads:J[:,:] = [dfi1/dx1, ….., dfi1/dxn….df in/dx1, ….., dfin/dxn]
- *args is there only for compatibility reason at call time
-
SpectralRadius
(Var, FluidIndex, n, x, *argv)¶ Computes the spectral radius associated to the flux.
- Parameters
Var (2D numpy array) – the variables of the problem (NbVars x NbGivenPoints)
FluidIndex (integer numpy array) – the fluid present at each given point (NbGivenPoints)
n (2D numpy array) – the x-y values of the normal at each given point (NbGivenPoints x 2)
x (2D numpy array) – (generally optional, required for this problem) the x-y locations at which the variables are given (NbGivenPoints x 2)
- Returns
the spectral radius computed at each given point
- Return type
Lmb (numpy array)
Note
- *args is there only for compatibility reason at call time
-
EigenValues
(Var, FluidIndex, n, x, *args)¶ Computes the eigenvalues associated to the flux.
- Parameters
Var (2D numpy array) – the variables of the problem (NbVars x NbGivenPoints)
FluidIndex (integer numpy array) – the fluid present at each given point (NbGivenPoints)
n (2D numpy array) – the x-y values of the normal at each given point (NbGivenPoints x 2)
x (2D numpy array) – (generally optional, required for this problem) the x-y locations at which the variables are given (NbGivenPoints x 2)
- Returns
(NbGivenPoints x NbEigs) the eigenvalues at each given point
- Return type
lbd (numpy array)
Note
- *args is there only for compatibility reason at call time
-
RightEigenVectors
(Var, FluidIndex, n, x, *args)¶ Computes the right eigenvectors associated to the eigenvalues.
- Parameters
Var (2D numpy array) – the variables of the problem (NbVars x NbGivenPoints)
FluidIndex (integer numpy array) – the fluid present at each given point (NbGivenPoints)
n (2D numpy array) – the x-y values of the normal at each given point (NbGivenPoints x 2)
x (2D numpy array) – (generally optional, required for this problem) the x-y locations at which the variables are given (NbGivenPoints x 2)
- Returns
(NbVars x MbVars x NbGivenPoints) the matrix of eigenvectors for each given point
- Return type
reg (numpy array)
Note
- *args is there only for compatibility reason at call time
-
LeftEigenVectors
(Var, FluidIndex, n, x, *args)¶ Computes the left eigenvectors associated to the eigenvalues.
- Parameters
Var (2D numpy array) – the variables of the problem (NbVars x NbGivenPoints)
FluidIndex (integer numpy array) – the fluid present at each given point (NbGivenPoints)
n (2D numpy array) – the x-y values of the normal at each given point (NbGivenPoints x 2)
x (2D numpy array) – (generally optional, required for this problem) the x-y locations at which the variables are given (NbGivenPoints x 2)
- Returns
(NbVars x MbVars x NbGivenPoints) the matrix of eigenvectors for each given point
- Return type
reg (numpy array)
Note
- *args is there only for compatibility reason at call time
Equations of state¶
-
class
3_LESARD.SourceCode.ProblemDefinition.GoverningEquations.LinearAdvectionRotation.
EOS
(Id)¶ Class furnishing the methods giving the Equation of State that can be used in combination with the
SimpleFluid
definition. Upon creation, fills the fieldEOS
with the function callback corresponding to the wished initialisation function and the fieldEOSName
with the associated label. The implemented EoS are linked as follows. See their respective documentation for more information.- Dummy
Note
After having been initialised the instance, the EOS is accessible from the field EOS within this class.
Warning
For advection, this class is void, just here for compatibility properties in the code
- Parameters
Id (integer) – the index corresponding to the equation of fluid the user wants when considering the governing equations given in this module and the associated fluid.
Methods
-
Dummy
(Var, FluidProp, Type)¶ Dummy EOS returning the identity for the sake of the code compatibility.
- Parameters
Type (string) – desired output (here inactive)
Var (float numpy array) – the value of the variables (NbVariables x NbGivenPoints)
FluidProp (FluidModel list) – the list of the fluid properties associated to each given Point
- Returns
the variable retrieved from the EoS (here identity as dummy function)
- Return type
val
Initial conditions¶
-
class
3_LESARD.SourceCode.ProblemDefinition.GoverningEquations.LinearAdvectionRotation.
InitialConditions
(Id, *params)¶ Class furnishing the solution initialisation routines that are suitable to study the the Linear Advection Rotation Equations, defined for working on a subdomain. Access the routine computing the initial conditions through the field
IC
, filled upon creation with the function callback corresponding to the index of the wished initialisation function. The implemented initialisation method are linked as follows. See their respective documentation for more information.- ConstantState
- Not implemented
- Bump
- Parameters
Id (integer) – the index corresponding to the initialisation method the user wants when considering the governing equations given in this module.
params (list of arguments) – the (fixed) arguments to pass to the selected function
Methods
-
ConstantState
(PointsID, Mesh, EOS=[], FluidProp=[], subdomain=[], *argv)¶ Initialising the solution to a constant state on the given points, all belonging to a same subdomain.
- Parameters
PointsID (integer array-like) – the index of the points to consider
Mesh (MeshStructure) – the considered mesh
EOS (function callback) – (optional), the equation of state given by the model of the fluid that is present at the given points
FluidProp (FluidSpecs) – (optional, not used here) the the properties of the fluid present where the given points are
Subdomain (shapely multipolyon) – (optional, not used here) the shapely polygon to which the given points belong
- Returns
The initialised values at the considered points (NbVariables x NbGivenPoints)
- Return type
Init (float numpy array)
Note
*argv is here only for compatibility on call
There is usually one initialisation routine per subdomain, and the PointsID are not necessarily contigous, be careful when assigning the values back in the regular solution.
-
Bump
(PointsID, Mesh, EOS=[], FluidProp=[], subdomain=[], *args)¶ IInitialising the solution to a bump centred at the center of mass of the subdmain. (The given points should all belonging to a same subdomain).
- Parameters
PointsID (integer array-like) – the index of the points to consider
Mesh (MeshStructure) – the considered mesh
EOS (function callback) – (optional), the equation of state given by the model of the fluid that is present at the given points
FluidProp (FluidSpecs) – (optional, not used here) the the properties of the fluid present where the given points are
Subdomain (shapely multipolyon) – (optional, not used here) the shapely polygon to which the given points belong
- Returns
The initialised values at the considered points (NbVariables x NbGivenPoints)
- Return type
Init (float numpy array)
Note
*args is here only for compatibility on call
- There is usually one initialisation routine per subdomain, and the PointsID are not necessarily contigous, be careful when assigning the values back in the regular solution.
- For running smothly the inner circle diameter of the given subdomain should be at least 2