#include <VentilationProblem.hpp>

Public Member Functions | |
| VentilationProblem () | |
| VentilationProblem (const std::string &rMeshDirFilePath, unsigned rootIndex=0u) | |
| ~VentilationProblem () | |
| void | SetOutflowPressure (double pressure) |
| void | SetConstantInflowPressures (double pressure) |
| void | SetConstantInflowFluxes (double flux) |
| void | SetPressureAtBoundaryNode (const Node< 3 > &rNode, double pressure) |
| double | GetPressureAtBoundaryNode (const Node< 3 > &rNode) |
| double | GetFluxAtOutflow () |
| void | SetFluxAtBoundaryNode (const Node< 3 > &rNode, double flux) |
| void | Solve () |
| void | GetSolutionAsFluxesAndPressures (std::vector< double > &rFluxesOnEdges, std::vector< double > &rPressuresOnNodes) |
| void | AddDataToVtk (VtkMeshWriter< 1, 3 > &rVtkWriter, const std::string &rSuffix) |
| void | WriteVtk (const std::string &rDirName, const std::string &rFileBaseName) |
| void | SetRadiusOnEdge (bool isOnEdges=true) |
| TetrahedralMesh< 1, 3 > & | rGetMesh () |
| void | Solve (TimeStepper &rTimeStepper, void(*pBoundaryConditionFunction)(VentilationProblem *, TimeStepper &rTimeStepper, const Node< 3 > &), const std::string &rDirName, const std::string &rFileBaseName) |
| void | SolveProblemFromFile (const std::string &rInFilePath, const std::string &rOutFileDir, const std::string &rOutFileName) |
| double | GetViscosity () const |
| void | SetViscosity (double viscosity) |
| double | GetDensity () const |
| void | SetDensity (double density) |
| void | SetDynamicResistance (bool dynamicResistance=true) |
| Swan2012AcinarUnit * | GetAcinus (const Node< 3 > &rNode) |
Private Member Functions | |
| void | SolveDirectFromFlux () |
| void | SetupIterativeSolver () |
| void | SolveIterativelyFromPressure () |
| double | CalculateResistance (Element< 1, 3 > &rElement, bool usePedley=false, double flux=DBL_MAX) |
Private Attributes | |
| TetrahedralMesh< 1, 3 > | mMesh |
| unsigned | mOutletNodeIndex |
| bool | mDynamicResistance |
| bool | mRadiusOnEdge |
| double | mViscosity |
| double | mDensity |
| std::map< unsigned, Swan2012AcinarUnit * > | mAcinarUnits |
| std::vector< double > | mFlux |
| std::vector< double > | mPressure |
| std::map< unsigned, double > | mPressureCondition |
| bool | mFluxGivenAtInflow |
| Mat | mTerminalInteractionMatrix |
| std::map< unsigned, unsigned > | mTerminalToNodeIndex |
| std::map< unsigned, unsigned > | mTerminalToEdgeIndex |
| Vec | mTerminalFluxChangeVector |
| Vec | mTerminalPressureChangeVector |
| KSP | mTerminalKspSolver |
A class for solving one-dimensional flow in pipe problems on branching trees.
At graph edges: each pipe models Poiseuille flow (flux is linearly proportional to the pressure drop) At graph nodes: the flux is balanced so that mass is conserved.
Works in 3D <1,3> Current functionality: pressure boundary conditions are set on each of the boundary nodes Solves for pressure at internal nodes and flux on edges
Definition at line 56 of file VentilationProblem.hpp.
| VentilationProblem::VentilationProblem | ( | ) |
Default constructor Attempts to read all parameters from a hard-coded file Loads a mesh from file(s) Identifies the outlet node (a.k.a root of tree or the mouth end) A check is made that it is a boundary node. We could also check that on trees with more than one bifurcation there are no boundary nodes in its 2nd generation successors. Creates a linear system of the appropriate size to match the mesh
| VentilationProblem::VentilationProblem | ( | const std::string & | rMeshDirFilePath, | |
| unsigned | rootIndex = 0u | |||
| ) |
Constructor Loads a mesh from file(s) Identifies the outlet node (a.k.a root of tree or the mouth end) A check is made that it is a boundary node. We could also check that on trees with more than one bifurcation there are no boundary nodes in its 2nd generation successors. Creates a linear system of the appropriate size to match the mesh
| rMeshDirFilePath | the path and root name of the .node and .edge files for the mesh | |
| rootIndex | the global index of the root/outlet node in the mesh (defaults to node zero). |
Definition at line 41 of file VentilationProblem.cpp.
References TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ConstructFromMeshReader(), EXCEPTION, AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryNodeIteratorBegin(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryNodeIteratorEnd(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNode(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumElements(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetWeightedDirectionForElement(), mAcinarUnits, mFlux, mMesh, mOutletNodeIndex, mPressure, mViscosity, Swan2012AcinarUnit::SetAirwayPressure(), Swan2012AcinarUnit::SetPleuralPressure(), Swan2012AcinarUnit::SetStretchRatio(), Swan2012AcinarUnit::SetTerminalBronchioleResistance(), and Swan2012AcinarUnit::SetUndeformedVolume().
| VentilationProblem::~VentilationProblem | ( | ) |
Destructor destroys the linear system
Definition at line 106 of file VentilationProblem.cpp.
References PetscTools::Destroy(), mAcinarUnits, mTerminalFluxChangeVector, mTerminalInteractionMatrix, mTerminalKspSolver, mTerminalPressureChangeVector, and PETSC_DESTROY_PARAM.
| void VentilationProblem::AddDataToVtk | ( | VtkMeshWriter< 1, 3 > & | rVtkWriter, | |
| const std::string & | rSuffix | |||
| ) |
Add flux and pressure data to a VtkMeshWriter.
| rVtkWriter | the mesh writer ready for the data | |
| rSuffix | Suffix with which to annotate e.g. pressure_001 |
Definition at line 522 of file VentilationProblem.cpp.
References VtkMeshWriter< ELEMENT_DIM, SPACE_DIM >::AddCellData(), VtkMeshWriter< ELEMENT_DIM, SPACE_DIM >::AddPointData(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryNodeIteratorBegin(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryNodeIteratorEnd(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), GetSolutionAsFluxesAndPressures(), mAcinarUnits, mMesh, and mOutletNodeIndex.
Referenced by Solve(), and WriteVtk().
| double VentilationProblem::CalculateResistance | ( | Element< 1, 3 > & | rElement, | |
| bool | usePedley = false, |
|||
| double | flux = DBL_MAX | |||
| ) | [private] |
Get the resistance of an edge. This defaults to Poiseuille resistance (in which only the geometry is used. Otherwise, Pedley's correction is calculated, which requires a flux to be given.
| rElement | The edge on which to perform this calculation | |
| usePedley | Whether to add Pedley's increasing correction term. Here the resistance increases which the sqrt of Reynold's number (dependent on flux). | |
| flux | The flux in the edge (used for Pedley correction). |
Definition at line 351 of file VentilationProblem.cpp.
References AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetAttribute(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNode(), TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetWeightedDirectionForElement(), mDensity, mMesh, mRadiusOnEdge, and mViscosity.
Referenced by SetupIterativeSolver(), and SolveDirectFromFlux().
| Swan2012AcinarUnit* VentilationProblem::GetAcinus | ( | const Node< 3 > & | rNode | ) | [inline] |
| rNode | The node to get the acinus at. Must be a boundary node! |
Definition at line 339 of file VentilationProblem.hpp.
References Node< SPACE_DIM >::GetIndex(), and mAcinarUnits.
| double VentilationProblem::GetDensity | ( | ) | const [inline] |
Definition at line 314 of file VentilationProblem.hpp.
References mDensity.
| double VentilationProblem::GetFluxAtOutflow | ( | ) |
Gets the most recent flux at the outflow/outlet (mouth)
Definition at line 457 of file VentilationProblem.cpp.
References mFlux, and mOutletNodeIndex.
Gets the most recent pressure at a boundary node
| rNode | The node to get the pressure for. |
Definition at line 447 of file VentilationProblem.cpp.
References Node< SPACE_DIM >::GetIndex(), and mPressure.
| void VentilationProblem::GetSolutionAsFluxesAndPressures | ( | std::vector< double > & | rFluxesOnEdges, | |
| std::vector< double > & | rPressuresOnNodes | |||
| ) |
| rFluxesOnEdges | The fluxes ordered by edge index (this vector is resized) | |
| rPressuresOnNodes | The pressures ordered by node index (this vector is resized) |
Definition at line 493 of file VentilationProblem.cpp.
References mFlux, and mPressure.
Referenced by AddDataToVtk().
| double VentilationProblem::GetViscosity | ( | ) | const [inline] |
Definition at line 298 of file VentilationProblem.hpp.
References mViscosity.
| TetrahedralMesh< 1, 3 > & VentilationProblem::rGetMesh | ( | ) |
Definition at line 506 of file VentilationProblem.cpp.
References mMesh.
| void VentilationProblem::SetConstantInflowFluxes | ( | double | flux | ) |
Sets the flux at each inflow/leaf-edge of the tree Flux is "volumetric flow rate"
| flux | The flux value in (mm^3)/s |
Definition at line 422 of file VentilationProblem.cpp.
References AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryNodeIteratorBegin(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryNodeIteratorEnd(), mMesh, mOutletNodeIndex, and SetFluxAtBoundaryNode().
| void VentilationProblem::SetConstantInflowPressures | ( | double | pressure | ) |
Sets the pressure at each inflow/leaf of the tree
| pressure | The pressure value in Pascals |
Definition at line 408 of file VentilationProblem.cpp.
References AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryNodeIteratorBegin(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryNodeIteratorEnd(), mMesh, mOutletNodeIndex, and SetPressureAtBoundaryNode().
Referenced by SolveProblemFromFile().
| void VentilationProblem::SetDensity | ( | double | density | ) | [inline] |
| density | the density in kg/(m^3) |
Definition at line 322 of file VentilationProblem.hpp.
References mDensity.
Referenced by SolveProblemFromFile().
| void VentilationProblem::SetDynamicResistance | ( | bool | dynamicResistance = true |
) | [inline] |
| dynamicResistance |
Definition at line 330 of file VentilationProblem.hpp.
References mDynamicResistance.
Sets a Dirichlet flux boundary condition for a given node.
The given boundary condition will be applied at the next time step and persist through time unless overwritten.
| rNode | The node to set the boundary condition for | |
| flux | The flux boundary condition in (mm^3)/s |
Definition at line 462 of file VentilationProblem.cpp.
References Node< SPACE_DIM >::ContainingElementsBegin(), EXCEPTION, Node< SPACE_DIM >::IsBoundaryNode(), mFlux, and mFluxGivenAtInflow.
Referenced by SetConstantInflowFluxes().
| void VentilationProblem::SetOutflowPressure | ( | double | pressure | ) |
Sets the pressure at the outflow/root of the tree
| pressure | The pressure value in Pascals |
Definition at line 402 of file VentilationProblem.cpp.
References AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNode(), mMesh, mOutletNodeIndex, mPressure, and SetPressureAtBoundaryNode().
Referenced by SolveProblemFromFile().
Sets a Dirichlet pressure boundary condition for a given node.
The given boundary condition will be applied at the next time step and persist through time unless overwritten.
| rNode | The node to set the boundary condition for | |
| pressure | The pressure boundary condition in Pascals |
Definition at line 435 of file VentilationProblem.cpp.
References EXCEPTION, Node< SPACE_DIM >::GetIndex(), Node< SPACE_DIM >::IsBoundaryNode(), mFluxGivenAtInflow, and mPressureCondition.
Referenced by SetConstantInflowPressures(), and SetOutflowPressure().
| void VentilationProblem::SetRadiusOnEdge | ( | bool | isOnEdges = true |
) |
Used to set mRadiusOnEdge flag. This is false by default in the constructor (conic pipes with radius defined at nodes). When true pipes are cylindrical.
| isOnEdges | The new value of mRadiusOnEdge |
Definition at line 501 of file VentilationProblem.cpp.
References mRadiusOnEdge.
| void VentilationProblem::SetupIterativeSolver | ( | ) | [private] |
Set up the PETSc machinery for solving the iterative problem: given pressure conditions at the terminals guess and refine flux conditions which match them.
This creates and fills a PETSc Mat. It also creates two PETSc Vecs and a KSP solver.
Definition at line 162 of file VentilationProblem.cpp.
References CalculateResistance(), Node< SPACE_DIM >::ContainingElementsBegin(), Node< SPACE_DIM >::ContainingElementsEnd(), PetscMatTools::Finalise(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryNodeIteratorBegin(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryNodeIteratorEnd(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorBegin(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorEnd(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNode(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumBoundaryNodes(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumElements(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), Node< SPACE_DIM >::IsBoundaryNode(), mMesh, mOutletNodeIndex, mTerminalFluxChangeVector, mTerminalInteractionMatrix, mTerminalKspSolver, mTerminalPressureChangeVector, mTerminalToEdgeIndex, mTerminalToNodeIndex, and PetscMatTools::SetOption().
Referenced by SolveIterativelyFromPressure().
| void VentilationProblem::SetViscosity | ( | double | viscosity | ) | [inline] |
| viscosity | the viscosity in kg/(mm*sec) |
Definition at line 306 of file VentilationProblem.hpp.
References mViscosity.
Referenced by SolveProblemFromFile().
| void VentilationProblem::Solve | ( | TimeStepper & | rTimeStepper, | |
| void(*)(VentilationProblem *, TimeStepper &rTimeStepper, const Node< 3 > &) | pBoundaryConditionFunction, | |||
| const std::string & | rDirName, | |||
| const std::string & | rFileBaseName | |||
| ) |
Assemble the linear system by writing in * flux balance at the nodes * Poiseuille flow in the edges
Solve the linear system repeatedly
| rTimeStepper | The start, end and time-step | |
| pBoundaryConditionFunction | setter | |
| rDirName | A directory name relative to CHASTE_TEST_OUTPUT. | |
| rFileBaseName | The base name of the new VTK file. |
Definition at line 550 of file VentilationProblem.cpp.
References AddDataToVtk(), TimeStepper::AdvanceOneTimeStep(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryNodeIteratorBegin(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryNodeIteratorEnd(), TimeStepper::GetTotalTimeStepsTaken(), TimeStepper::IsTimeAtEnd(), mMesh, mOutletNodeIndex, Solve(), and VtkMeshWriter< ELEMENT_DIM, SPACE_DIM >::WriteFilesUsingMesh().
| void VentilationProblem::Solve | ( | ) |
Solve the linear system either * directly from fluxes * iteratively from pressures
Definition at line 479 of file VentilationProblem.cpp.
References mFluxGivenAtInflow, SolveDirectFromFlux(), and SolveIterativelyFromPressure().
Referenced by Solve(), and SolveProblemFromFile().
| void VentilationProblem::SolveDirectFromFlux | ( | ) | [private] |
Use flux boundary conditions at leaves (and pressure condition at root) to perform a direct solve. This involves * solving directly for parent flux up the tree (using flux balance at each node) * solving for child pressure (using Poiseuille or Pedley resistance) down the tree
Definition at line 122 of file VentilationProblem.cpp.
References CalculateResistance(), Node< SPACE_DIM >::ContainingElementsBegin(), Node< SPACE_DIM >::ContainingElementsEnd(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorBegin(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElementIteratorEnd(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNode(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), Node< SPACE_DIM >::IsBoundaryNode(), mDynamicResistance, mFlux, mMesh, and mPressure.
Referenced by Solve(), and SolveIterativelyFromPressure().
| void VentilationProblem::SolveIterativelyFromPressure | ( | ) | [private] |
Use pressure boundary conditions at leaves (and pressure condition at root) to perform convert to flux boundary conditions (assuming Poiseuille flow) and then perform a direct solve with SolveDirectFromFlux(). Note that pressure to flux conversion requires accumulation of the resistance down the tree since weighted_sum(resistance)*flux = delta pressure.
Note in the mDynamicResistance case we ignore dynamic resistance when setting up the weighted sum of resistances. However, since we use dynamic resistance in SolveDirectFromFlux() the solution will converge to the one which uses dynamic resistance, despite the flux corrections being slightly off.
Definition at line 237 of file VentilationProblem.cpp.
References AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumBoundaryNodes(), mFlux, mMesh, mOutletNodeIndex, mPressure, mPressureCondition, mTerminalFluxChangeVector, mTerminalInteractionMatrix, mTerminalKspSolver, mTerminalPressureChangeVector, mTerminalToEdgeIndex, mTerminalToNodeIndex, NEVER_REACHED, SetupIterativeSolver(), and SolveDirectFromFlux().
Referenced by Solve().
| void VentilationProblem::SolveProblemFromFile | ( | const std::string & | rInFilePath, | |
| const std::string & | rOutFileDir, | |||
| const std::string & | rOutFileName | |||
| ) |
Read a problem definition from a file and use then solve that problem
| rInFilePath | Path to file which contains the problem definition | |
| rOutFileDir | Path to folder for output (relative to CHASTE_TEST_OUTPUT) | |
| rOutFileName | Name for VTK output |
Definition at line 596 of file VentilationProblem.cpp.
References EXCEPTION, SetConstantInflowPressures(), SetDensity(), SetOutflowPressure(), SetViscosity(), Solve(), and WriteVtk().
| void VentilationProblem::WriteVtk | ( | const std::string & | rDirName, | |
| const std::string & | rFileBaseName | |||
| ) |
Output the solution to a Vtk file
| rDirName | A directory name relative to CHASTE_TEST_OUTPUT. | |
| rFileBaseName | The base name of the new VTK file. |
Definition at line 514 of file VentilationProblem.cpp.
References AddDataToVtk(), mMesh, and VtkMeshWriter< ELEMENT_DIM, SPACE_DIM >::WriteFilesUsingMesh().
Referenced by SolveProblemFromFile().
std::map<unsigned, Swan2012AcinarUnit*> VentilationProblem::mAcinarUnits [private] |
One acinar unit for each terminal node.
Definition at line 85 of file VentilationProblem.hpp.
Referenced by AddDataToVtk(), GetAcinus(), VentilationProblem(), and ~VentilationProblem().
double VentilationProblem::mDensity [private] |
Density in kg/(mm^3). rho (for dry air) ~ 1.2041 kg/m^3 = 1.2e-9 kg/mm^3 Default to Swan (page 224) rho = 1.51e-6 g/mm^3 <-- USED HERE = 1.51e-9 kg/mm^3 = 1.51e-6 kg/(m s^2) * s^2 / (mm)^2
Definition at line 83 of file VentilationProblem.hpp.
Referenced by CalculateResistance(), GetDensity(), and SetDensity().
bool VentilationProblem::mDynamicResistance [private] |
Use dynamic (flux related) resistance and a nonlinear solver
Definition at line 61 of file VentilationProblem.hpp.
Referenced by SetDynamicResistance(), and SolveDirectFromFlux().
std::vector<double> VentilationProblem::mFlux [private] |
Used to hold the flux solution (and boundary conditions) in edge index ordering
Definition at line 86 of file VentilationProblem.hpp.
Referenced by GetFluxAtOutflow(), GetSolutionAsFluxesAndPressures(), SetFluxAtBoundaryNode(), SolveDirectFromFlux(), SolveIterativelyFromPressure(), and VentilationProblem().
bool VentilationProblem::mFluxGivenAtInflow [private] |
Used to switch solution methods. If the flux is given at the boundary then the entire system can be solved directly by back substitution.
Definition at line 89 of file VentilationProblem.hpp.
Referenced by SetFluxAtBoundaryNode(), SetPressureAtBoundaryNode(), and Solve().
TetrahedralMesh<1,3> VentilationProblem::mMesh [private] |
The 1d in 3d branching tree mesh
Definition at line 59 of file VentilationProblem.hpp.
Referenced by AddDataToVtk(), CalculateResistance(), rGetMesh(), SetConstantInflowFluxes(), SetConstantInflowPressures(), SetOutflowPressure(), SetupIterativeSolver(), Solve(), SolveDirectFromFlux(), SolveIterativelyFromPressure(), VentilationProblem(), and WriteVtk().
unsigned VentilationProblem::mOutletNodeIndex [private] |
The outlet node is the root of the branching tree structure
Definition at line 60 of file VentilationProblem.hpp.
Referenced by AddDataToVtk(), GetFluxAtOutflow(), SetConstantInflowFluxes(), SetConstantInflowPressures(), SetOutflowPressure(), SetupIterativeSolver(), Solve(), SolveIterativelyFromPressure(), and VentilationProblem().
std::vector<double> VentilationProblem::mPressure [private] |
Used to hold the pressure solution (and outlet boundary pressure) in node index ordering).
Definition at line 87 of file VentilationProblem.hpp.
Referenced by GetPressureAtBoundaryNode(), GetSolutionAsFluxesAndPressures(), SetOutflowPressure(), SolveDirectFromFlux(), SolveIterativelyFromPressure(), and VentilationProblem().
std::map<unsigned, double> VentilationProblem::mPressureCondition [private] |
Pressure boundary conditions at terminal nodes.
Definition at line 88 of file VentilationProblem.hpp.
Referenced by SetPressureAtBoundaryNode(), and SolveIterativelyFromPressure().
bool VentilationProblem::mRadiusOnEdge [private] |
False by default (conical pipes with radius defined at nodes). When true pipes are cylindrical.
Definition at line 62 of file VentilationProblem.hpp.
Referenced by CalculateResistance(), and SetRadiusOnEdge().
Used as the output of the mTerminalInteractionMatrix terminal pressure to flux solver
Definition at line 103 of file VentilationProblem.hpp.
Referenced by SetupIterativeSolver(), SolveIterativelyFromPressure(), and ~VentilationProblem().
The symmetric matrix is an estimate of the dense matrix system which determines how flux changes at terminal nodes are reflected in pressure changes at terminals in the form P_{terminal} = A Q_{terminal}. Each entry in the dense matrix (for Poiseuille flow) is the sum of the resistances of all pipes which are common ancestors of a pair of terminals. In order to iteratively match pressure conditions we must invert this equation.
Definition at line 99 of file VentilationProblem.hpp.
Referenced by SetupIterativeSolver(), SolveIterativelyFromPressure(), and ~VentilationProblem().
KSP VentilationProblem::mTerminalKspSolver [private] |
The linear solver for the mTerminalInteractionMatrix terminal pressure to flux solver
Definition at line 105 of file VentilationProblem.hpp.
Referenced by SetupIterativeSolver(), SolveIterativelyFromPressure(), and ~VentilationProblem().
Used as the input of the mTerminalInteractionMatrix terminal pressure to flux solver
Definition at line 104 of file VentilationProblem.hpp.
Referenced by SetupIterativeSolver(), SolveIterativelyFromPressure(), and ~VentilationProblem().
std::map<unsigned, unsigned> VentilationProblem::mTerminalToEdgeIndex [private] |
A mapping from the indexing scheme used in the mTerminalInteractionMatrix to the full mesh element indexing
Definition at line 101 of file VentilationProblem.hpp.
Referenced by SetupIterativeSolver(), and SolveIterativelyFromPressure().
std::map<unsigned, unsigned> VentilationProblem::mTerminalToNodeIndex [private] |
A mapping from the indexing scheme used in the mTerminalInteractionMatrix to the full mesh node indexing
Definition at line 100 of file VentilationProblem.hpp.
Referenced by SetupIterativeSolver(), and SolveIterativelyFromPressure().
double VentilationProblem::mViscosity [private] |
(Dynamic) viscosity in kg/(mm*second). Default to value from Swan et al. 2012. 10.1016/j.jtbi.2012.01.042 (page 224) mu = 1.92e-5 Pa*s <-- USED HERE = 1.92e-5 kg/(m*s) = 1.92e-8 kg/(mm*s) or kPa.s Consider http://en.wikipedia.org/wiki/Viscosity#Air which gives mu = 1.81x10^(-5) kg/(m*s) -- 1.86x10^(-5) kg/(m*s)
Definition at line 72 of file VentilationProblem.hpp.
Referenced by CalculateResistance(), GetViscosity(), SetViscosity(), and VentilationProblem().
1.6.2