LifeV
|
AbstractNumericalFlux Gives a common interface for hyperbolic's flux function. More...
#include <HyperbolicFluxNumerical.hpp>
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_Type > | vectorPtr_Type |
typedef GetPot | dataFile_Type |
typedef KN< Real > | normal_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 §ion="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... | |
AbstractNumericalFlux Gives a common interface for hyperbolic's flux function.
This class gives a common interface for hyperbolic's flux functions , in particoular it computes
The class can handle the dependeces of the flux function of external vector or scalar fields.
Definition at line 137 of file HyperbolicFluxNumerical.hpp.
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.
typedef std::function< Real ( const Real& ) > scalarFunction_Type |
Definition at line 150 of file HyperbolicFluxNumerical.hpp.
typedef SolverType::vector_type vector_Type |
Definition at line 152 of file HyperbolicFluxNumerical.hpp.
typedef std::shared_ptr< vector_Type > vectorPtr_Type |
Definition at line 153 of file HyperbolicFluxNumerical.hpp.
typedef GetPot dataFile_Type |
Definition at line 155 of file HyperbolicFluxNumerical.hpp.
typedef KN< Real > normal_Type |
Definition at line 156 of file HyperbolicFluxNumerical.hpp.
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.
physicalFlux | Physical flux for the problem. |
firstDerivativePhysicalFlux | First derivative, respect the the unknown, of the physical flux. |
fESpace | Finite element space of the hyperbolic problem. |
data | Data for the problem. |
section | Section for read the data from GetPot file. |
Definition at line 400 of file HyperbolicFluxNumerical.hpp.
|
virtual |
Virtual destructor.
Definition at line 418 of file HyperbolicFluxNumerical.hpp.
|
pure virtual |
Computes the face contribution of the flux.
Given a face it evaluates
between to elements sharing the face.
leftState | Left value of the unknown respect to the face. |
rightValue | Right value of the unknown respect to the face. |
normal | Normal of the face. |
iElem | The ID of the current element in the mesh. |
t | Current time. |
x | Abscissa. |
y | Ordinate. |
z | Quota. |
Implemented in GodunovNumericalFlux< Mesh, SolverType >.
|
inline |
Add one external field.
Add one extra field for the dependece from .
field | The filed to be added. |
Definition at line 223 of file HyperbolicFluxNumerical.hpp.
|
inline |
Return the physical flux.
Definition at line 237 of file HyperbolicFluxNumerical.hpp.
|
inline |
Return the first derivative, respect to the unknown, of the physical flux.
Definition at line 246 of file HyperbolicFluxNumerical.hpp.
|
inline |
Evaluate the flux dot normal in a given point of a face.
normal | The normal of the face. |
iElem | The ID of the current element in the mesh. |
t | Current time. |
x | Abscissa. |
y | Ordinate. |
z | Quota. |
unknown | The value of the unknown. |
Definition at line 262 of file HyperbolicFluxNumerical.hpp.
|
inline |
Evaluate the first derivative of the flux dot normal in a given point of a face.
normal | The normal of the face. |
iElem | The ID of the current element in the mesh. |
t | Current time. |
x | Abscissa. |
y | Ordinate. |
z | Quota. |
unknown | The value of the unknown. |
Definition at line 284 of file HyperbolicFluxNumerical.hpp.
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 it evaluates
between to elements sharing the face.
leftState | Left value of the unknown respect to the face. |
rightValue | Right value of the unknown respect to the face. |
normal | Normal of the face. |
iElem | The ID of the current element in the mesh. |
t | Current time. |
x | Abscissa. |
y | Ordinate. |
z | Quota. |
Definition at line 431 of file HyperbolicFluxNumerical.hpp.
|
protected |
Return a scalar function from a general vector function dot normal in a given point of a face.
function | Function to be reduced. |
normal | The normal of the face. |
iElem | The ID of the current element in the mesh. |
t | Current time. |
x | Abscissa. |
y | Ordinate. |
z | Quota. |
plusMinus |
Definition at line 469 of file HyperbolicFluxNumerical.hpp.
|
protected |
Return a scalar function from the absolute value of a general vector function dot normal in a given point of a face.
function | Function to be reduced. |
normal | The normal of the face. |
iElem | The ID of the current element in the mesh. |
t | Current time. |
x | Abscissa. |
y | Ordinate. |
z | Quota. |
plusMinus |
|
protected |
Physical flux function.
Definition at line 374 of file HyperbolicFluxNumerical.hpp.
|
protected |
First derivative, respect to unknown, of physical flux function.
Definition at line 377 of file HyperbolicFluxNumerical.hpp.
|
protected |
Tollerance for the Brent algorithm for computing the CFL condition.
Definition at line 380 of file HyperbolicFluxNumerical.hpp.
|
protected |
Maximum of iterations for the Brent algorithm for computing the CFL condition.
Definition at line 383 of file HyperbolicFluxNumerical.hpp.
Finite element space of the hyperbolic solver.
Definition at line 386 of file HyperbolicFluxNumerical.hpp.
|
protected |
Vector of pointers for the dependences of the permeability to external vector fields.
Definition at line 389 of file HyperbolicFluxNumerical.hpp.