LifeV
OneDFSISolver Class Reference

OneDFSISolver - Solver class for the 1D model. More...

#include <OneDFSISolver.hpp>

+ Collaboration diagram for OneDFSISolver:

Private Attributes

physicsPtr_Type M_physicsPtr
 L2 Projection of the second derivative of Q over P1 space. More...
 
fluxPtr_Type M_fluxPtr
 
sourcePtr_Type M_sourcePtr
 
feSpacePtr_Type M_feSpacePtr
 
commPtr_Type M_commPtr
 
Displayer M_displayer
 
std::shared_ptr< MatrixElementalM_elementalMassMatrixPtr
 element mass matrix More...
 
std::shared_ptr< MatrixElementalM_elementalStiffnessMatrixPtr
 element stiffness matrix More...
 
std::shared_ptr< MatrixElementalM_elementalGradientMatrixPtr
 element gradient matrix More...
 
std::shared_ptr< MatrixElementalM_elementalDivergenceMatrixPtr
 element divergence matrix More...
 
vectorPtrContainer_Type M_rhs
 Right hand sides of the linear system i: "mass * M_Ui = M_rhsi". More...
 
vectorPtrContainer_Type M_residual
 Residual of the linear system. More...
 
vectorPtrContainer_Type M_fluxVector
 Flux F(U) (in P1) More...
 
vectorPtrContainer_Type M_sourceVector
 Source term S (in P1) More...
 
scalarVectorContainer_Type M_dFdUVector
 diffFlux = dF(U)/dU (in P0) More...
 
scalarVectorContainer_Type M_dSdUVector
 diffSrc = dSource(U)/dU (in P0) More...
 
matrixPtr_Type M_homogeneousMassMatrixPtr
 tridiagonal mass matrix More...
 
matrixPtr_Type M_homogeneousGradientMatrixPtr
 tridiagonal gradient matrix More...
 
matrixPtrContainer_Type M_dSdUMassMatrixPtr
 tridiagonal mass matrices multiplied by diffSrcij More...
 
matrixPtrContainer_Type M_dFdUStiffnessMatrixPtr
 tridiagonal stiffness matrices multiplied by diffFluxij More...
 
matrixPtrContainer_Type M_dFdUGradientMatrixPtr
 tridiagonal gradient matrices multiplied by diffFluxij More...
 
matrixPtrContainer_Type M_dSdUDivergenceMatrixPtr
 tridiagonal divergence matrices multiplied by diffSrcij More...
 
linearSolverPtr_Type M_linearSolverPtr
 The linear solver. More...
 
linearSolverPtr_Type M_linearViscoelasticSolverPtr
 

Typedef & Enumerator

typedef OneDFSIPhysics physics_Type
 
typedef std::shared_ptr< physics_TypephysicsPtr_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 data_Type::scalarVector_Type scalarVector_Type
 
typedef std::array< scalarVector_Type, 4 > scalarVectorContainer_Type
 
typedef FESpace< mesh_Type, MapEpetrafeSpace_Type
 
typedef std::shared_ptr< feSpace_TypefeSpacePtr_Type
 
typedef Epetra_Comm comm_Type
 
typedef std::shared_ptr< comm_TypecommPtr_Type
 
typedef SolverAmesos linearSolver_Type
 
typedef std::shared_ptr< linearSolver_TypelinearSolverPtr_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::shared_ptr< matrix_TypematrixPtr_Type
 
typedef std::array< matrixPtr_Type, 4 > matrixPtrContainer_Type
 
typedef std::map< std::string, vectorPtr_Typesolution_Type
 
typedef std::shared_ptr< solution_TypesolutionPtr_Type
 
typedef solution_Type::const_iterator solutionConstIterator_Type
 
typedef OneDFSI::bcLine_Type bcLine_Type
 
typedef OneDFSI::bcSide_Type bcSide_Type
 
typedef OneDFSI::bcType_Type bcType_Type
 

Constructors & Destructor

 OneDFSISolver ()
 Empty Constructor. More...
 
virtual ~OneDFSISolver ()
 Destructor. More...
 

Methods

void buildConstantMatrices ()
 Build constant matrices (mass and grad) More...
 
void setupSolution (solution_Type &solution)
 Setup the solution using the default FESpace map. More...
 
void setupSolution (solution_Type &solution, const MapEpetra &map, const bool &onlyMainQuantities=false)
 Setup the solution using user defined FESpace map. More...
 
void initialize (solution_Type &solution)
 Initialize all the variables of the solution to a reference condition with $Q=0$, $A=A^0$, and $P=P_\mathrm{ext}$. More...
 
void computeW1W2 (solution_Type &solution)
 Update the Riemann variables. More...
 
void computePressure (solution_Type &solution, const Real &timeStep)
 Update the pressure. More...
 
void computeAreaRatio (solution_Type &solution)
 Update the ratio between $A$ and $A^0$. More...
 
void computeArea (solution_Type &solution)
 Compute A from the area ratio: $\displaystyle\frac{A}{A^0}-1$. More...
 
void updateRHS (const solution_Type &solution, const Real &timeStep)
 Compute the right hand side. More...
 
void iterate (OneDFSIBCHandler &bcH, solution_Type &solution, const Real &time, const Real &timeStep)
 Update convective term and BC. Then solve the linearized system. More...
 
vector_Type viscoelasticFlowRateCorrection (const vector_Type &newArea, const vector_Type &newElasticFlowRate, const vector_Type &oldViscoelasticFlowRate, const Real &timeStep, OneDFSIBCHandler &bcHandler, const bool &updateSystemMatrix=true)
 Apply the viscoelastic flow rate correction. More...
 
Real computeCFL (const solution_Type &solution, const Real &timeStep) const
 CFL computation (correct for constant mesh) More...
 
void resetOutput (const solution_Type &solution)
 Reset the output files. More...
 
void postProcess (const solution_Type &solution, const Real &time)
 Save results on output files. More...
 

Set Methods

void setProblem (const physicsPtr_Type &physicsPtr, const fluxPtr_Type &fluxPtr, const sourcePtr_Type &sourcePtr)
 Set problem classes. More...
 
void setCommunicator (const commPtr_Type &commPtr)
 Set the communicator. More...
 
void setFESpace (const feSpacePtr_Type &feSpacePtr)
 Set the FEspace. More...
 
void setLinearSolver (const linearSolverPtr_Type &linearSolverPtr)
 Set the linear solver. More...
 
void setLinearViscoelasticSolver (const linearSolverPtr_Type &linearViscoelasticSolverPtr)
 Set the viscoelastic linear solver. More...
 

Get Methods

const physicsPtr_Typephysics () const
 Get the physics class. More...
 
const fluxPtr_Typeflux () const
 Get the flux class. More...
 
const sourcePtr_Typesource () const
 Get the source class. More...
 
UInt boundaryDOF (const bcSide_Type &bcSide) const
 Return the ID of the boundary node given a side. More...
 
Real boundaryValue (const solution_Type &solution, const bcType_Type &bcType, const bcSide_Type &bcSide) const
 Return the value of a quantity ( $P$, $A$, $Q$, $W_1$, $W_2$) on a specified boundary. More...
 
void boundaryEigenValuesEigenVectors (const bcSide_Type &bcSide, const solution_Type &solution, container2D_Type &eigenvalues, container2D_Type &leftEigenvector1, container2D_Type &leftEigenvector2)
 Return the value of the eigenvalues and eigenvectors on a specified boundary. More...
 
const vectorPtrContainer_Typeresidual () const
 Get the residual container. More...
 
const matrixPtr_TypemassMatrix () const
 Get the system matrix without BC. More...
 

Private Methods

void updateFlux (const solution_Type &solution)
 Update the P1 flux vector from U: M_fluxi = F_h(Un) i=1,2 (works only for P1Seg elements) More...
 
void updatedFdU (const solution_Type &solution)
 Call _updateFlux and update the P0 derivative of flux vector from U: More...
 
void updateSource (const solution_Type &solution)
 Update the P1 source vector from U: M_sourcei = S_h(Un) i=1,2 (works only for P1Seg elements) More...
 
void updatedSdU (const solution_Type &solution)
 Call _updateSource and update the P0 derivative of source vector from U: More...
 
void updateMatrices ()
 Update the matrices. More...
 
void updateElementalMatrices (const Real &dFdU, const Real &dSdU)
 Update the element matrices with the current element. More...
 
void matrixAssemble (const UInt &ii, const UInt &jj)
 Assemble the matrices. More...
 
void applyDirichletBCToMatrix (matrix_Type &matrix)
 Update the matrices to take into account Dirichlet BC. More...
 
vector_Type inertialFlowRateCorrection (const vector_Type &)
 Apply the inertial Flux correction: More...
 
vector_Type longitudinalFlowRateCorrection ()
 Apply the longitudinal Flux correction: More...
 

Unimplemented Methods

 OneDFSISolver (const OneDFSISolver &solver)
 
OneDFSISolveroperator= (const OneDFSISolver &solver)
 

Detailed Description

OneDFSISolver - Solver class for the 1D model.

Author
Vincent Martin, Tiziano Passerini, Lucia Mirabella, Gilles Fourestey, Cristiano Malossi
See also
Equations and networks of 1-D models [7]
Geometrical multiscale coupling of 1-D models [12] [13] [3]

EQUATIONS:
The conservative form of the generic hyperbolic problem is

\[ \frac{\partial \mathbf U}{\partial t} + \frac{\partial \mathbf F(\mathbf U)}{\partial z} + \mathbf S(\mathbf U) = 0, \]

where $\mathbf U$ are the conservative variables, $\mathbf F$ the corresponding fluxes, and $\mathbf S$ represents the source terms.

NUMERICAL DISCRETIZATION:
We discretize the problem by a second-order Taylor–Galerkin scheme. Let us consider the time interval $[t^n, t^{n+1}]$, for $n=0,1,2,\dots$, with $t^n=n \Delta t$, $\Delta t$ being the time step. Given $ \mathbf U_h^n$ we compute $\mathbf U_h^{n+1}$ by using the following scheme

\[ (\mathbf U^{n+1}_h,\varphi_h) = (\mathbf U^n_h,\varphi_h) + \Delta t \left( \mathbf F(\mathbf U^n_h)- \displaystyle\frac{\Delta t}{2} \displaystyle\frac{\partial \mathbf F(\mathbf U^n_h)}{\partial \mathbf U}\left(\mathbf S(\mathbf U^n_h)+ \displaystyle\frac{\partial \mathbf F(\mathbf U^n_h)}{\partial z}\right), \displaystyle\frac{\partial \varphi_h}{\partial z}\right) + \Delta t \left( -\mathbf S(\mathbf U^n_h)+\displaystyle\frac{\Delta t}{2} \displaystyle\frac{\partial \mathbf S(\mathbf U^n_h)}{\partial \mathbf U}\left( \mathbf S(\mathbf U^n_h)+ \displaystyle\frac{\partial \mathbf F(\mathbf U^n_h)}{\partial z}\right), \varphi_h\right). \]

IMPLEMENTATION:
The implementation of the Taylor-Galerkin scheme is the following:

(Un+1, phi) = //! massFactor^{-1} * Un+1
(Un, phi) //! mass * U
+dt * ( Fh(Un), dphi/dz ) //! grad * F(U)
-dt^2/2 * (diffFh(Un) Sh(Un), dphi/dz ) //! gradDiffFlux(U) * S(U)
+dt^2/2 * (diffSh(Un) dFh/dz(Un), phi ) //! divDiffSrc(U) * F(U)
-dt^2/2 * (diffFh(Un) dFh/dz(Un), dphi/dz ) //!stiffDiffFlux(U) * F(U)
-dt * ( Sh(Un), phi ) //! mass * S(U)
+dt^2/2 * (diffSh(Un) Sh(Un), phi ) //! massDiffSrc(U) * S(U)

Let's define:

  1. ()_{i in nodes} is the basis of P1 (the "hat" functions)
  2. (1_{i+1/2})_{i+1/2 in elements} is the basis of P0 (constant per element). The vertices of the element "i+1/2" are the nodes "i" and "i+1".

Then:

  1. Uh is in P1 : U = sum_{i in nodes} U_i phi_i
  2. Fh(U) is in P1 : F(U) = sum_{i in nodes} F(U_i) phi_i
  3. diffFh(U) is in P0 : diffFlux(U) = sum_{i+1/2 in elements} 1/2 { dF/dU(U_i) + dF/dU(U_i+1) } 1_{i+1/2} (means of the two extremal values of the cell)
  4. dF/dz(U) = sum_{i in nodes} F(U_i) d(phi_i)/dz
  5. Sh(U) is in P1 : S(U) = sum_{i in nodes} S(U_i) phi_i
  6. diffSh(U) is in P0 : diffSrc(U) = sum_{i+1/2 in elements} 1/2 { dS/dU(U_i) + dS/dU(U_i+1) } 1_{i+1/2} (means of the two extremal values of the cell)

DEVELOPMENT NOTES:
The option taken here is to define the different tridiagonal matrix operators (div, grad, mass, stiff) and reconstruct them at each time step (as they depend on diffFlux and diffSrc). They are thus rebuilt at the element level and reassembled. Afterwards, there remains to do only some tridiagonal matrix vector products to obtain the right hand side. This procedure might appear a bit memory consuming (there are 18 tridiagonal matrices stored), but it has the advantage of being very clear. If it is too costly, it should be quite easy to improve it.

Definition at line 149 of file OneDFSISolver.hpp.

Member Typedef Documentation

◆ physics_Type

Definition at line 156 of file OneDFSISolver.hpp.

◆ physicsPtr_Type

typedef std::shared_ptr< physics_Type > physicsPtr_Type

Definition at line 157 of file OneDFSISolver.hpp.

◆ flux_Type

Definition at line 159 of file OneDFSISolver.hpp.

◆ fluxPtr_Type

typedef std::shared_ptr< flux_Type > fluxPtr_Type

Definition at line 160 of file OneDFSISolver.hpp.

◆ source_Type

Definition at line 162 of file OneDFSISolver.hpp.

◆ sourcePtr_Type

typedef std::shared_ptr< source_Type > sourcePtr_Type

Definition at line 163 of file OneDFSISolver.hpp.

◆ data_Type

Definition at line 165 of file OneDFSISolver.hpp.

◆ mesh_Type

Definition at line 166 of file OneDFSISolver.hpp.

◆ container2D_Type

◆ scalarVector_Type

◆ scalarVectorContainer_Type

Definition at line 170 of file OneDFSISolver.hpp.

◆ feSpace_Type

Definition at line 172 of file OneDFSISolver.hpp.

◆ feSpacePtr_Type

typedef std::shared_ptr< feSpace_Type > feSpacePtr_Type

Definition at line 173 of file OneDFSISolver.hpp.

◆ comm_Type

typedef Epetra_Comm comm_Type

Definition at line 175 of file OneDFSISolver.hpp.

◆ commPtr_Type

typedef std::shared_ptr< comm_Type > commPtr_Type

Definition at line 176 of file OneDFSISolver.hpp.

◆ linearSolver_Type

Definition at line 178 of file OneDFSISolver.hpp.

◆ linearSolverPtr_Type

typedef std::shared_ptr< linearSolver_Type > linearSolverPtr_Type

Definition at line 179 of file OneDFSISolver.hpp.

◆ vector_Type

Definition at line 181 of file OneDFSISolver.hpp.

◆ vectorPtr_Type

typedef std::shared_ptr< vector_Type > vectorPtr_Type

Definition at line 182 of file OneDFSISolver.hpp.

◆ vectorPtrContainer_Type

typedef std::array< vectorPtr_Type, 2 > vectorPtrContainer_Type

Definition at line 183 of file OneDFSISolver.hpp.

◆ matrix_Type

Definition at line 185 of file OneDFSISolver.hpp.

◆ matrixPtr_Type

typedef std::shared_ptr<matrix_Type> matrixPtr_Type

Definition at line 186 of file OneDFSISolver.hpp.

◆ matrixPtrContainer_Type

typedef std::array<matrixPtr_Type, 4 > matrixPtrContainer_Type

Definition at line 187 of file OneDFSISolver.hpp.

◆ solution_Type

typedef std::map< std::string, vectorPtr_Type > solution_Type

Definition at line 189 of file OneDFSISolver.hpp.

◆ solutionPtr_Type

typedef std::shared_ptr< solution_Type > solutionPtr_Type

Definition at line 190 of file OneDFSISolver.hpp.

◆ solutionConstIterator_Type

typedef solution_Type::const_iterator solutionConstIterator_Type

Definition at line 191 of file OneDFSISolver.hpp.

◆ bcLine_Type

Definition at line 193 of file OneDFSISolver.hpp.

◆ bcSide_Type

Definition at line 194 of file OneDFSISolver.hpp.

◆ bcType_Type

Definition at line 195 of file OneDFSISolver.hpp.

Constructor & Destructor Documentation

◆ OneDFSISolver() [1/2]

OneDFSISolver ( )
explicit

Empty Constructor.

Need a call to: setCommunicator(), setProblem(), setFESpace()

Definition at line 61 of file OneDFSISolver.cpp.

+ Here is the caller graph for this function:

◆ ~OneDFSISolver()

virtual ~OneDFSISolver ( )
inlinevirtual

Destructor.

Definition at line 210 of file OneDFSISolver.hpp.

◆ OneDFSISolver() [2/2]

OneDFSISolver ( const OneDFSISolver solver)
explicitprivate

Member Function Documentation

◆ buildConstantMatrices()

void buildConstantMatrices ( )

Build constant matrices (mass and grad)

Definition at line 93 of file OneDFSISolver.cpp.

◆ setupSolution() [1/2]

void setupSolution ( solution_Type solution)
inline

Setup the solution using the default FESpace map.

Parameters
solutionsolution container

Definition at line 225 of file OneDFSISolver.hpp.

◆ setupSolution() [2/2]

void setupSolution ( solution_Type solution,
const MapEpetra map,
const bool &  onlyMainQuantities = false 
)

Setup the solution using user defined FESpace map.

Parameters
solutionsolution container
mapmap for initializing the solution vectors
onlyMainQuantitiesif true setup only $Q$, $P$, and $\displaystyle\frac{A}{A^0}-1$

Definition at line 138 of file OneDFSISolver.cpp.

◆ initialize()

void initialize ( solution_Type solution)

Initialize all the variables of the solution to a reference condition with $Q=0$, $A=A^0$, and $P=P_\mathrm{ext}$.

Parameters
solutionthe solution container

Definition at line 180 of file OneDFSISolver.cpp.

◆ computeW1W2()

void computeW1W2 ( solution_Type solution)

Update the Riemann variables.

Parameters
solutionthe solution container is passed with $A^n$, $Q^n$ and it is updated with $W_1^n$, $W_2^n$

Definition at line 200 of file OneDFSISolver.cpp.

+ Here is the caller graph for this function:

◆ computePressure()

void computePressure ( solution_Type solution,
const Real timeStep 
)

Update the pressure.

This method compute the value of the pressure (elastic and if necessary also viscoelastic) adding it to the solution.

Parameters
solutionthe solution container is passed with $A^n$, $Q^n$, $W_1^n$, $W_2^n$ and is updated with $P^n$
timeSteptime step

Definition at line 210 of file OneDFSISolver.cpp.

+ Here is the caller graph for this function:

◆ computeAreaRatio()

void computeAreaRatio ( solution_Type solution)

Update the ratio between $A$ and $A^0$.

Parameters
solutionthe solution container is passed with $A^n$, is updated with $\displaystyle\frac{A}{A^0}-1$

Definition at line 226 of file OneDFSISolver.cpp.

+ Here is the caller graph for this function:

◆ computeArea()

void computeArea ( solution_Type solution)

Compute A from the area ratio: $\displaystyle\frac{A}{A^0}-1$.

Parameters
solutionthe solution container is passed with $\displaystyle\frac{A}{A^0}-1$ and is updated with $A^n$

Definition at line 235 of file OneDFSISolver.cpp.

◆ updateRHS()

void updateRHS ( const solution_Type solution,
const Real timeStep 
)

Compute the right hand side.

Parameters
solutionthe solution container
timeStepthe time step.

Definition at line 244 of file OneDFSISolver.cpp.

◆ iterate()

void iterate ( OneDFSIBCHandler bcH,
solution_Type solution,
const Real time,
const Real timeStep 
)

Update convective term and BC. Then solve the linearized system.

Parameters
bcHthe BC handler
timethe time
timeStepthe time step

Definition at line 304 of file OneDFSISolver.cpp.

◆ viscoelasticFlowRateCorrection()

OneDFSISolver::vector_Type viscoelasticFlowRateCorrection ( const vector_Type newArea,
const vector_Type newElasticFlowRate,
const vector_Type oldViscoelasticFlowRate,
const Real timeStep,
OneDFSIBCHandler bcHandler,
const bool &  updateSystemMatrix = true 
)

Apply the viscoelastic flow rate correction.

To introduce the viscoelastic component of the wall within the formulation, we use an operator-splitting technique, where the flow rate is splitted into two components such that $Q = \hat{Q} + \tilde{Q}$, where $\hat{Q}$ is the solution of the pure elastic problem and $\tilde{Q}$ is the viscoelastic correction. On each time interval $[t^n, t^{n+1}]$ with $n \ge 0$, firstly we solve the elastic part for $\hat{Q}^{n+1}$, using $Q^n$ to compute the contributions in such equation to the right hand side, and then we correct the flow rate by solving the following equation By using the mass conservation equation, we remove the time dependence from the viscoelastic wall term. The resulting problem is

\[ \displaystyle\frac{1}{A}\displaystyle\frac{\partial \tilde{Q}}{\partial t} - \displaystyle\frac{\partial}{\partial z}\left(\displaystyle\frac{\gamma}{\rho A^{3/2}}\displaystyle\frac{\partial Q}{\partial z}\right) = 0, \]

which is closed by a proper set of homogeneous boundary conditions for $\tilde{Q}$. The corresponding finite element formulation reads: given $(A_h^{n+1},\hat{Q}_h^{n+1})$, find $\tilde{Q}_h^{n+1}$ such that

\[ \left(\displaystyle\frac{\tilde{Q}^{n+1}_h}{A^{n+1}_h},\varphi_h\right) + \Delta t \left(\displaystyle\frac{\gamma}{\rho \left(A^{n+1}_h\right)^{3/2}}\displaystyle\frac{\partial \tilde{Q}^{n+1}_h}{\partial z}, \displaystyle\frac{\partial \varphi_h}{\partial z}\right) = \left(\displaystyle\frac{\tilde{Q}^{n}_h}{A^{n+1}_h},\varphi_h\right) - \Delta t \left(\displaystyle\frac{\gamma}{\rho \left(A^{n+1}_h\right)^{3/2}} \displaystyle\frac{\partial \hat{Q}^{n+1}_h}{\partial z},\displaystyle\frac{\partial \varphi_h}{\partial z}\right)+ \Delta t \left[\displaystyle\frac{\gamma}{\rho \left(A^{n+1}_h\right)^{3/2}}\displaystyle\frac{\partial\hat{Q}^{n+1}_h}{\partial z}\,\varphi_h\right]^L_0. \]

Parameters
areaarea
flowRateflow rate
timeStepthe time step
bcHthe BC handler
updateSystemMatrixflag for the recomputation of the system matrix
Returns
the viscoelastic flow rate correction $\tilde{Q}$

Definition at line 351 of file OneDFSISolver.cpp.

◆ computeCFL()

Real computeCFL ( const solution_Type solution,
const Real timeStep 
) const

CFL computation (correct for constant mesh)

Parameters
solutionthe solution container
timeStepthe time step
Returns
CFL

Definition at line 430 of file OneDFSISolver.cpp.

◆ resetOutput()

void resetOutput ( const solution_Type solution)

Reset the output files.

Parameters
solutionthe solution container

Definition at line 453 of file OneDFSISolver.cpp.

◆ postProcess()

void postProcess ( const solution_Type solution,
const Real time 
)

Save results on output files.

Parameters
solutionsolution container
timesolution time

Definition at line 465 of file OneDFSISolver.cpp.

◆ setProblem()

void setProblem ( const physicsPtr_Type physicsPtr,
const fluxPtr_Type fluxPtr,
const sourcePtr_Type sourcePtr 
)

Set problem classes.

Parameters
physicsPtrpointer to the physics class.
fluxPtrpointer to the flux class.
sourcePtrpointer to the source class.

Definition at line 489 of file OneDFSISolver.cpp.

◆ setCommunicator()

void setCommunicator ( const commPtr_Type commPtr)

Set the communicator.

Parameters
commPtrpointer to the Epetra MPI communicator

Definition at line 499 of file OneDFSISolver.cpp.

◆ setFESpace()

void setFESpace ( const feSpacePtr_Type feSpacePtr)

Set the FEspace.

Parameters
feSpacePtrpointer to the FE space

Definition at line 506 of file OneDFSISolver.cpp.

◆ setLinearSolver()

void setLinearSolver ( const linearSolverPtr_Type linearSolverPtr)

Set the linear solver.

Parameters
linearSolverPtrpointer to the linear solver for the hyperbolic problem

Definition at line 531 of file OneDFSISolver.cpp.

◆ setLinearViscoelasticSolver()

void setLinearViscoelasticSolver ( const linearSolverPtr_Type linearViscoelasticSolverPtr)

Set the viscoelastic linear solver.

Parameters
linearViscoelasticSolverPtrpointer to the linear solver for the viscoelastic problem

Definition at line 537 of file OneDFSISolver.cpp.

◆ physics()

const physicsPtr_Type& physics ( ) const
inline

Get the physics class.

Returns
shared pointer to the physics class.

Definition at line 395 of file OneDFSISolver.hpp.

◆ flux()

const fluxPtr_Type& flux ( ) const
inline

Get the flux class.

Returns
shared pointer to the flux class.

Definition at line 404 of file OneDFSISolver.hpp.

◆ source()

const sourcePtr_Type& source ( ) const
inline

Get the source class.

Returns
shared pointer to the source class.

Definition at line 413 of file OneDFSISolver.hpp.

◆ boundaryDOF()

UInt boundaryDOF ( const bcSide_Type bcSide) const

Return the ID of the boundary node given a side.

Parameters
bcSideSide of the boundary.
Returns
ID of the boundary node.

Definition at line 546 of file OneDFSISolver.cpp.

+ Here is the caller graph for this function:

◆ boundaryValue()

Real boundaryValue ( const solution_Type solution,
const bcType_Type bcType,
const bcSide_Type bcSide 
) const

Return the value of a quantity ( $P$, $A$, $Q$, $W_1$, $W_2$) on a specified boundary.

Given a bcType and a bcSide it return the value of the quantity.

Parameters
bcTypeType of the asked boundary value.
bcSideSide of the boundary.
Returns
value of the quantity on the specified side.

Definition at line 571 of file OneDFSISolver.cpp.

◆ boundaryEigenValuesEigenVectors()

void boundaryEigenValuesEigenVectors ( const bcSide_Type bcSide,
const solution_Type solution,
container2D_Type eigenvalues,
container2D_Type leftEigenvector1,
container2D_Type leftEigenvector2 
)

Return the value of the eigenvalues and eigenvectors on a specified boundary.

Parameters
bcSideSide of the boundary.
solutionsolution container.
eigenvaluesoutput eigenvalues.
leftEigenvector1output left eigenvector associated to the first eigenvalue.
leftEigenvector1output left eigenvector associated to the second eigenvalue.

Definition at line 623 of file OneDFSISolver.cpp.

◆ residual()

const vectorPtrContainer_Type& residual ( ) const
inline

Get the residual container.

Returns
System residual container

Definition at line 451 of file OneDFSISolver.hpp.

◆ massMatrix()

const matrixPtr_Type& massMatrix ( ) const
inline

Get the system matrix without BC.

Returns
shared pointer to the system matrix without BC

Definition at line 460 of file OneDFSISolver.hpp.

◆ updateFlux()

void updateFlux ( const solution_Type solution)
private

Update the P1 flux vector from U: M_fluxi = F_h(Un) i=1,2 (works only for P1Seg elements)

Definition at line 641 of file OneDFSISolver.cpp.

+ Here is the caller graph for this function:

◆ updatedFdU()

void updatedFdU ( const solution_Type solution)
private

Call _updateFlux and update the P0 derivative of flux vector from U:

M_diffFluxij = dF_h/dU(Un) i,j=1,2 M_diffFluxij(elem) = 1/2 [ dF/dU(U(node1(elem))) + dF/dU(U(node2(elem))) ]

(mean value of the two extremal values of dF/dU) BEWARE: works only for P1Seg elements

Definition at line 656 of file OneDFSISolver.cpp.

+ Here is the caller graph for this function:

◆ updateSource()

void updateSource ( const solution_Type solution)
private

Update the P1 source vector from U: M_sourcei = S_h(Un) i=1,2 (works only for P1Seg elements)

Definition at line 689 of file OneDFSISolver.cpp.

+ Here is the caller graph for this function:

◆ updatedSdU()

void updatedSdU ( const solution_Type solution)
private

Call _updateSource and update the P0 derivative of source vector from U:

M_diffSrcij = dS_h/dU(Un) i,j=1,2 M_diffSrcij(elem) = 1/2 [ dS/dU(U(node1(elem))) + dS/dU(U(node2(elem))) ]

(mean value of the two extremal values of dS/dU) BEWARE: works only for P1Seg elements

Definition at line 704 of file OneDFSISolver.cpp.

+ Here is the caller graph for this function:

◆ updateMatrices()

void updateMatrices ( )
private

Update the matrices.

M_massMatrixDiffSrcij, M_stiffMatrixDiffFluxij M_gradMatrixDiffFluxij, and M_divMatrixDiffSrcij (i,j=1,2)

from the values of diffFlux(Un) and diffSrc(Un) that are computed with _updateMatrixCoefficients.

call of _updateMatrixCoefficients, _updateMatrixElementalrices and _assemble_matrices.

Definition at line 736 of file OneDFSISolver.cpp.

+ Here is the caller graph for this function:

◆ updateElementalMatrices()

void updateElementalMatrices ( const Real dFdU,
const Real dSdU 
)
private

Update the element matrices with the current element.

Update the gradient matrix gradient operator: grad_{ij} = {fe} coeff {d }{d x}

BEWARE :

Parameters
0the first argument "0" corresponds to the first and only coordinate (1D!), and HERE it starts from 0... (Damm'!)
-M_coeffGrad: the sign "-" in the second argument is added to correspond to the described operator. (There is a minus in the elemOper implementation).

update the divergence matrix divergence operator: (transpose of the gradient) div_{ij} = {fe} coeff {d }{d x}

Note
formally this M_elementalDivergenceMatrixPtr is not necessary as it is the transpose of the M_elementalGradientMatrixPtr. But for the sake of clarity, I prefer to keep it. (low cost!)

BEWARE : same remarks as grad (see above).

Definition at line 768 of file OneDFSISolver.cpp.

◆ matrixAssemble()

void matrixAssemble ( const UInt ii,
const UInt jj 
)
private

Assemble the matrices.

Definition at line 819 of file OneDFSISolver.cpp.

◆ applyDirichletBCToMatrix()

void applyDirichletBCToMatrix ( matrix_Type matrix)
private

Update the matrices to take into account Dirichlet BC.

Modify the matrix to take into account the Dirichlet boundary conditions (works for P1Seg and canonic numbering!)

Definition at line 835 of file OneDFSISolver.cpp.

◆ inertialFlowRateCorrection()

OneDFSISolver::vector_Type inertialFlowRateCorrection ( const vector_Type flux)
private

Apply the inertial Flux correction:

We use a finite element scheme for the correction term: given the solution of Taylor-Galerkin scheme, solve ( 1/Ah(n+1) Qtildeh(n+1), phi) + //! 1/A * massFactor^{-1} * Un+1 ( m / rho ) * ( dQtildeh(n+1)/dz, dphi/dz ) //! stiff * Qtilde(U) = ( m / rho ) * ( dQhath(n+1)/dz, dphi/dz ) //! stiff * Qhat(U)

m = rho_w h0 / ( 2 sqrt(pi) sqrt(A0) )

Definition at line 846 of file OneDFSISolver.cpp.

◆ longitudinalFlowRateCorrection()

OneDFSISolver::vector_Type longitudinalFlowRateCorrection ( )
private

Apply the longitudinal Flux correction:

We use a finite element scheme for the correction term: given the solution of Taylor-Galerkin scheme, solve ( 1/Ah(n+1) Qtildeh(n+1), phi) + //! 1/A * massFactor^{-1} * Un+1 = ( 1/Ah(n+1) Qtildeh(n), phi) + //! 1/A * massFactor^{-1} * Un+1

  • ( a / rho ) * ( d3Ahath(n+1)/dz3, phi ) //! mass * d3Ahat(U)/dz

Definition at line 940 of file OneDFSISolver.cpp.

◆ operator=()

OneDFSISolver& operator= ( const OneDFSISolver solver)
private

Field Documentation

◆ M_physicsPtr

physicsPtr_Type M_physicsPtr
private

L2 Projection of the second derivative of Q over P1 space.

Definition at line 572 of file OneDFSISolver.hpp.

◆ M_fluxPtr

fluxPtr_Type M_fluxPtr
private

Definition at line 573 of file OneDFSISolver.hpp.

◆ M_sourcePtr

sourcePtr_Type M_sourcePtr
private

Definition at line 574 of file OneDFSISolver.hpp.

◆ M_feSpacePtr

feSpacePtr_Type M_feSpacePtr
private

Definition at line 575 of file OneDFSISolver.hpp.

◆ M_commPtr

commPtr_Type M_commPtr
private

Definition at line 576 of file OneDFSISolver.hpp.

◆ M_displayer

Displayer M_displayer
private

Definition at line 577 of file OneDFSISolver.hpp.

◆ M_elementalMassMatrixPtr

std::shared_ptr< MatrixElemental > M_elementalMassMatrixPtr
private

element mass matrix

Definition at line 579 of file OneDFSISolver.hpp.

◆ M_elementalStiffnessMatrixPtr

std::shared_ptr< MatrixElemental > M_elementalStiffnessMatrixPtr
private

element stiffness matrix

Definition at line 580 of file OneDFSISolver.hpp.

◆ M_elementalGradientMatrixPtr

std::shared_ptr< MatrixElemental > M_elementalGradientMatrixPtr
private

element gradient matrix

Definition at line 581 of file OneDFSISolver.hpp.

◆ M_elementalDivergenceMatrixPtr

std::shared_ptr< MatrixElemental > M_elementalDivergenceMatrixPtr
private

element divergence matrix

Definition at line 582 of file OneDFSISolver.hpp.

◆ M_rhs

vectorPtrContainer_Type M_rhs
private

Right hand sides of the linear system i: "mass * M_Ui = M_rhsi".

Definition at line 585 of file OneDFSISolver.hpp.

◆ M_residual

vectorPtrContainer_Type M_residual
private

Residual of the linear system.

Definition at line 588 of file OneDFSISolver.hpp.

◆ M_fluxVector

vectorPtrContainer_Type M_fluxVector
private

Flux F(U) (in P1)

Definition at line 591 of file OneDFSISolver.hpp.

◆ M_sourceVector

vectorPtrContainer_Type M_sourceVector
private

Source term S (in P1)

Definition at line 594 of file OneDFSISolver.hpp.

◆ M_dFdUVector

scalarVectorContainer_Type M_dFdUVector
private

diffFlux = dF(U)/dU (in P0)

Definition at line 597 of file OneDFSISolver.hpp.

◆ M_dSdUVector

scalarVectorContainer_Type M_dSdUVector
private

diffSrc = dSource(U)/dU (in P0)

Definition at line 600 of file OneDFSISolver.hpp.

◆ M_homogeneousMassMatrixPtr

matrixPtr_Type M_homogeneousMassMatrixPtr
private

tridiagonal mass matrix

Definition at line 603 of file OneDFSISolver.hpp.

◆ M_homogeneousGradientMatrixPtr

matrixPtr_Type M_homogeneousGradientMatrixPtr
private

tridiagonal gradient matrix

Definition at line 606 of file OneDFSISolver.hpp.

◆ M_dSdUMassMatrixPtr

matrixPtrContainer_Type M_dSdUMassMatrixPtr
private

tridiagonal mass matrices multiplied by diffSrcij

Definition at line 609 of file OneDFSISolver.hpp.

◆ M_dFdUStiffnessMatrixPtr

matrixPtrContainer_Type M_dFdUStiffnessMatrixPtr
private

tridiagonal stiffness matrices multiplied by diffFluxij

Definition at line 612 of file OneDFSISolver.hpp.

◆ M_dFdUGradientMatrixPtr

matrixPtrContainer_Type M_dFdUGradientMatrixPtr
private

tridiagonal gradient matrices multiplied by diffFluxij

Definition at line 615 of file OneDFSISolver.hpp.

◆ M_dSdUDivergenceMatrixPtr

matrixPtrContainer_Type M_dSdUDivergenceMatrixPtr
private

tridiagonal divergence matrices multiplied by diffSrcij

Definition at line 618 of file OneDFSISolver.hpp.

◆ M_linearSolverPtr

linearSolverPtr_Type M_linearSolverPtr
private

The linear solver.

Definition at line 621 of file OneDFSISolver.hpp.

◆ M_linearViscoelasticSolverPtr

linearSolverPtr_Type M_linearViscoelasticSolverPtr
private

Definition at line 622 of file OneDFSISolver.hpp.


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