LifeV
OneDFSIFunctionSolverDefinedCompatibility Class Reference

OneDFSIFunctionSolverDefinedCompatibility - Class which implements Compatibility boundary conditions for the 1D segment. More...

#include <OneDFSIFunctionSolverDefined.hpp>

+ Inheritance diagram for OneDFSIFunctionSolverDefinedCompatibility:
+ Collaboration diagram for OneDFSIFunctionSolverDefinedCompatibility:

Protected Attributes

UInt M_bcElement
 ID of the boundary edge. More...
 
UInt M_bcInternalNode
 Dof of the internal node adjacent to the boundary. More...
 
container2D_Type M_eigenvalues
 Eigen values of the jacobian diffFlux (= dF/dU = H) More...
 
container2D_Type M_deltaEigenvalues
 
container2D_Type M_leftEigenvector1
 Left eigen vectors for the two eigen values. More...
 
container2D_Type M_leftEigenvector2
 
container2D_Type M_deltaLeftEigenvector1
 
container2D_Type M_deltaLeftEigenvector2
 
- Protected Attributes inherited from OneDFSIFunctionSolverDefinedRiemann
container2D_Type M_bcU
 Value of U at the boundary. More...
 
container2D_Type M_bcW
 Value of W at the boundary. More...
 
- Protected Attributes inherited from OneDFSIFunctionSolverDefined
fluxPtr_Type M_fluxPtr
 
sourcePtr_Type M_sourcePtr
 
solutionPtr_Type M_solutionPtr
 
UInt M_bcNode
 
bcSide_Type M_bcSide
 
bcType_Type M_bcType
 

Type definitions

typedef OneDFSIFunctionSolverDefinedRiemann super
 
typedef super::fluxPtr_Type fluxPtr_Type
 
typedef super::sourcePtr_Type sourcePtr_Type
 
typedef super::solutionPtr_Type solutionPtr_Type
 
typedef super::mesh_Type mesh_Type
 
typedef super::container2D_Type container2D_Type
 

Constructors & Destructor

 OneDFSIFunctionSolverDefinedCompatibility (const bcSide_Type &bcSide, const bcType_Type &bcType)
 Constructor. More...
 
 OneDFSIFunctionSolverDefinedCompatibility (const OneDFSIFunctionSolverDefinedCompatibility &bcFunctionCompatibility)
 Copy constructor. More...
 
virtual ~OneDFSIFunctionSolverDefinedCompatibility ()
 Destructor. More...
 

Methods

virtual Real operator() (const Real &, const Real &timeStep)
 Operator() More...
 

Protected Methods

void setupNode ()
 Automatically identify the boundary node. More...
 
Real computeRHS (const Real &timeStep)
 Compute the rhs. More...
 
void computeEigenValuesVectors ()
 Compute the current eigenvalues and eigenvectors. More...
 
Real evaluateRHS (const Real &eigenvalue, const container2D_Type &eigenvector, const container2D_Type &deltaEigenvector, const Real &timeStep)
 Compute the rhs. More...
 
Real computeCFL (const Real &eigenvalue, const Real &timeStep) const
 Compute the current CFL. More...
 
Real scalarProduct (const container2D_Type &vector1, const container2D_Type &vector2)
 Scalar product between 2 2D vectors. More...
 

Additional Inherited Members

- Public Types inherited from OneDFSIFunctionSolverDefinedRiemann
typedef OneDFSIFunctionSolverDefined super
 
typedef super::container2D_Type container2D_Type
 
- Public Types inherited from OneDFSIFunctionSolverDefined
typedef OneDFSIFunction bcFunction_Type
 
typedef std::shared_ptr< bcFunction_TypebcFunctionPtr_Type
 
typedef OneDFSIFlux flux_Type
 
typedef std::shared_ptr< flux_TypefluxPtr_Type
 
typedef OneDFSISource source_Type
 
typedef std::shared_ptr< source_TypesourcePtr_Type
 
typedef OneDFSIData data_Type
 
typedef data_Type::mesh_Type mesh_Type
 
typedef data_Type::container2D_Type container2D_Type
 
typedef SolverAmesos linearSolver_Type
 
typedef linearSolver_Type::vector_type vector_Type
 
typedef std::shared_ptr< vector_TypevectorPtr_Type
 
typedef std::array< vectorPtr_Type, 2 > vectorPtrContainer_Type
 
typedef linearSolver_Type::matrix_type matrix_Type
 
typedef std::map< std::string, vectorPtr_Typesolution_Type
 
typedef std::shared_ptr< solution_TypesolutionPtr_Type
 
typedef OneDFSI::bcLine_Type bcLine_Type
 
typedef OneDFSI::bcSide_Type bcSide_Type
 
typedef OneDFSI::bcType_Type bcType_Type
 
- Public Member Functions inherited from OneDFSIFunctionSolverDefinedRiemann
 OneDFSIFunctionSolverDefinedRiemann (const bcSide_Type &bcSide, const bcType_Type &bcType)
 Constructor. More...
 
 OneDFSIFunctionSolverDefinedRiemann (const OneDFSIFunctionSolverDefinedRiemann &bcFunctionRiemann)
 Copy constructor. More...
 
virtual ~OneDFSIFunctionSolverDefinedRiemann ()
 Destructor. More...
 
- Public Member Functions inherited from OneDFSIFunctionSolverDefined
 OneDFSIFunctionSolverDefined (const bcSide_Type &bcSide, const bcType_Type &bcType)
 Constructor. More...
 
 OneDFSIFunctionSolverDefined (const OneDFSIFunctionSolverDefined &bcFunctionDefault)
 Copy constructor. More...
 
virtual ~OneDFSIFunctionSolverDefined ()
 Destructor. More...
 
void setFluxSource (const fluxPtr_Type &fluxPtr, const sourcePtr_Type &sourcePtr)
 Set the flux and the source classes for the problem. More...
 
void setSolution (const solutionPtr_Type &solutionPtr)
 Set the solution of the problem. More...
 
- Protected Member Functions inherited from OneDFSIFunctionSolverDefinedRiemann
void updateBCVariables ()
 Update the boundary condition variables. More...
 
- Protected Member Functions inherited from OneDFSIFunctionSolverDefined

Detailed Description

OneDFSIFunctionSolverDefinedCompatibility - Class which implements Compatibility boundary conditions for the 1D segment.

Author
Lucia Mirabella, Tiziano Passerini, Cristiano Malossi

The compatibility equations are derived using the pseudo-characteristic teory:

\[ \mathbf L(\mathbf U^{n+1}-\mathbf U^{n}_\star + \mathbf U^0 - \mathbf U^0_\star) = \Delta t \left( \mathbf \Lambda \displaystyle\frac{\partial \mathbf L}{\partial z}(\mathbf U^n_\star - \mathbf U^0_\star) - \mathbf L \mathbf B(\mathbf U^n_\star) + \mathbf L \mathbf B(\mathbf U^0_\star) \right). \]

Definition at line 269 of file OneDFSIFunctionSolverDefined.hpp.

Member Typedef Documentation

◆ super

◆ fluxPtr_Type

◆ sourcePtr_Type

◆ solutionPtr_Type

◆ mesh_Type

Definition at line 282 of file OneDFSIFunctionSolverDefined.hpp.

◆ container2D_Type

Constructor & Destructor Documentation

◆ OneDFSIFunctionSolverDefinedCompatibility() [1/2]

OneDFSIFunctionSolverDefinedCompatibility ( const bcSide_Type bcSide,
const bcType_Type bcType 
)
explicit

Constructor.

Parameters
bcLinethe line of the boundary condition (first or second).
bcTypethe type of the boundary condition ( $Q$, $A$, $P$, $S$, $W_1$, $W_2$).

Definition at line 147 of file OneDFSIFunctionSolverDefined.cpp.

+ Here is the caller graph for this function:

◆ OneDFSIFunctionSolverDefinedCompatibility() [2/2]

Copy constructor.

Parameters
bcFunctionCompatibilityOneDFSIFunctionSolverDefinedCompatibility

Definition at line 160 of file OneDFSIFunctionSolverDefined.cpp.

+ Here is the caller graph for this function:

◆ ~OneDFSIFunctionSolverDefinedCompatibility()

virtual ~OneDFSIFunctionSolverDefinedCompatibility ( )
inlinevirtual

Destructor.

Definition at line 305 of file OneDFSIFunctionSolverDefined.hpp.

Member Function Documentation

◆ operator()()

virtual Real operator() ( const Real ,
const Real timeStep 
)
inlinevirtual

Operator()

Evaluate the function.

Parameters
timethe current time.
timeStepthe time step.
Returns
the value of the function.

Reimplemented from OneDFSIFunctionSolverDefinedRiemann.

Reimplemented in OneDFSIFunctionSolverDefinedWindkessel3, and OneDFSIFunctionSolverDefinedAbsorbing.

Definition at line 321 of file OneDFSIFunctionSolverDefined.hpp.

◆ setupNode()

void setupNode ( )
protectedvirtual

Automatically identify the boundary node.

Reimplemented from OneDFSIFunctionSolverDefined.

Definition at line 177 of file OneDFSIFunctionSolverDefined.cpp.

◆ computeRHS()

Real computeRHS ( const Real timeStep)
protected

Compute the rhs.

Parameters
timeStepthe time step.
Returns
rhs of the problem.

Definition at line 201 of file OneDFSIFunctionSolverDefined.cpp.

+ Here is the caller graph for this function:

◆ computeEigenValuesVectors()

void computeEigenValuesVectors ( )
protected

Compute the current eigenvalues and eigenvectors.

Definition at line 223 of file OneDFSIFunctionSolverDefined.cpp.

+ Here is the caller graph for this function:

◆ evaluateRHS()

Real evaluateRHS ( const Real eigenvalue,
const container2D_Type eigenvector,
const container2D_Type deltaEigenvector,
const Real timeStep 
)
protected

Compute the rhs.

Parameters
eigenvalueeigenvalue
eigenvectoreigenvector
deltaEigenvectorderivative of the eigenvector
timeStepthe time step.
Returns
rhs of the problem

Definition at line 235 of file OneDFSIFunctionSolverDefined.cpp.

◆ computeCFL()

Real computeCFL ( const Real eigenvalue,
const Real timeStep 
) const
protected

Compute the current CFL.

Parameters
eigenvalueeigenvalue
timeStepthe time step.
Returns
CFL

Definition at line 287 of file OneDFSIFunctionSolverDefined.cpp.

+ Here is the caller graph for this function:

◆ scalarProduct()

Real scalarProduct ( const container2D_Type vector1,
const container2D_Type vector2 
)
inlineprotected

Scalar product between 2 2D vectors.

vector1 first vector vector2 second vector

Returns
scalar product

Definition at line 371 of file OneDFSIFunctionSolverDefined.hpp.

+ Here is the caller graph for this function:

Field Documentation

◆ M_bcElement

UInt M_bcElement
protected

ID of the boundary edge.

Definition at line 379 of file OneDFSIFunctionSolverDefined.hpp.

◆ M_bcInternalNode

UInt M_bcInternalNode
protected

Dof of the internal node adjacent to the boundary.

Definition at line 381 of file OneDFSIFunctionSolverDefined.hpp.

◆ M_eigenvalues

container2D_Type M_eigenvalues
protected

Eigen values of the jacobian diffFlux (= dF/dU = H)

Definition at line 384 of file OneDFSIFunctionSolverDefined.hpp.

◆ M_deltaEigenvalues

container2D_Type M_deltaEigenvalues
protected

Definition at line 385 of file OneDFSIFunctionSolverDefined.hpp.

◆ M_leftEigenvector1

container2D_Type M_leftEigenvector1
protected

Left eigen vectors for the two eigen values.

Definition at line 388 of file OneDFSIFunctionSolverDefined.hpp.

◆ M_leftEigenvector2

container2D_Type M_leftEigenvector2
protected

Definition at line 389 of file OneDFSIFunctionSolverDefined.hpp.

◆ M_deltaLeftEigenvector1

container2D_Type M_deltaLeftEigenvector1
protected

Definition at line 390 of file OneDFSIFunctionSolverDefined.hpp.

◆ M_deltaLeftEigenvector2

container2D_Type M_deltaLeftEigenvector2
protected

Definition at line 391 of file OneDFSIFunctionSolverDefined.hpp.


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