Fluid Selectors


The modules defining the fluid selectors are located in the folder _Solver/FluidSelectors.

  • Each module defines a fluid selector and should comport a class FS

  • The spatial scheme and motion equation associated to the fluid selector should be implemented within the class as well.

  • The class should at least contain the methods Initialiser, Iteration, UpdateFSVAttributes (see below for an example of a format).


Already implemented fluid selectors




Naive Level Set + Continuous galerkin + Lax-Friedrichs


class 3_LESARD.SourceCode._Solver.FluidSelectors.NaiveLS_CG_LFX.FS(Problem, Mesh, QuadratureRules, *Params)

Classical level set approach, initialised according to a signed distance and iterated with a continuous galerkin - lfx scheme, without redistancing.

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 also contains the further field

  •   FluidSelectorName (string) – the name of the level set

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

  •   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)


Parameters
  • Problem (Problem) – the considered problem

  • Mesh (MeshStructure) – the considered mesh instance

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


Methods

Initialiser(Solution)

External initialisation routine.

Parameters

Solution (Solution structure) – the solution instance the scheme is working on

Returns

the level set values at the MeshVertex and Dofs, to be copied to the Solution structure (NbSubdomains x (NbMeshPoints+NbDofs))

Return type

LSValues (float numpy array)

Flux(Var, Motion, *args)

Flux function of the (conservative) EulerEquations

Parameters
  • Var (float numpy array) – the value of the Level Set values (NbFluids x NbGivenPoints)

  • Motion (float numpy array) – the velocity array (2xNbPoints)

Returns

the obtained flux (spatialdim x NbVars x NbPoints)

Return type

Flux (3D numpy array)

Note

*args is there for compatibility reason when e.g. a location xy is given when called

Jacobian(Var, Motion, *args)

Computes the Jacobian of the flux for the advection of the LevelSet

Parameters
  • Var (float numpy array) – the value of the Level Set values (NbFluids x NbGivenPoints)

  • Motion (float numpy array) – the velocity array (2xNbPoints)

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 for compatibility reason when e.g. a location xy is given when called the x-y locations at which the variables are given (NbGivenPoints x 2)

SpectralRadius(Var, Motion, n, *args)

Computes the spectral radius associated to the motion flux.

Parameters
  • Var (2D numpy array) – the variables of the problem (NbVars x NbGivenPoints)

  • Motion (float numpy array) – the velocity array (2xNbPoints)

  • n (2D numpy array) – the x-y values of the normal at each given point (NbGivenPoints x 2)

Returns

the array containing the spectral radius at each given point

Return type

Lmb (numpy array)

Note

*args is there for compatibility reason when e.g. a the xy locations at which the variables are given (NbGivenPoints x 2)

FlagDofs(Solution)

Flags the degrees of freedom according to the different values of the associated level set.

Parameters

Solution (Solution structure) – the solution instance the scheme is working on

Returns

updating directly the value of FluidFlag in the solution instance

Return type

None

FlagPoints(LSValues, CellId=- 1, x=[], FaceId=- 1)

Flags the degrees of freedom according to the different values of the associated level set.

Parameters
  • LSValues (float numpy multidarray) – the vector within which the values of the level set at dofs are stored (full vector in the same format as stored in Solution.LSValues)

  • CellId (integer) – the id of the primal cell the points under consideration lie on (important esp. for points on the cell’s boundary)

  • x (2D numpy array) – (optional) barycentric coordinates of the point under consideration, ordered as the physical vertices of the element (NbPoints x NbBarycoor) If not given the FlagDofs will be simply computed at the Dofs of the CellId.

  • FaceId (integer) – When CellId and x are given, the face where the given point lie (important for code speedup only)

Returns

If CellId and x are given, the list of indices corresponding to the fluid located at the positions given in x.

If CellId is given but not x, the list of indices corresponding to the fluid located at the Dofs (according to their FV representation) If no CellId neither x is given, the vector within which the values of the level set at dofs are stored (full vector in the same format as stored in LSValues)

Return type

Flags (integer list)

ReconstructFSAtVertices(Solution)

Routine that maps the values from the Dofs to the Physical vertices of the mesh of the FluidSpotter and flags.

Parameters

Solution (Solution) – the currently being computed solution

Returns

fills direclty the RSol and RFluidFlag values in the data structure

Return type

None

Note

This is only for later plotting purposes.

UpdateSubdomains(Solution)

Method modifying the subdomains associated to the Solution upon the values of LSValues.

Parameters

Solution (Solution) – the solution instance the scheme is working on

Returns

updating directly the subdomains in the solution instance

Return type

None

Warning

Not safe yet if the first path given out is internal

Redistancing(Solution)

Redistances the Level set values upon the computed subdomains.

Parameters

Solution (Solution) – the solution instance the scheme is working on

Returns

updating directly the LSValues in the solution instance

Return type

None

ComputeFlux(U, LSValues)
Emulates a flux computation as done in the iteration routine in order

to be used in the DeC subtimestep evaluation.

Parameters
  • U (numpy float array) – buffer containing the current state values at each dof (NbVars x NbDofs)

  • LS (numpy float array) – buffer containing the current FluidSpotters values at each dof(NbFluids x NbDofs)

Returns

fluxes in the format requires by the UpdateFSValues routine

Return type

fluxes (numpy float (multiD)array)

Iteration(Solution, fluxes, i, du=0, dt=1)

Main iteration of the scheme, implementing the continuous galerkin scheme.

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

  • 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 computed residuals (NbVars x NbDofs)

Return type

Resu (float numpy array)

UpdateFSVAttributes(Solution)

Wrapper routine to update the quantitites that depend on the UpdateFSValues fresh output

Parameters

Solution (Solution) – the solution instance the scheme is working on

Returns

updating directly the subvalues in the solution instance

Return type

None

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 InnerBFValues and FaceWiseBFValues in the scheme structure, containing

the values of each basis function at each quadrature point for each element and face.

Return type

None