LifeV
AbstractNumericalFlux< Mesh, SolverType > Class Template Referenceabstract

AbstractNumericalFlux Gives a common interface for hyperbolic's flux function. More...

#include <HyperbolicFluxNumerical.hpp>

+ Inheritance diagram for AbstractNumericalFlux< Mesh, SolverType >:
+ Collaboration diagram for AbstractNumericalFlux< Mesh, SolverType >:

Protected Attributes

vectorFunction_Type M_physicalFlux
 Physical flux function. More...
 
vectorFunction_Type M_firstDerivativePhysicalFlux
 First derivative, respect to unknown, of physical flux function. More...
 
Real M_CFLBrentToll
 Tollerance for the Brent algorithm for computing the CFL condition. More...
 
UInt M_CFLBrentMaxIter
 Maximum of iterations for the Brent algorithm for computing the CFL condition. More...
 
const FESpace< Mesh, MapEpetra > & M_fESpace
 Finite element space of the hyperbolic solver. More...
 
std::vector< const vectorPtr_Type *> M_fields
 Vector of pointers for the dependences of the permeability to external vector fields. More...
 

Public Types

typedef std::function< Vector(const Real &, const Real &, const Real &, const Real &, const std::vector< Real > &) > vectorFunction_Type
 
typedef std::function< Real(const Real &) > scalarFunction_Type
 
typedef SolverType::vector_type vector_Type
 
typedef std::shared_ptr< vector_TypevectorPtr_Type
 
typedef GetPot dataFile_Type
 
typedef KN< Realnormal_Type
 

Constructors and destructor

 AbstractNumericalFlux (const vectorFunction_Type &physicalFlux, const vectorFunction_Type &firstDerivativePhysicalFlux, const FESpace< Mesh, MapEpetra > &fESpace, const dataFile_Type &data, const std::string &section="numerical_flux/")
 Constructor for the class. More...
 
virtual ~AbstractNumericalFlux ()
 Virtual destructor. More...
 

Operators

virtual Real operator() (const Real &leftState, const Real &rightState, const normal_Type &normal, const UInt &iElem, const Real &t=0, const Real &x=0, const Real &y=0, const Real &z=0) const =0
 Computes the face contribution of the flux. More...
 

Set Methods

void setExternalField (const vectorPtr_Type &field)
 Add one external field. More...
 

Get Methods

vectorFunction_Type physicalFlux () const
 Return the physical flux. More...
 
vectorFunction_Type firstDerivativePhysicalFlux () const
 Return the first derivative, respect to the unknown, of the physical flux. More...
 
Real physicalFluxDotNormal (const normal_Type &normal, const UInt &iElem, const Real &t, const Real &x, const Real &y, const Real &z, const Real &unknown) const
 Evaluate the flux dot normal in a given point of a face. More...
 
Real firstDerivativePhysicalFluxDotNormal (const normal_Type &normal, const UInt &iElem, const Real &t, const Real &x, const Real &y, const Real &z, const Real &unknown) const
 Evaluate the first derivative of the flux dot normal in a given point of a face. More...
 
Real normInfinity (const Real &leftState, const Real &rightState, const normal_Type &normal, const UInt &iElem, const Real &t=0, const Real &x=0, const Real &y=0, const Real &z=0) const
 Computes the local infinity norm of the first derivative of the flux dot normal. More...
 

Protected Methods

scalarFunction_Type computeFunctionDotNormal (const vectorFunction_Type &function, const normal_Type &normal, const UInt &iElem, const Real &t, const Real &x, const Real &y, const Real &z, const Real &plusMinus) const
 Return a scalar function from a general vector function dot normal in a given point of a face. More...
 
scalarFunction_Type computeAbsFunctionDotNormal (const vectorFunction_Type &function, const normal_Type &normal, const UInt &iElem, const Real &t, const Real &x, const Real &y, const Real &z, const Real &plusMinus) const
 Return a scalar function from the absolute value of a general vector function dot normal in a given point of a face. More...
 

Detailed Description

template<typename Mesh, typename SolverType = LifeV::SolverAztecOO>
class LifeV::AbstractNumericalFlux< Mesh, SolverType >

AbstractNumericalFlux Gives a common interface for hyperbolic's flux function.

Author
Alessio Fumagalli aless.nosp@m.io.f.nosp@m.umaga.nosp@m.lli@.nosp@m.mail..nosp@m.poli.nosp@m.mi.it
Michel Kern miche.nosp@m.l.ke.nosp@m.rn@in.nosp@m.ria..nosp@m.fr

This class gives a common interface for hyperbolic's flux functions $ \mathbf{F}( u ) $, in particoular it computes

  1. the flux function dot product the outward unit normal $ \mathbf{n} $;
  2. the first derivative of the flux function, respect to $ u $, dot product the outward unit normal;
  3. the maximum value of $ \hat{\mathbf{F}} \cdot \mathbf{n} $ between to adjacent elements;
  4. the value of $ \Vert \mathbf{F}^\prime \cdot \mathbf{n}_{e, K} \Vert_{L^\infty(a_0, b_0) } $ between to adjacent elements.

The class can handle the dependeces of the flux function $ \mathbf{F} $ of external vector or scalar fields.

Note
In the implementation of the physical flux $ \mathbf{F} $ we suppose that the first parameter of the vector is the unknown. See the test case for an example.

Definition at line 137 of file HyperbolicFluxNumerical.hpp.

Member Typedef Documentation

◆ vectorFunction_Type

typedef std::function< Vector ( const Real&, const Real&, const Real&, const Real&, const std::vector<Real>& ) > vectorFunction_Type

Definition at line 148 of file HyperbolicFluxNumerical.hpp.

◆ scalarFunction_Type

typedef std::function< Real ( const Real& ) > scalarFunction_Type

Definition at line 150 of file HyperbolicFluxNumerical.hpp.

◆ vector_Type

◆ vectorPtr_Type

typedef std::shared_ptr< vector_Type > vectorPtr_Type

Definition at line 153 of file HyperbolicFluxNumerical.hpp.

◆ dataFile_Type

Definition at line 155 of file HyperbolicFluxNumerical.hpp.

◆ normal_Type

typedef KN< Real > normal_Type

Definition at line 156 of file HyperbolicFluxNumerical.hpp.

Constructor & Destructor Documentation

◆ AbstractNumericalFlux()

AbstractNumericalFlux ( const vectorFunction_Type physicalFlux,
const vectorFunction_Type firstDerivativePhysicalFlux,
const FESpace< Mesh, MapEpetra > &  fESpace,
const dataFile_Type data,
const std::string &  section = "numerical_flux/" 
)

Constructor for the class.

Parameters
physicalFluxPhysical flux for the problem.
firstDerivativePhysicalFluxFirst derivative, respect the the unknown, of the physical flux.
fESpaceFinite element space of the hyperbolic problem.
dataData for the problem.
sectionSection for read the data from GetPot file.

Definition at line 400 of file HyperbolicFluxNumerical.hpp.

+ Here is the caller graph for this function:

◆ ~AbstractNumericalFlux()

~AbstractNumericalFlux ( )
virtual

Virtual destructor.

Definition at line 418 of file HyperbolicFluxNumerical.hpp.

Member Function Documentation

◆ operator()()

virtual Real operator() ( const Real leftState,
const Real rightState,
const normal_Type normal,
const UInt iElem,
const Real t = 0,
const Real x = 0,
const Real y = 0,
const Real z = 0 
) const
pure virtual

Computes the face contribution of the flux.

Given a face $ e \in \partial K $ it evaluates

\[ \max \hat{\mathbf{F}} \cdot \mathbf{n} \]

between to elements sharing the face.

Parameters
leftStateLeft value of the unknown respect to the face.
rightValueRight value of the unknown respect to the face.
normalNormal of the face.
iElemThe ID of the current element in the mesh.
tCurrent time.
xAbscissa.
yOrdinate.
zQuota.
Note
We assume left and right side of $ e $ is given by: the normal direction goes from the left side to the right side.

Implemented in GodunovNumericalFlux< Mesh, SolverType >.

◆ setExternalField()

void setExternalField ( const vectorPtr_Type field)
inline

Add one external field.

Add one extra field for the dependece from $ \mathbf{F} $.

Parameters
fieldThe filed to be added.

Definition at line 223 of file HyperbolicFluxNumerical.hpp.

◆ physicalFlux()

vectorFunction_Type physicalFlux ( ) const
inline

Return the physical flux.

Returns
The physical flux.

Definition at line 237 of file HyperbolicFluxNumerical.hpp.

◆ firstDerivativePhysicalFlux()

vectorFunction_Type firstDerivativePhysicalFlux ( ) const
inline

Return the first derivative, respect to the unknown, of the physical flux.

Returns
The first derivative of the physical flux.

Definition at line 246 of file HyperbolicFluxNumerical.hpp.

◆ physicalFluxDotNormal()

Real physicalFluxDotNormal ( const normal_Type normal,
const UInt iElem,
const Real t,
const Real x,
const Real y,
const Real z,
const Real unknown 
) const
inline

Evaluate the flux dot normal in a given point of a face.

Parameters
normalThe normal of the face.
iElemThe ID of the current element in the mesh.
tCurrent time.
xAbscissa.
yOrdinate.
zQuota.
unknownThe value of the unknown.
Returns
The value of $ \mathbf{ \hat{ F } } \cdot \mathbf{ n }$ in the point $ (x, y, z, t, u) $.

Definition at line 262 of file HyperbolicFluxNumerical.hpp.

◆ firstDerivativePhysicalFluxDotNormal()

Real firstDerivativePhysicalFluxDotNormal ( const normal_Type normal,
const UInt iElem,
const Real t,
const Real x,
const Real y,
const Real z,
const Real unknown 
) const
inline

Evaluate the first derivative of the flux dot normal in a given point of a face.

Parameters
normalThe normal of the face.
iElemThe ID of the current element in the mesh.
tCurrent time.
xAbscissa.
yOrdinate.
zQuota.
unknownThe value of the unknown.
Returns
The value of $ \mathbf{ \hat{ F^\prime } } \cdot \mathbf{ n }$ in the point $ (x, y, z, t, u) $.

Definition at line 284 of file HyperbolicFluxNumerical.hpp.

◆ normInfinity()

Real normInfinity ( const Real leftState,
const Real rightState,
const normal_Type normal,
const UInt iElem,
const Real t = 0,
const Real x = 0,
const Real y = 0,
const Real z = 0 
) const

Computes the local infinity norm of the first derivative of the flux dot normal.

Given a face $ e \in \partial K $ it evaluates

\[ \displaystyle \max_{e \in \partial K^- \cap \partial K^+} \vert \hat{\mathbf{F^\prime}} \cdot \mathbf{n} \vert \]

between to elements sharing the face.

Parameters
leftStateLeft value of the unknown respect to the face.
rightValueRight value of the unknown respect to the face.
normalNormal of the face.
iElemThe ID of the current element in the mesh.
tCurrent time.
xAbscissa.
yOrdinate.
zQuota.
Note
We assume left and right side of $ e $ is given by: the normal direction goes from the left side to the right side.

Definition at line 431 of file HyperbolicFluxNumerical.hpp.

◆ computeFunctionDotNormal()

std::function< Real(const Real &) > computeFunctionDotNormal ( const vectorFunction_Type function,
const normal_Type normal,
const UInt iElem,
const Real t,
const Real x,
const Real y,
const Real z,
const Real plusMinus 
) const
protected

Return a scalar function from a general vector function dot normal in a given point of a face.

Parameters
functionFunction to be reduced.
normalThe normal of the face.
iElemThe ID of the current element in the mesh.
tCurrent time.
xAbscissa.
yOrdinate.
zQuota.
plusMinus
Returns
The function $ g(u) = \pm \mathbf{f}(x,y,z,t,u) \cdot \mathbf{n} $.

Definition at line 469 of file HyperbolicFluxNumerical.hpp.

+ Here is the caller graph for this function:

◆ computeAbsFunctionDotNormal()

scalarFunction_Type computeAbsFunctionDotNormal ( const vectorFunction_Type function,
const normal_Type normal,
const UInt iElem,
const Real t,
const Real x,
const Real y,
const Real z,
const Real plusMinus 
) const
protected

Return a scalar function from the absolute value of a general vector function dot normal in a given point of a face.

Parameters
functionFunction to be reduced.
normalThe normal of the face.
iElemThe ID of the current element in the mesh.
tCurrent time.
xAbscissa.
yOrdinate.
zQuota.
plusMinus
Returns
The function $ g(u) = \pm \vert \mathbf{f}(x,y,z,t,u) \cdot \mathbf{n} \vert $.

Field Documentation

◆ M_physicalFlux

vectorFunction_Type M_physicalFlux
protected

Physical flux function.

Definition at line 374 of file HyperbolicFluxNumerical.hpp.

◆ M_firstDerivativePhysicalFlux

vectorFunction_Type M_firstDerivativePhysicalFlux
protected

First derivative, respect to unknown, of physical flux function.

Definition at line 377 of file HyperbolicFluxNumerical.hpp.

◆ M_CFLBrentToll

Real M_CFLBrentToll
protected

Tollerance for the Brent algorithm for computing the CFL condition.

Definition at line 380 of file HyperbolicFluxNumerical.hpp.

◆ M_CFLBrentMaxIter

UInt M_CFLBrentMaxIter
protected

Maximum of iterations for the Brent algorithm for computing the CFL condition.

Definition at line 383 of file HyperbolicFluxNumerical.hpp.

◆ M_fESpace

const FESpace<Mesh, MapEpetra>& M_fESpace
protected

Finite element space of the hyperbolic solver.

Definition at line 386 of file HyperbolicFluxNumerical.hpp.

◆ M_fields

std::vector< const vectorPtr_Type* > M_fields
protected

Vector of pointers for the dependences of the permeability to external vector fields.

Definition at line 389 of file HyperbolicFluxNumerical.hpp.


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