Spatial Modifiers


Modifications applied to the result of a classical scheme such as limiter or filters have to be implemented in _Solver/SpatialModifiers as a module. The modifier has to be implemented within a class whose name should match its type (e.g “Limiter” or “Filtering”), and the main method that is intended for modifying the residuals should be named Mollify in any case.


From a user perspective, it is possible to apply several mollifiers to the original residual by specifying the wished sequence as a litteral string. The input will then be interpreted by the class Amendements upon the instance creation. In practice, the input is parsed by the method ExtractSequence and converted to the actual modifyer by StringToResu.



Page navigation




Usage


Assuming one has computed a residual resu at some time setp n, and likes to modify it by the means of some limiter or filter/streamlines dissipator, the wish of the amendement sequence to use and the related parameter has to be given in the settings file.


Specifying the sequence of wished modification of the residual

In the settings file (see the User Manual) indicate for the amendments the composition of modifications in the natural way showed below, where $Limit and $Filter are the indices of the desired mollifiers in the mapping file Mappers.py, the values between the symbols <> are the arguments to give to the desired mollifiers (see below). Check each mollifier’s definition to know the required ones, if any. If no parameter is required, writing <> is not neccessary.

# No limiter
0

# A single modification
Limiter_$Limit

# A single modification that should be added to the original residual
Original+Filtering_$Filter

# A single modification with argument to give to the mollifier
Limiter_$Limit<listof,parameters>

# Or a combination
Limiter_$Limit+Filtering_$Filter<some,parameters,topass>

More complicated composition can be performed in a similar way. In those cases, only the most nested mollifiers will be applied to the original residuals (those are the terms not followed by parenthesis).

Limiter_3<2>(Limiter_1<1.2>(Original+Filtering_2)+Filtering_1<p>(Limiter_1<2.2>))*Limiter_2

Warning

  • There should be NO spurious spaces in the definition of the modification strings.



Specifying arguments

The list of arguments to pass to the mollifier should be given between the symbols <>, in the order required by the selected mollifier (see below for their definition). The types of the arguments (float, integer, string, etc) may be different, but will be passed to the related function as strings.


By example, if one limiter takes a float and a string as arguments, one would write in the parameter file

Limiter_1<14.2,thestring>

Note

When the required argument is of type Quadrature, a specific writing format should be used. The mapping from the argument to the asked quadrature instance is made automatically during the argument parsing in the StringToResu routine, (converting the string to quadrature there rather than in the modifier class avoids an import loop at the technical level).


To allow this automatic mapping, a specific format for the quadrature argument has to be kept: it should be

#InnerQuadratureSettings##BoundaryQuadratureSettings#

where the parameters InnerQuadratureSettings and BoundaryQuadratureSettings have to be given as in Step 3: Design your settings, but with underscores separating the given properties instead of spaces. By example, if one would like to select the 5th scheme of Williams Shunn Jameson from the Quadpy library within the elements and a hand base quadrature on the boundary, the respective quadrature settings for the filtering method would read

# InnerQuadratureSettings
2_15_5

# BoundaryQuadratureSettings
1

and one would write in the parameter file

Filtering_1<#2_15_5##1#,other,parameters>

Note

The full list of the available quadrature is given at Quadrature Rules.





Modifiers sequence generation


The module SpatialAmendements located in the file _Solver/SpatialModifiers/SpatialAmendements.py`furnishes a single class :code:`Amendements that interfaces the user instructions for limiting and amendments of the schemes to the actuals amemdement routine, respecting the sequence and order of mollifications.

Upon instance creation, it converts the string-user input defining the sequence of modifications to a unique equivalent modification routine that is to apply on any freshly computed residual. To do so, it converts the string to a sequence of unitary application process and combines them according to the given syntax. The method Mollify is therefore the generic interface to apply any sequence of mollifyer to a spatial scheme.



class 3_LESARD.SourceCode._Solver.SpatialModifiers.SpatialAmendements.Amendements(AmendementsSequence, Problem, Mesh)

Class furnishing all the schemes tools and registers the scheme’s wished parameters. The main iteration is accessible with the method “Mollify”.

Fields

Upon instance creation, the given parameters are repeated within the structure with identical names (see below). It also has for further fields:

  •   MolliProcess (list of callbacks) – the list of amendements to apply to the designated sequential intermediate residuals


Parameters
  • AmendementsSequence (string) – the user-understandable list of amendements to apply to the designated sequential intermediate residuals

  • Problem (Problem) – the considered problem

  • Mesh (MeshStructure) – the considered mesh instance


Methods

ExtractSequence()

Converting the string user input into a tabular giving the sequence of modification to apply

Parameters

None – self.AmendementsSequence of the structure should be filled

Returns

the unitary modification list as strings, each in the format ModifyierType_Index(ApplicationTarget)

Return type

Tab (string list)

StringToResu(stringList)

Interfacing the user wishes with the actual routines, according to the unitary modification sequence list.

Parameters

StringList (list of strings) – the unitary modification list, each in the format ModifyierType_Index(ApplicationTarget)

Returns

fills directly the self.MolliProcess of the structure, a tabular containing each step of modification [modifier, ApplicationTarget]

Return type

None

Mollify(FlagPoints, SolBuff, resu, fluxes, i, du=0, dt=1)

Limting routine according to the given sequence stored in the structure

Parameters
  • FlagPoints (function callback) – the handle of a function flagging the points to the relevant fluid

  • Solbuff (solution structure) – structure containing the current solution’s values to iterate

  • resu (float numpy array) – previoulsy computed residuals that have to be modified (NbVars x NbElementDofs)

  • fluxes (numpy (multiD)array) – pre-computed fluxes at the points of interest. For this scheme, access with fluxes[element, face, coordinate(fx or fy), variable, pointindex]

  • i (integer) – the index of the considered element within which the partial residuals will be computed

  • du (float numpy array) – (optional) when using DeC, the difference in the time iteration

  • dt (float) – (optional) when using DeC, the full time step

Returns

the modified residuals (NbVars x NbElementDofs)

Return type

Lresu (float numpy array)

Warning

the fluxes are not changing through the intermediary steps




Available limiters


class 3_LESARD.SourceCode._Solver.SpatialModifiers.Limiter_Psi.Limiter(Problem, Mesh, *Params)

Class furnishing all the schemes tools and registers the scheme’s wished parameters. The main iteration is accessible with the method “Mollify”.

Fields

All the fields given as arguments at the instance’s creation are repeated within the structure with the same names (see the parameters below).


Parameters
  • Problem (Problem) – the considered problem

  • Mesh (MeshStructure) – the considered mesh instance

  • Params (string list) – the limiters’s parameters as wished by the user


Methods

Mollify(FlagPoints, SolBuff, resu, fluxes, i, du=0, dt=1)

Limting routine according to the Psi rule

Parameters
  • FlagPoints (function callback) – the handle of a function flagging the points to the relevant fluid (not used here)

  • Solution (solution structure) – structure containing the current solution’s values to iterate

  • resu (float numpy array) – previoulsy computed residuals that have to be modified (NbVars x NbElementDofs)

  • fluxes (numpy (multiD)array) – pre-computed fluxes at the points of interest. For this scheme, access with fluxes[element, face, coordinate(fx or fy), variable, pointindex]

  • i (integer) – the index of the considered element within which the partial residuals will be computed

  • du (float numpy array) – (optional) when using DeC, the difference in the time iteration

  • dt (float) – (optional) when using DeC, the full time step

Returns

the limited residuals (NbVars x NbElementDofs)

Return type

Lresu (float numpy array)




Available filters


class 3_LESARD.SourceCode._Solver.SpatialModifiers.Filtering_Streamline.Filtering(Problem, Mesh, QuadratureRules, *Params)

Class furnishing all the schemes tools and registers the scheme’s wished parameters. The main iteration is accessible with the method “Mollify”.

Fields

All the fields given as arguments at the instance’s creation are repeated within the structure with the same names (see the parameters below). It additionaly has upon instance creation the fields:

  •   Order (integer) – The order of the used basis functions (matching the order of the mesh)

  •   InnerQBFValues (float numpy array) – Value of the basis function at the quadrature points of each element (NbInnerElementsx NbQuadraturePoints x NbBasisFunc)

  •   InnerQWeigths (float numpy array) – Weights of the quadrature points of each element (NbInnerElements x NbQuadraturePoints)

  •   InnerQGFValues (float numpy array) – Value of the gradient of the basis function at the quadrature points of each element (NbInnerElements x NbQuadraturePoints x NbBasisFunc x dim)

  •   InnerQPoints (float numpy array) – Quadrature points of each element (NbInnerElements x NbQuadraturePoints)

  •   FaceWiseQWeigths (float numpy array) – Weights of the quadrature points of each element’s face (NbInnerElements x NbFace x NbQuadraturePoints)

  •   FaceWiseQPoints (float numpy array) – Quadrature points of each element’s face (NbInnerElements x NbFace x NbQuadraturePoints)

  •   FaceWiseQBFValues (float numpy array) – Value of the basis function at the quadrature points of each element’s face (NbInnerElements x NbFace x NbQuadraturePoints x NbBasisFunc)


Warning

The streamline stabilisation only seems to work for advection so far, it crashes on euler


Parameters
  • Problem (Problem) – the considered problem

  • Mesh (MeshStructure) – the considered mesh instance

  • QuadratureRules (Quadratures) – the quadrature instance with which the scheme will be performed (madatory paramter)

  • Params (string list) – the limiters’s parameters as wished by the user (useless here)


Methods

Mollify(FlagPoints, Solution, resu, fluxes, i, du=0, dt=1)

Streamline dissipation, that has to be ADDED to the residual. In the input line, do not simply write Filtering_1, but write 0+Filtering_1 so that it gets added to the initial residual

Parameters
  • FlagPoints (function callback) – the handle of a function flagging the points to the relevant fluid

  • Solution (solution structure) – structure containing the current solution’s values to iterate

  • resu (float numpy array) – previoulsy computed residuals that have to be modified (NbVars x NbElementDofs)

  • fluxes (numpy (multiD)array) – pre-computed fluxes at the points of interest. For this scheme, access with fluxes[element, face, coordinate(fx or fy), variable, pointindex]

  • i (integer) – the index of the considered element within which the partial residuals will be computed

  • du (float numpy array) – (optional) when using DeC, the difference in the time iteration

  • dt (float) – (optional) when using DeC, the full time step

Returns

the limited residuals (NbVars x NbElementDofs)

Return type

Lresu (float numpy array)

PreEvaluateQuadratureQuantities()

Precomputing and storing the values of the basis functions at the quadrature points of each element

Parameters

None – Considers the given parameters of the class directly

Returns

Creates the fields InnerQBFValues, InnerQPoints, InnerQWeigths and InnerQGFValues in the modifier structure, containing the values of each basis function at each the pointsm, the wigths and the basis functions values there for each element.

Return type

None