LifeV
OneDFSIData Class Reference

OneDFSIData - Class which read and holds all the data for the One Dimensional Model Solver. More...

#include <OneDFSIData.hpp>

+ Collaboration diagram for OneDFSIData:

Private Attributes

OneDFSI::physicsType_Type M_physicsType
 Model. More...
 
OneDFSI::fluxTerm_Type M_fluxType
 
OneDFSI::sourceTerm_Type M_sourceType
 
timePtr_Type M_timeDataPtr
 Data containers for time and mesh. More...
 
meshPtr_Type M_meshPtr
 
bool M_viscoelasticWall
 Physical Wall Model. More...
 
Real M_viscoelasticAngle
 
Real M_viscoelasticPeriod
 
scalarVector_Type M_viscoelasticCoefficient
 
bool M_inertialWall
 
Real M_densityWall
 
Real M_inertialModulus
 
bool M_longitudinalWall
 
std::string M_postprocessingDirectory
 Miscellaneous. More...
 
std::string M_postprocessingFile
 full directory name (including path) More...
 
Real M_CFLmax
 output file name More...
 
Real M_jacobianPerturbationArea
 Jacobian perturbation. More...
 
Real M_jacobianPerturbationFlowRate
 
Real M_jacobianPerturbationStress
 
bool M_computeCoefficients
 Physical Parameters. More...
 
Int M_powerLawCoefficient
 
Real M_density
 
Real M_viscosity
 
bool M_thickVessel
 
Real M_young
 
Real M_poisson
 
Real M_externalPressure
 
Real M_venousPressure
 
Real M_robertsonCorrection
 
scalarVector_Type M_thickness
 
Real M_friction
 
scalarVector_Type M_area0
 
scalarVector_Type M_alpha
 
scalarVector_Type M_beta0
 
scalarVector_Type M_beta1
 
scalarVector_Type M_dArea0dz
 Derivatives of main coefficients. More...
 
scalarVector_Type M_dAlphadz
 
scalarVector_Type M_dBeta0dz
 
scalarVector_Type M_dBeta1dz
 
scalarVector_Type M_flux11
 Flux matrix. More...
 
scalarVector_Type M_flux12
 
scalarVector_Type M_flux21
 
scalarVector_Type M_flux22
 
scalarVector_Type M_celerity1
 Celerities of the linear problem (eigenvalues of the flux matrix) More...
 
scalarVector_Type M_celerity2
 
scalarVector_Type M_celerity1LeftEigenvector1
 Eigenvector for first and second eigenvalue. More...
 
scalarVector_Type M_celerity1LeftEigenvector2
 
scalarVector_Type M_celerity2LeftEigenvector1
 
scalarVector_Type M_celerity2LeftEigenvector2
 
scalarVector_Type M_source10
 Source matrix. More...
 
scalarVector_Type M_source20
 
scalarVector_Type M_source11
 
scalarVector_Type M_source12
 
scalarVector_Type M_source21
 
scalarVector_Type M_source22
 

Type definitions

enum  OneD_distributionLaw { uniform, linear, pointwise }
 
typedef TimeData time_Type
 
typedef std::shared_ptr< time_TypetimePtr_Type
 
typedef RegionMesh< LinearLinemesh_Type
 
typedef std::shared_ptr< mesh_TypemeshPtr_Type
 
typedef ublas::vector< RealscalarVector_Type
 
typedef std::array< Real, 2 > container2D_Type
 

Constructors & Destructor

 OneDFSIData ()
 Empty constructor. More...
 
virtual ~OneDFSIData ()
 Destructor. More...
 

Methods

void setup (const GetPot &dataFile, const std::string &section="1D_Model")
 Setup method. More...
 
void LIFEV_DEPRECATED (oldStyleSetup(const GetPot &dataFile, const std::string &section="1dnetwork"))
 Deprecated setup method. (DO NOT USE THIS) More...
 
void updateCoefficients ()
 Update all the physical coefficients. More...
 
template<typename VectorType >
Real computeSpatialDerivativeAtNode (const VectorType &vector, const UInt &iNode, const UInt &bcFiniteDifferenceOrder=2)
 Compute the spatial derivative of a quantity at a node. More...
 
void showMe (std::ostream &output=std::cout) const
 Initialize linear parameters (NOT WORKING) More...
 

Set Methods

void setPostprocessingDirectory (const std::string &directory)
 Set the post-processing directory. More...
 
void setPostprocessingFile (const std::string &file)
 Set the post-processing file name. More...
 
void setJacobianPerturbationArea (const Real &jacobianPerturbationArea)
 Set the area perturbation parameter to compute the Jacobian matrix (in the Multiscale framework) More...
 
void setJacobianPerturbationFlowRate (const Real &jacobianPerturbationFlowRate)
 Set the flow rate perturbation parameter to compute the Jacobian matrix (in the Multiscale framework) More...
 
void setJacobianPerturbationStress (const Real &jacobianPerturbationStress)
 Set the stress perturbation parameter to compute the Jacobian matrix (in the Multiscale framework) More...
 
void setTimeData (const timePtr_Type timeDataPtr)
 Set data time container. More...
 
void setDensity (const Real &density)
 Set the fluid density. More...
 
void setViscosity (const Real &viscosity)
 Set the fluid viscosity. More...
 
void setDensityWall (const Real &densityWall)
 Set the wall density. More...
 
void setThickness (const Real &thickness, const UInt &i)
 Set the wall thickness. More...
 
void setYoung (const Real &young)
 Set the wall Young modulus. More...
 
void setPoisson (const Real &poisson)
 Set the wall Poisson number. More...
 
void setExternalPressure (const Real &externalPressure)
 Set the wall external pressure. More...
 
void setVenousPressure (const Real &venousPressure)
 Set the venous pressure at the terminals. More...
 
void setArea0 (const Real &area0, const UInt &i)
 Set the reference area. More...
 
void setBeta0 (const Real &beta0, const UInt &i)
 Set the wall $\beta_0$. More...
 
void setdBeta0dz (const Real &dBeta0dz, const UInt &i)
 Set the wall $\frac{d\beta_0}{dz}$. More...
 

Get Methods

const OneDFSI::physicsType_TypephysicsType () const
 Get the physics type. More...
 
const OneDFSI::fluxTerm_TypefluxType () const
 Get the flux type. More...
 
const OneDFSI::sourceTerm_TypesourceType () const
 Get the source type. More...
 
timePtr_Type dataTime () const
 Get data time container. More...
 
meshPtr_Type mesh () const
 Get the mesh container. More...
 
Real length () const
 Get the length of the 1D segment. More...
 
UInt numberOfElements () const
 Get the number of elements in the 1D segment. More...
 
UInt numberOfNodes () const
 Get the number of nodes in the 1D segment. More...
 
const bool & viscoelasticWall () const
 Get the flag identifying if the wall is viscoelastic. More...
 
const RealviscoelasticCoefficient (const UInt &i) const
 Get the viscoelastic coefficient $\gamma$. More...
 
const bool & inertialWall () const
 Get the flag identifying if the wall has inertia. More...
 
const RealdensityWall () const
 Get the density of the wall. More...
 
const RealinertialModulus () const
 Get the inertial coefficient (to be defined) More...
 
const bool & longitudinalWall () const
 Get the flag identifying if the wall has a longitudinal pre-stress. More...
 
const std::string & postprocessingDirectory () const
 Get the post-processing directory. More...
 
const std::string & postprocessingFile () const
 Get the post-processing file. More...
 
const RealCFLmax () const
 Get the imposed CFL condition. More...
 
const RealjacobianPerturbationArea () const
 Get the area perturbation parameter to compute the Jacobian matrix (in the Multiscale framework) More...
 
const RealjacobianPerturbationFlowRate () const
 Get the flow rate perturbation parameter to compute the Jacobian matrix (in the Multiscale framework) More...
 
const RealjacobianPerturbationStress () const
 Get the stress perturbation parameter to compute the Jacobian matrix (in the Multiscale framework) More...
 
const RealdensityRho () const
 Get the fluid density. More...
 
const Realviscosity () const
 Get the fluid viscosity. More...
 
const Realyoung () const
 Get the wall Young modulus. More...
 
const Realpoisson () const
 Get the wall Poisson number. More...
 
const RealexternalPressure () const
 Get the wall external pressure. More...
 
const RealvenousPressure () const
 Get the venous pressure. More...
 
const RealrobertsonCorrection () const
 Get the Robertson correction coefficient (Not tested: maybe wrong in the code) More...
 
const Realthickness (const UInt &i) const
 Get the wall thickness. More...
 
const Realfriction () const
 Get the wall friction coefficient $K_r$. More...
 
const Realarea0 (const UInt &i) const
 Get the reference area $A^0$. More...
 
const Realalpha (const UInt &i) const
 Get the Coriolis coefficient $\alpha$. More...
 
const Realbeta0 (const UInt &i) const
 Get the $\beta_0$ coefficient. More...
 
const Realbeta1 (const UInt &i) const
 Get the $\beta_1$ coefficient. More...
 
const RealdArea0dz (const UInt &i) const
 Get $\frac{dA^0}{dz}$. More...
 
const RealdAlphadz (const UInt &i) const
 Get $\frac{d\alpha}{dz}$. More...
 
const RealdBeta0dz (const UInt &i) const
 Get $\frac{d\beta_0}{dz}$. More...
 
const RealdBeta1dz (const UInt &i) const
 Get $\frac{d\beta_1}{dz}$. More...
 
const Realflux11 (const UInt &i) const
 Get the flux coefficient $F_{11}$ (used only in the linear problem) More...
 
const Realflux12 (const UInt &i) const
 Get the flux coefficient $F_{12}$ (used only in the linear problem) More...
 
const Realflux21 (const UInt &i) const
 Get the flux coefficient $F_{21}$ (used only in the linear problem) More...
 
const Realflux22 (const UInt &i) const
 Get the flux coefficient $F_{22}$ (used only in the linear problem) More...
 
const Realcelerity1 (const UInt &i) const
 Get the first eigenvector $\lambda_{1}$ (used only in the linear problem) More...
 
const Realcelerity2 (const UInt &i) const
 Get the second eigenvector $\lambda_{2}$ (used only in the linear problem) More...
 
const RealleftEigenVector11 (const UInt &i) const
 Get the left eigenvector coefficient $L_{11}$ (used only in the linear problem) More...
 
const RealleftEigenVector12 (const UInt &i) const
 Get the left eigenvector coefficient $L_{12}$ (used only in the linear problem) More...
 
const RealleftEigenVector21 (const UInt &i) const
 Get the left eigenvector coefficient $L_{21}$ (used only in the linear problem) More...
 
const RealleftEigenVector22 (const UInt &i) const
 Get the left eigenvector coefficient $L_{22}$ (used only in the linear problem) More...
 
const Realsource10 (const UInt &i) const
 Get the source coefficient $S_{10}$ (used only in the linear problem) More...
 
const Realsource20 (const UInt &i) const
 Get the source coefficient $S_{20}$ (used only in the linear problem) More...
 
const Realsource11 (const UInt &i) const
 Get the source coefficient $S_{11}$ (used only in the linear problem) More...
 
const Realsource12 (const UInt &i) const
 Get the source coefficient $S_{12}$ (used only in the linear problem) More...
 
const Realsource21 (const UInt &i) const
 Get the source coefficient $S_{21}$ (used only in the linear problem) More...
 
const Realsource22 (const UInt &i) const
 Get the source coefficient $S_{22}$ (used only in the linear problem) More...
 

Private Methods

void linearInterpolation (scalarVector_Type &vector, const GetPot &dataFile, const std::string &quantity, const Real &defaultValue, const bool &isArea=false)
 Compute the linear interpolation of a general quantity. More...
 
void computeSpatialDerivatives ()
 Compute the derivatives of alpha, area0, beta0, and beta1 using centered differences. More...
 
void resetContainers ()
 Reset all the containers. More...
 

Detailed Description

OneDFSIData - Class which read and holds all the data for the One Dimensional Model Solver.

Authors
Vincent Martin, Cristiano Malossi
See also
Equations and networks of 1-D models [7]
Geometrical multiscale coupling of 1-D models [12] [13] [3]

Physical Parameters

Main parameters: $A^0, \alpha, \beta_0, \beta_1, \gamma, K_r, \rho$.

Euler equations:

\[ \left\{\begin{array}{l} \displaystyle\frac{\partial A}{\partial t} + \frac{\partial Q}{\partial z} = 0, \\[2ex] \displaystyle\frac{\partial Q}{\partial t} + \alpha \frac{\partial}{\partial z}\left(\frac{Q^2}{A}\right) + \frac{A}{\rho}\frac{\partial P}{\partial z} + K_r \frac{Q}{A} = 0 \end{array}\right. \]

with

\[ P-P_\mathrm{ext} = \psi(A,A^0,\beta_0, \beta_1, \gamma) = \underbrace{\sqrt{\frac{\pi}{A^0}}\frac{h E}{1-\nu^2}}_{\beta_0} \left(\left(\frac{A}{A^0}\right)^{\beta_1}-1\right) + \underbrace{\frac{T \tan\phi}{4 \sqrt{\pi}}\frac{h E}{1-\nu^2}}_{\displaystyle\gamma} \frac{1}{A\sqrt{A}} \frac{\partial A}{\partial t}, \]

Linear Parameters

Parameters: $F_{11}, F_{12}, F_{21}, F_{22}, \lambda_1, \lambda_2$

Equations:

\[ \left\{\begin{array}{l} \displaystyle\frac{dU_1}{dt} + F_{11} \frac{dU_1}{dz} + F_{12} \frac{dU_2}{dz} = 0\\[2ex] \displaystyle\frac{dU_2}{dt} + F_{21} \frac{dU_1}{dz} + F_{22} \frac{dU_2}{dz} = 0 \end{array}\right. \]

The flux matrix $\mathbf F = [F_{11}, F_{12}; F_{21}, F_{22}]$ has the eigenvalues $\lambda_1, \lambda_2$.

Definition at line 106 of file OneDFSIData.hpp.

Member Typedef Documentation

◆ time_Type

Definition at line 113 of file OneDFSIData.hpp.

◆ timePtr_Type

typedef std::shared_ptr< time_Type > timePtr_Type

Definition at line 114 of file OneDFSIData.hpp.

◆ mesh_Type

Definition at line 116 of file OneDFSIData.hpp.

◆ meshPtr_Type

typedef std::shared_ptr< mesh_Type > meshPtr_Type

Definition at line 117 of file OneDFSIData.hpp.

◆ scalarVector_Type

typedef ublas::vector< Real > scalarVector_Type

Definition at line 120 of file OneDFSIData.hpp.

◆ container2D_Type

typedef std::array< Real, 2 > container2D_Type

Definition at line 121 of file OneDFSIData.hpp.

Member Enumeration Documentation

◆ OneD_distributionLaw

Enumerator
uniform 
linear 
pointwise 

Definition at line 123 of file OneDFSIData.hpp.

Constructor & Destructor Documentation

◆ OneDFSIData()

OneDFSIData ( )
explicit

Empty constructor.

Definition at line 51 of file OneDFSIData.cpp.

+ Here is the caller graph for this function:

◆ ~OneDFSIData()

virtual ~OneDFSIData ( )
inlinevirtual

Destructor.

Definition at line 140 of file OneDFSIData.hpp.

Member Function Documentation

◆ setup()

void setup ( const GetPot dataFile,
const std::string &  section = "1D_Model" 
)

Setup method.

Parameters
dataFileGetPot dataFile
sectionsection in the dataFile

Definition at line 114 of file OneDFSIData.cpp.

◆ LIFEV_DEPRECATED()

void LIFEV_DEPRECATED ( oldStyleSetup(const GetPot &dataFile, const std::string &section="1dnetwork")  )

Deprecated setup method. (DO NOT USE THIS)

This method has been implemented for compatibility reason with some old applications. Please, for any new application use the setup method in place of this one. Future compatibility for this method is not guaranteed.

Parameters
dataFileGetPot dataFile
Returns
section section in the dataFile

◆ updateCoefficients()

void updateCoefficients ( )

Update all the physical coefficients.

This method can be called after any update in the physical coefficient. It recomputes the main coefficients $\alpha, \beta_0, \beta_1, \gamma, K_r$

Definition at line 416 of file OneDFSIData.cpp.

+ Here is the caller graph for this function:

◆ computeSpatialDerivativeAtNode()

Real computeSpatialDerivativeAtNode ( const VectorType &  vector,
const UInt iNode,
const UInt bcFiniteDifferenceOrder = 2 
)
inline

Compute the spatial derivative of a quantity at a node.

Note: works only for homogeneous discretizations.

Parameters
vectorthe quantity vector
iNodenode
Returns
spatial derivative

Definition at line 1022 of file OneDFSIData.hpp.

◆ showMe()

void showMe ( std::ostream &  output = std::cout) const

Initialize linear parameters (NOT WORKING)

The linearization of the Euler model yields

\[ \begin{array}{l} \displaystyle F = [ Q; A c^2]\\[2ex] \displaystyle B = [ 0; k_R / A^0]\\[2ex] \displaystyle c = \sqrt{\frac{\beta_0\beta_1}{\rho} } \end{array} \]

the vessel stiffer on the left side of interval [xl, xr] These routines change the elastic modulus along the vessel

When x < alpha - delta/2, the Young modulus is E * factor When x > alpha + delta/2, the Young modulus is E When alpha - delta/2 < x < alpha + delta/2, the Young modulus changes smoothly from the larger to the smaller value, according to a polynomial law of order n.

The grid size can be adapted (yesadaptive=1) in the nieghborhood of alpha, where the spatial derivative of the parameter will be maximum. However, the grid size is not allowed to be smaller than min_deltaxMake the vessel stiffer on the right side of interval [xl, xr]

See also
stiffenVesselLeftDisplay some information about the model.
Parameters
outputStream where the informations must be printed

Definition at line 691 of file OneDFSIData.cpp.

◆ setPostprocessingDirectory()

void setPostprocessingDirectory ( const std::string &  directory)
inline

Set the post-processing directory.

Parameters
directorypost-processing directory

Definition at line 245 of file OneDFSIData.hpp.

◆ setPostprocessingFile()

void setPostprocessingFile ( const std::string &  file)
inline

Set the post-processing file name.

Parameters
filepost-processing file name

Definition at line 254 of file OneDFSIData.hpp.

◆ setJacobianPerturbationArea()

void setJacobianPerturbationArea ( const Real jacobianPerturbationArea)
inline

Set the area perturbation parameter to compute the Jacobian matrix (in the Multiscale framework)

Parameters
jacobianPerturbationAreaarea perturbation parameter

Definition at line 263 of file OneDFSIData.hpp.

◆ setJacobianPerturbationFlowRate()

void setJacobianPerturbationFlowRate ( const Real jacobianPerturbationFlowRate)
inline

Set the flow rate perturbation parameter to compute the Jacobian matrix (in the Multiscale framework)

Parameters
jacobianPerturbationFlowRateflow rate perturbation parameter

Definition at line 272 of file OneDFSIData.hpp.

◆ setJacobianPerturbationStress()

void setJacobianPerturbationStress ( const Real jacobianPerturbationStress)
inline

Set the stress perturbation parameter to compute the Jacobian matrix (in the Multiscale framework)

Parameters
jacobianPerturbationStressstress perturbation parameter

Definition at line 281 of file OneDFSIData.hpp.

◆ setTimeData()

void setTimeData ( const timePtr_Type  timeDataPtr)
inline

Set data time container.

Parameters
timeDataPtrshared_ptr to TimeData container

Definition at line 290 of file OneDFSIData.hpp.

◆ setDensity()

void setDensity ( const Real density)
inline

Set the fluid density.

Parameters
densityfluid density

Definition at line 299 of file OneDFSIData.hpp.

◆ setViscosity()

void setViscosity ( const Real viscosity)
inline

Set the fluid viscosity.

Parameters
viscosityfluid viscosity

Definition at line 308 of file OneDFSIData.hpp.

◆ setDensityWall()

void setDensityWall ( const Real densityWall)
inline

Set the wall density.

Parameters
densityWallwall density

Definition at line 317 of file OneDFSIData.hpp.

◆ setThickness()

void setThickness ( const Real thickness,
const UInt i 
)
inline

Set the wall thickness.

Parameters
inode id
thicknesswall thickness

Definition at line 327 of file OneDFSIData.hpp.

◆ setYoung()

void setYoung ( const Real young)
inline

Set the wall Young modulus.

Parameters
youngwall Young modulus

Definition at line 336 of file OneDFSIData.hpp.

◆ setPoisson()

void setPoisson ( const Real poisson)
inline

Set the wall Poisson number.

Parameters
poissonwall Poisson number

Definition at line 345 of file OneDFSIData.hpp.

◆ setExternalPressure()

void setExternalPressure ( const Real externalPressure)
inline

Set the wall external pressure.

Parameters
externalPressurewall external pressure

Definition at line 354 of file OneDFSIData.hpp.

◆ setVenousPressure()

void setVenousPressure ( const Real venousPressure)
inline

Set the venous pressure at the terminals.

Parameters
venousPressurevenous pressure at the terminals

Definition at line 363 of file OneDFSIData.hpp.

◆ setArea0()

void setArea0 ( const Real area0,
const UInt i 
)
inline

Set the reference area.

Parameters
inode id
area0reference area

Definition at line 373 of file OneDFSIData.hpp.

◆ setBeta0()

void setBeta0 ( const Real beta0,
const UInt i 
)
inline

Set the wall $\beta_0$.

Parameters
inode id
beta0wall $\beta_0$

Definition at line 383 of file OneDFSIData.hpp.

◆ setdBeta0dz()

void setdBeta0dz ( const Real dBeta0dz,
const UInt i 
)
inline

Set the wall $\frac{d\beta_0}{dz}$.

Parameters
inode id
thicknesswall $\frac{d\beta_0}{dz}$

Definition at line 393 of file OneDFSIData.hpp.

◆ physicsType()

const OneDFSI::physicsType_Type& physicsType ( ) const
inline

Get the physics type.

Returns
Physics type

Definition at line 408 of file OneDFSIData.hpp.

◆ fluxType()

const OneDFSI::fluxTerm_Type& fluxType ( ) const
inline

Get the flux type.

Returns
Flux type

Definition at line 417 of file OneDFSIData.hpp.

◆ sourceType()

const OneDFSI::sourceTerm_Type& sourceType ( ) const
inline

Get the source type.

Returns
Source type

Definition at line 426 of file OneDFSIData.hpp.

◆ dataTime()

timePtr_Type dataTime ( ) const
inline

Get data time container.

Returns
shared_ptr to TimeData container

Definition at line 435 of file OneDFSIData.hpp.

◆ mesh()

meshPtr_Type mesh ( ) const
inline

Get the mesh container.

Returns
shared_ptr to the mesh

Definition at line 444 of file OneDFSIData.hpp.

◆ length()

Real length ( ) const
inline

Get the length of the 1D segment.

Returns
length of the 1D segment

Definition at line 453 of file OneDFSIData.hpp.

◆ numberOfElements()

UInt numberOfElements ( ) const
inline

Get the number of elements in the 1D segment.

Returns
number of elements in the 1D segment

Definition at line 462 of file OneDFSIData.hpp.

◆ numberOfNodes()

UInt numberOfNodes ( ) const
inline

Get the number of nodes in the 1D segment.

Returns
number of nodes in the 1D segment

Definition at line 471 of file OneDFSIData.hpp.

◆ viscoelasticWall()

const bool& viscoelasticWall ( ) const
inline

Get the flag identifying if the wall is viscoelastic.

Returns
true if the wall is viscoelastic, false otherwise

Definition at line 480 of file OneDFSIData.hpp.

◆ viscoelasticCoefficient()

const Real& viscoelasticCoefficient ( const UInt i) const
inline

Get the viscoelastic coefficient $\gamma$.

Returns
viscoelastic coefficient $\gamma$

Definition at line 489 of file OneDFSIData.hpp.

◆ inertialWall()

const bool& inertialWall ( ) const
inline

Get the flag identifying if the wall has inertia.

Returns
true if the wall has intertia, false otherwise

Definition at line 498 of file OneDFSIData.hpp.

◆ densityWall()

const Real& densityWall ( ) const
inline

Get the density of the wall.

Returns
density of the wall

Definition at line 507 of file OneDFSIData.hpp.

◆ inertialModulus()

const Real& inertialModulus ( ) const
inline

Get the inertial coefficient (to be defined)

Returns
inertial coefficient (to be defined)

Definition at line 516 of file OneDFSIData.hpp.

◆ longitudinalWall()

const bool& longitudinalWall ( ) const
inline

Get the flag identifying if the wall has a longitudinal pre-stress.

Returns
true if the wall has a longitudinal pre-stress, false otherwise

Definition at line 525 of file OneDFSIData.hpp.

◆ postprocessingDirectory()

const std::string& postprocessingDirectory ( ) const
inline

Get the post-processing directory.

Returns
post-processing directory

Definition at line 534 of file OneDFSIData.hpp.

◆ postprocessingFile()

const std::string& postprocessingFile ( ) const
inline

Get the post-processing file.

Returns
post-processing file

Definition at line 543 of file OneDFSIData.hpp.

◆ CFLmax()

const Real& CFLmax ( ) const
inline

Get the imposed CFL condition.

Returns
imposed CFL condition

Definition at line 552 of file OneDFSIData.hpp.

◆ jacobianPerturbationArea()

const Real& jacobianPerturbationArea ( ) const
inline

Get the area perturbation parameter to compute the Jacobian matrix (in the Multiscale framework)

Returns
area perturbation parameter

Definition at line 561 of file OneDFSIData.hpp.

◆ jacobianPerturbationFlowRate()

const Real& jacobianPerturbationFlowRate ( ) const
inline

Get the flow rate perturbation parameter to compute the Jacobian matrix (in the Multiscale framework)

Returns
flow rate perturbation parameter

Definition at line 570 of file OneDFSIData.hpp.

◆ jacobianPerturbationStress()

const Real& jacobianPerturbationStress ( ) const
inline

Get the stress perturbation parameter to compute the Jacobian matrix (in the Multiscale framework)

Returns
stress perturbation parameter

Definition at line 579 of file OneDFSIData.hpp.

◆ densityRho()

const Real& densityRho ( ) const
inline

Get the fluid density.

Returns
fluid density

Definition at line 588 of file OneDFSIData.hpp.

◆ viscosity()

const Real& viscosity ( ) const
inline

Get the fluid viscosity.

Returns
fluid viscosity

Definition at line 597 of file OneDFSIData.hpp.

◆ young()

const Real& young ( ) const
inline

Get the wall Young modulus.

Returns
wall Young modulus

Definition at line 606 of file OneDFSIData.hpp.

◆ poisson()

const Real& poisson ( ) const
inline

Get the wall Poisson number.

Returns
wall Poisson number

Definition at line 615 of file OneDFSIData.hpp.

◆ externalPressure()

const Real& externalPressure ( ) const
inline

Get the wall external pressure.

Returns
wall external pressure

Definition at line 624 of file OneDFSIData.hpp.

◆ venousPressure()

const Real& venousPressure ( ) const
inline

Get the venous pressure.

Returns
venous pressure.

Definition at line 633 of file OneDFSIData.hpp.

◆ robertsonCorrection()

const Real & robertsonCorrection ( ) const

Get the Robertson correction coefficient (Not tested: maybe wrong in the code)

Returns
Robertson correction coefficient

Definition at line 782 of file OneDFSIData.cpp.

◆ thickness()

const Real& thickness ( const UInt i) const
inline

Get the wall thickness.

Returns
wall thickness

Definition at line 648 of file OneDFSIData.hpp.

◆ friction()

const Real& friction ( ) const
inline

Get the wall friction coefficient $K_r$.

Returns
wall friction coefficient $K_r$

Definition at line 657 of file OneDFSIData.hpp.

◆ area0()

const Real& area0 ( const UInt i) const
inline

Get the reference area $A^0$.

Parameters
inode id
Returns
reference area $A^0$

Definition at line 667 of file OneDFSIData.hpp.

◆ alpha()

const Real& alpha ( const UInt i) const
inline

Get the Coriolis coefficient $\alpha$.

Parameters
inode id
Returns
Coriolis coefficient $\alpha$

Definition at line 677 of file OneDFSIData.hpp.

◆ beta0()

const Real& beta0 ( const UInt i) const
inline

Get the $\beta_0$ coefficient.

Parameters
inode id
Returns
$\beta_0$ coefficient

Definition at line 687 of file OneDFSIData.hpp.

◆ beta1()

const Real& beta1 ( const UInt i) const
inline

Get the $\beta_1$ coefficient.

Parameters
inode id
Returns
$\beta_1$ coefficient

Definition at line 697 of file OneDFSIData.hpp.

◆ dArea0dz()

const Real& dArea0dz ( const UInt i) const
inline

Get $\frac{dA^0}{dz}$.

Parameters
inode id
Returns
$\frac{dA^0}{dz}$

Definition at line 707 of file OneDFSIData.hpp.

◆ dAlphadz()

const Real& dAlphadz ( const UInt i) const
inline

Get $\frac{d\alpha}{dz}$.

Parameters
inode id
Returns
$\frac{d\alpha}{dz}$

Definition at line 717 of file OneDFSIData.hpp.

◆ dBeta0dz()

const Real& dBeta0dz ( const UInt i) const
inline

Get $\frac{d\beta_0}{dz}$.

Parameters
inode id
Returns
$\frac{d\beta_0}{dz}$

Definition at line 727 of file OneDFSIData.hpp.

◆ dBeta1dz()

const Real& dBeta1dz ( const UInt i) const
inline

Get $\frac{d\beta_1}{dz}$.

Parameters
inode id
Returns
$\frac{d\beta_1}{dz}$

Definition at line 737 of file OneDFSIData.hpp.

◆ flux11()

const Real& flux11 ( const UInt i) const
inline

Get the flux coefficient $F_{11}$ (used only in the linear problem)

Parameters
inode id
Returns
$F_{11}$

Definition at line 747 of file OneDFSIData.hpp.

◆ flux12()

const Real& flux12 ( const UInt i) const
inline

Get the flux coefficient $F_{12}$ (used only in the linear problem)

Parameters
inode id
Returns
$F_{12}$

Definition at line 757 of file OneDFSIData.hpp.

◆ flux21()

const Real& flux21 ( const UInt i) const
inline

Get the flux coefficient $F_{21}$ (used only in the linear problem)

Parameters
inode id
Returns
$F_{21}$

Definition at line 767 of file OneDFSIData.hpp.

◆ flux22()

const Real& flux22 ( const UInt i) const
inline

Get the flux coefficient $F_{22}$ (used only in the linear problem)

Parameters
inode id
Returns
$F_{22}$

Definition at line 777 of file OneDFSIData.hpp.

◆ celerity1()

const Real& celerity1 ( const UInt i) const
inline

Get the first eigenvector $\lambda_{1}$ (used only in the linear problem)

Parameters
inode id
Returns
$\lambda_{1}$

Definition at line 787 of file OneDFSIData.hpp.

◆ celerity2()

const Real& celerity2 ( const UInt i) const
inline

Get the second eigenvector $\lambda_{2}$ (used only in the linear problem)

Parameters
inode id
Returns
$\lambda_{2}$

Definition at line 797 of file OneDFSIData.hpp.

◆ leftEigenVector11()

const Real& leftEigenVector11 ( const UInt i) const
inline

Get the left eigenvector coefficient $L_{11}$ (used only in the linear problem)

Parameters
inode id
Returns
$L_{11}$

Definition at line 807 of file OneDFSIData.hpp.

◆ leftEigenVector12()

const Real& leftEigenVector12 ( const UInt i) const
inline

Get the left eigenvector coefficient $L_{12}$ (used only in the linear problem)

Parameters
inode id
Returns
$L_{12}$

Definition at line 817 of file OneDFSIData.hpp.

◆ leftEigenVector21()

const Real& leftEigenVector21 ( const UInt i) const
inline

Get the left eigenvector coefficient $L_{21}$ (used only in the linear problem)

Parameters
inode id
Returns
$L_{21}$

Definition at line 827 of file OneDFSIData.hpp.

◆ leftEigenVector22()

const Real& leftEigenVector22 ( const UInt i) const
inline

Get the left eigenvector coefficient $L_{22}$ (used only in the linear problem)

Parameters
inode id
Returns
$L_{22}$

Definition at line 837 of file OneDFSIData.hpp.

◆ source10()

const Real& source10 ( const UInt i) const
inline

Get the source coefficient $S_{10}$ (used only in the linear problem)

Parameters
inode id
Returns
$S_{10}$

Definition at line 847 of file OneDFSIData.hpp.

◆ source20()

const Real& source20 ( const UInt i) const
inline

Get the source coefficient $S_{20}$ (used only in the linear problem)

Parameters
inode id
Returns
$S_{20}$

Definition at line 857 of file OneDFSIData.hpp.

◆ source11()

const Real& source11 ( const UInt i) const
inline

Get the source coefficient $S_{11}$ (used only in the linear problem)

Parameters
inode id
Returns
$S_{11}$

Definition at line 867 of file OneDFSIData.hpp.

◆ source12()

const Real& source12 ( const UInt i) const
inline

Get the source coefficient $S_{12}$ (used only in the linear problem)

Parameters
inode id
Returns
$S_{12}$

Definition at line 877 of file OneDFSIData.hpp.

◆ source21()

const Real& source21 ( const UInt i) const
inline

Get the source coefficient $S_{21}$ (used only in the linear problem)

Parameters
inode id
Returns
$S_{21}$

Definition at line 887 of file OneDFSIData.hpp.

◆ source22()

const Real& source22 ( const UInt i) const
inline

Get the source coefficient $S_{22}$ (used only in the linear problem)

Parameters
inode id
Returns
$S_{22}$

Definition at line 897 of file OneDFSIData.hpp.

◆ linearInterpolation()

void linearInterpolation ( scalarVector_Type vector,
const GetPot dataFile,
const std::string &  quantity,
const Real defaultValue,
const bool &  isArea = false 
)
private

Compute the linear interpolation of a general quantity.

Very useful for tapering.

Parameters
vectorinterpolated vector
dataFiledata file
quantityquantity
defaultValuedefault value
isAreaflag identifying if the vector is an area (the linear interpolation is done on the radius)

Definition at line 795 of file OneDFSIData.cpp.

◆ computeSpatialDerivatives()

void computeSpatialDerivatives ( )
private

Compute the derivatives of alpha, area0, beta0, and beta1 using centered differences.

Definition at line 828 of file OneDFSIData.cpp.

+ Here is the caller graph for this function:

◆ resetContainers()

void resetContainers ( )
private

Reset all the containers.

Definition at line 840 of file OneDFSIData.cpp.

+ Here is the caller graph for this function:

Field Documentation

◆ M_physicsType

OneDFSI::physicsType_Type M_physicsType
private

Model.

Definition at line 933 of file OneDFSIData.hpp.

◆ M_fluxType

OneDFSI::fluxTerm_Type M_fluxType
private

Definition at line 934 of file OneDFSIData.hpp.

◆ M_sourceType

OneDFSI::sourceTerm_Type M_sourceType
private

Definition at line 935 of file OneDFSIData.hpp.

◆ M_timeDataPtr

timePtr_Type M_timeDataPtr
private

Data containers for time and mesh.

Definition at line 938 of file OneDFSIData.hpp.

◆ M_meshPtr

meshPtr_Type M_meshPtr
private

Definition at line 939 of file OneDFSIData.hpp.

◆ M_viscoelasticWall

bool M_viscoelasticWall
private

Physical Wall Model.

Definition at line 942 of file OneDFSIData.hpp.

◆ M_viscoelasticAngle

Real M_viscoelasticAngle
private

Definition at line 943 of file OneDFSIData.hpp.

◆ M_viscoelasticPeriod

Real M_viscoelasticPeriod
private

Definition at line 944 of file OneDFSIData.hpp.

◆ M_viscoelasticCoefficient

scalarVector_Type M_viscoelasticCoefficient
private

Definition at line 945 of file OneDFSIData.hpp.

◆ M_inertialWall

bool M_inertialWall
private

Definition at line 947 of file OneDFSIData.hpp.

◆ M_densityWall

Real M_densityWall
private

Definition at line 948 of file OneDFSIData.hpp.

◆ M_inertialModulus

Real M_inertialModulus
private

Definition at line 949 of file OneDFSIData.hpp.

◆ M_longitudinalWall

bool M_longitudinalWall
private

Definition at line 950 of file OneDFSIData.hpp.

◆ M_postprocessingDirectory

std::string M_postprocessingDirectory
private

Miscellaneous.

Definition at line 953 of file OneDFSIData.hpp.

◆ M_postprocessingFile

std::string M_postprocessingFile
private

full directory name (including path)

Definition at line 954 of file OneDFSIData.hpp.

◆ M_CFLmax

Real M_CFLmax
private

output file name

Definition at line 955 of file OneDFSIData.hpp.

◆ M_jacobianPerturbationArea

Real M_jacobianPerturbationArea
private

Jacobian perturbation.

Definition at line 958 of file OneDFSIData.hpp.

◆ M_jacobianPerturbationFlowRate

Real M_jacobianPerturbationFlowRate
private

Definition at line 959 of file OneDFSIData.hpp.

◆ M_jacobianPerturbationStress

Real M_jacobianPerturbationStress
private

Definition at line 960 of file OneDFSIData.hpp.

◆ M_computeCoefficients

bool M_computeCoefficients
private

Physical Parameters.

Definition at line 963 of file OneDFSIData.hpp.

◆ M_powerLawCoefficient

Int M_powerLawCoefficient
private

Definition at line 964 of file OneDFSIData.hpp.

◆ M_density

Real M_density
private

Definition at line 966 of file OneDFSIData.hpp.

◆ M_viscosity

Real M_viscosity
private

Definition at line 967 of file OneDFSIData.hpp.

◆ M_thickVessel

bool M_thickVessel
private

Definition at line 969 of file OneDFSIData.hpp.

◆ M_young

Real M_young
private

Definition at line 970 of file OneDFSIData.hpp.

◆ M_poisson

Real M_poisson
private

Definition at line 971 of file OneDFSIData.hpp.

◆ M_externalPressure

Real M_externalPressure
private

Definition at line 973 of file OneDFSIData.hpp.

◆ M_venousPressure

Real M_venousPressure
private

Definition at line 974 of file OneDFSIData.hpp.

◆ M_robertsonCorrection

Real M_robertsonCorrection
private

Definition at line 975 of file OneDFSIData.hpp.

◆ M_thickness

scalarVector_Type M_thickness
private

Definition at line 977 of file OneDFSIData.hpp.

◆ M_friction

Real M_friction
private

Definition at line 978 of file OneDFSIData.hpp.

◆ M_area0

scalarVector_Type M_area0
private

Definition at line 980 of file OneDFSIData.hpp.

◆ M_alpha

scalarVector_Type M_alpha
private

Definition at line 981 of file OneDFSIData.hpp.

◆ M_beta0

scalarVector_Type M_beta0
private

Definition at line 982 of file OneDFSIData.hpp.

◆ M_beta1

scalarVector_Type M_beta1
private

Definition at line 983 of file OneDFSIData.hpp.

◆ M_dArea0dz

scalarVector_Type M_dArea0dz
private

Derivatives of main coefficients.

Definition at line 986 of file OneDFSIData.hpp.

◆ M_dAlphadz

scalarVector_Type M_dAlphadz
private

Definition at line 987 of file OneDFSIData.hpp.

◆ M_dBeta0dz

scalarVector_Type M_dBeta0dz
private

Definition at line 988 of file OneDFSIData.hpp.

◆ M_dBeta1dz

scalarVector_Type M_dBeta1dz
private

Definition at line 989 of file OneDFSIData.hpp.

◆ M_flux11

scalarVector_Type M_flux11
private

Flux matrix.

Definition at line 992 of file OneDFSIData.hpp.

◆ M_flux12

scalarVector_Type M_flux12
private

Definition at line 993 of file OneDFSIData.hpp.

◆ M_flux21

scalarVector_Type M_flux21
private

Definition at line 994 of file OneDFSIData.hpp.

◆ M_flux22

scalarVector_Type M_flux22
private

Definition at line 995 of file OneDFSIData.hpp.

◆ M_celerity1

scalarVector_Type M_celerity1
private

Celerities of the linear problem (eigenvalues of the flux matrix)

Definition at line 998 of file OneDFSIData.hpp.

◆ M_celerity2

scalarVector_Type M_celerity2
private

Definition at line 999 of file OneDFSIData.hpp.

◆ M_celerity1LeftEigenvector1

scalarVector_Type M_celerity1LeftEigenvector1
private

Eigenvector for first and second eigenvalue.

Definition at line 1002 of file OneDFSIData.hpp.

◆ M_celerity1LeftEigenvector2

scalarVector_Type M_celerity1LeftEigenvector2
private

Definition at line 1003 of file OneDFSIData.hpp.

◆ M_celerity2LeftEigenvector1

scalarVector_Type M_celerity2LeftEigenvector1
private

Definition at line 1004 of file OneDFSIData.hpp.

◆ M_celerity2LeftEigenvector2

scalarVector_Type M_celerity2LeftEigenvector2
private

Definition at line 1005 of file OneDFSIData.hpp.

◆ M_source10

scalarVector_Type M_source10
private

Source matrix.

Definition at line 1008 of file OneDFSIData.hpp.

◆ M_source20

scalarVector_Type M_source20
private

Definition at line 1009 of file OneDFSIData.hpp.

◆ M_source11

scalarVector_Type M_source11
private

Definition at line 1010 of file OneDFSIData.hpp.

◆ M_source12

scalarVector_Type M_source12
private

Definition at line 1011 of file OneDFSIData.hpp.

◆ M_source21

scalarVector_Type M_source21
private

Definition at line 1012 of file OneDFSIData.hpp.

◆ M_source22

scalarVector_Type M_source22
private

Definition at line 1013 of file OneDFSIData.hpp.


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