v0.10.0
Classes | Public Types | Public Member Functions | Public Attributes | List of all members
MixTransport::UnsaturatedFlowElement Struct Reference

Implementation of operators, problem and finite elements for unsaturated flow. More...

#include <users_modules/tutorials/low/mix_transport/src/UnsaturatedFlow.hpp>

Inheritance diagram for MixTransport::UnsaturatedFlowElement:
[legend]
Collaboration diagram for MixTransport::UnsaturatedFlowElement:
[legend]

Classes

struct  BcData
 Class storing information about boundary condition. More...
 
struct  FaceRule
 Set integration rule to boundary elements. More...
 
struct  MonitorPostProc
 
struct  OpDivTauU_HdivL2
 
struct  OpEvaluateInitiallHead
 
struct  OpIntegrateFluxes
 
struct  OpPostProcMaterial
 
struct  OpResidualFlux
 Assemble flux residual. More...
 
struct  OpResidualMass
 
struct  OpRhsBcOnValues
 Evaluate boundary condition at the boundary. More...
 
struct  OpTauDotSigma_HdivHdiv
 
struct  OpVDivSigma_L2Hdiv
 
struct  OpVU_L2L2
 
struct  postProcessVol
 Post proces method for volume element Assemble vectors and matrices and apply essential boundary conditions. More...
 
struct  preProcessVol
 Pre-peprocessing Set head pressute rate and get inital essential boundary conditions. More...
 
struct  VolRule
 Set integration rule to volume elements. More...
 

Public Types

typedef std::map< int, boost::shared_ptr< GenericMaterial > > MaterialsDoubleMap
 
typedef map< int, boost::shared_ptr< BcData > > BcMap
 

Public Member Functions

 UnsaturatedFlowElement (MoFEM::Interface &m_field)
 
 ~UnsaturatedFlowElement ()
 
virtual MoFEMErrorCode getMaterial (const EntityHandle ent, int &block_id) const
 For given element handle get material block Id. More...
 
MoFEMErrorCode getBcOnValues (const EntityHandle ent, const int gg, const double x, const double y, const double z, double &value)
 Get value on boundary. More...
 
MoFEMErrorCode getBcOnFluxes (const EntityHandle ent, const double x, const double y, const double z, double &flux)
 essential (Neumann) boundary condition (set fluxes) More...
 
MoFEMErrorCode addFields (const std::string &values, const std::string &fluxes, const int order)
 add fields More...
 
MoFEMErrorCode addFiniteElements (const std::string &fluxes_name, const std::string &values_name)
 add finite elements More...
 
MoFEMErrorCode buildProblem (Range zero_flux_ents, BitRefLevel ref_level=BitRefLevel().set(0))
 Build problem. More...
 
MoFEMErrorCode setFiniteElements (ForcesAndSourcesCore::RuleHookFun vol_rule=VolRule(), ForcesAndSourcesCore::RuleHookFun face_rule=FaceRule())
 Create finite element instances. More...
 
MoFEMErrorCode createMatrices ()
 Create vectors and matrices. More...
 
MoFEMErrorCode destroyMatrices ()
 Delete matrices and vector when no longer needed. More...
 
MoFEMErrorCode calculateEssentialBc ()
 Calculate boundary conditions for fluxes. More...
 
MoFEMErrorCode calculateInitialPc ()
 Calculate inital pressure head distribution. More...
 
MoFEMErrorCode solveProblem (bool set_initial_pc=true)
 solve problem More...
 
- Public Member Functions inherited from MixTransport::MixTransportElement
 MixTransportElement (MoFEM::Interface &m_field)
 construction of this data structure More...
 
virtual ~MixTransportElement ()
 destructor More...
 
MoFEMErrorCode getDirichletBCIndices (IS *is)
 get dof indices where essential boundary conditions are applied More...
 
virtual MoFEMErrorCode getSource (const EntityHandle ent, const double x, const double y, const double z, double &flux)
 set source term More...
 
virtual MoFEMErrorCode getResistivity (const EntityHandle ent, const double x, const double y, const double z, MatrixDouble3by3 &inv_k)
 natural (Dirichlet) boundary conditions (set values) More...
 
MoFEMErrorCode addFields (const std::string &values, const std::string &fluxes, const int order)
 Add fields to database. More...
 
MoFEMErrorCode addFiniteElements (const std::string &fluxes_name, const std::string &values_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS")
 add finite elements More...
 
MoFEMErrorCode buildProblem (BitRefLevel &ref_level)
 Build problem. More...
 
MoFEMErrorCode postProc (const string out_file)
 Post process results. More...
 
MoFEMErrorCode createMatrices ()
 create matrices More...
 
MoFEMErrorCode solveLinearProblem ()
 solve problem More...
 
MoFEMErrorCode calculateResidual ()
 calculate residual More...
 
MoFEMErrorCode evaluateError ()
 Calculate error on elements. More...
 
MoFEMErrorCode destroyMatrices ()
 destroy matrices More...
 

Public Attributes

DM dM
 Discrete manager for unsaturated flow problem. More...
 
MaterialsDoubleMap dMatMap
 materials database More...
 
BcMap bcValueMap
 Store boundary condition for head capillary pressure. More...
 
EntityHandle lastEvalBcValEnt
 
int lastEvalBcBlockValId
 
BcMap bcFluxMap
 
EntityHandle lastEvalBcFluxEnt
 
int lastEvalBcBlockFluxId
 
boost::shared_ptr< ForcesAndSourcesCore > feFaceBc
 Elemnet to calculate essential bc. More...
 
boost::shared_ptr< ForcesAndSourcesCore > feFaceRhs
 Face element apply natural bc. More...
 
boost::shared_ptr< ForcesAndSourcesCore > feVolInitialPc
 Calculate inital boundary conditions. More...
 
boost::shared_ptr< ForcesAndSourcesCore > feVolRhs
 Assemble residual vector. More...
 
boost::shared_ptr< ForcesAndSourcesCore > feVolLhs
 Assemble tangent matrix. More...
 
boost::shared_ptr< MethodForForceScalingscaleMethodFlux
 Method scaling fluxes. More...
 
boost::shared_ptr< MethodForForceScalingscaleMethodValue
 Method scaling values. More...
 
boost::shared_ptr< FEMethodtsMonitor
 
boost::shared_ptr< VectorDouble > headRateAtGaussPts
 Vector keeps head rate. More...
 
std::vector< int > bcVecIds
 
VectorDouble bcVecVals
 
VectorDouble vecValsOnBc
 
Vec D1
 Vector with inital head capillary pressure. More...
 
Vec ghostFlux
 Ghost Vector of integrated fluxes. More...
 
- Public Attributes inherited from MixTransport::MixTransportElement
MoFEM::InterfacemField
 
MyVolumeFE feVol
 Instance of volume element. More...
 
MyTriFE feTri
 Instance of surface element. More...
 
VectorDouble valuesAtGaussPts
 values at integration points on element More...
 
MatrixDouble valuesGradientAtGaussPts
 gradients at integration points on element More...
 
VectorDouble divergenceAtGaussPts
 divergence at integration points on element More...
 
MatrixDouble fluxesAtGaussPts
 fluxes at integration points on element More...
 
set< int > bcIndices
 
std::map< int, BlockDatasetOfBlocks
 maps block set id with appropriate BlockData More...
 
Vec D
 
Vec D0
 
Vec F
 
Mat Aij
 
map< double, EntityHandleerrorMap
 
double sumErrorFlux
 
double sumErrorDiv
 
double sumErrorJump
 

Detailed Description

Implementation of operators, problem and finite elements for unsaturated flow.

Examples
unsaturated_transport.cpp.

Definition at line 115 of file UnsaturatedFlow.hpp.

Member Typedef Documentation

◆ BcMap

typedef map<int, boost::shared_ptr<BcData> > MixTransport::UnsaturatedFlowElement::BcMap
Examples
UnsaturatedFlow.hpp.

Definition at line 165 of file UnsaturatedFlow.hpp.

◆ MaterialsDoubleMap

Examples
UnsaturatedFlow.hpp.

Definition at line 131 of file UnsaturatedFlow.hpp.

Constructor & Destructor Documentation

◆ UnsaturatedFlowElement()

MixTransport::UnsaturatedFlowElement::UnsaturatedFlowElement ( MoFEM::Interface m_field)
Examples
UnsaturatedFlow.hpp.

Definition at line 119 of file UnsaturatedFlow.hpp.

120  : MixTransportElement(m_field), dM(PETSC_NULL), lastEvalBcValEnt(0),
122  lastEvalBcBlockFluxId(-1) {}
MixTransportElement(MoFEM::Interface &m_field)
construction of this data structure
DM dM
Discrete manager for unsaturated flow problem.

◆ ~UnsaturatedFlowElement()

MixTransport::UnsaturatedFlowElement::~UnsaturatedFlowElement ( )
Examples
UnsaturatedFlow.hpp.

Definition at line 124 of file UnsaturatedFlow.hpp.

124  {
125  if (dM != PETSC_NULL) {
126  CHKERR DMDestroy(&dM);
127  CHKERRABORT(PETSC_COMM_WORLD, ierr);
128  }
129  }
#define CHKERR
Inline error check.
Definition: definitions.h:604
DM dM
Discrete manager for unsaturated flow problem.

Member Function Documentation

◆ addFields()

MoFEMErrorCode MixTransport::UnsaturatedFlowElement::addFields ( const std::string &  values,
const std::string &  fluxes,
const int  order 
)

add fields

Examples
UnsaturatedFlow.hpp.

Definition at line 1128 of file UnsaturatedFlow.hpp.

1129  {
1131  // Fields
1134  CHKERR mField.add_field(values + "_t", L2, AINSWORTH_LEGENDRE_BASE, 1);
1135  // CHKERR mField.add_field(fluxes+"_residual",L2,AINSWORTH_LEGENDRE_BASE,1);
1136 
1137  // meshset consisting all entities in mesh
1138  EntityHandle root_set = mField.get_moab().get_root_set();
1139  // add entities to field
1140 
1142  if (it->getName().compare(0, 4, "SOIL") != 0)
1143  continue;
1144  CHKERR mField.add_ents_to_field_by_type(dMatMap[it->getMeshsetId()]->tEts,
1145  MBTET, fluxes);
1146  CHKERR mField.add_ents_to_field_by_type(dMatMap[it->getMeshsetId()]->tEts,
1147  MBTET, values);
1148  CHKERR mField.add_ents_to_field_by_type(dMatMap[it->getMeshsetId()]->tEts,
1149  MBTET, values + "_t");
1150  // CHKERR mField.add_ents_to_field_by_type(
1151  // dMatMap[it->getMeshsetId()]->tEts,MBTET,fluxes+"_residual"
1152  // );
1153  }
1154 
1155  CHKERR mField.set_field_order(root_set, MBTET, fluxes, order + 1);
1156  CHKERR mField.set_field_order(root_set, MBTRI, fluxes, order + 1);
1157  CHKERR mField.set_field_order(root_set, MBTET, values, order);
1158  CHKERR mField.set_field_order(root_set, MBTET, values + "_t", order);
1159  // CHKERR mField.set_field_order(root_set,MBTET,fluxes+"_residual",order);
1161  }
field with continuous normal traction
Definition: definitions.h:179
virtual moab::Interface & get_moab()=0
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
virtual MoFEMErrorCode add_field(const std::string &name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
MaterialsDoubleMap dMatMap
materials database
constexpr int order
virtual MoFEMErrorCode add_ents_to_field_by_type(const Range &ents, const EntityType type, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:152
#define _IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet having a particular BC meshset in a moFEM field.
#define CHKERR
Inline error check.
Definition: definitions.h:604
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415
field with C-1 continuity
Definition: definitions.h:180

◆ addFiniteElements()

MoFEMErrorCode MixTransport::UnsaturatedFlowElement::addFiniteElements ( const std::string &  fluxes_name,
const std::string &  values_name 
)

add finite elements

Examples
UnsaturatedFlow.hpp.

Definition at line 1164 of file UnsaturatedFlow.hpp.

1165  {
1167 
1168  // Define element "MIX". Note that this element will work with fluxes_name
1169  // and values_name. This reflect bilinear form for the problem
1178  values_name + "_t");
1179  // CHKERR
1180  // mField.modify_finite_element_add_field_data("MIX",fluxes_name+"_residual");
1181 
1183  if (it->getName().compare(0, 4, "SOIL") != 0)
1184  continue;
1186  dMatMap[it->getMeshsetId()]->tEts, MBTET, "MIX");
1187  }
1188 
1189  // Define element to integrate natural boundary conditions, i.e. set values.
1190  CHKERR mField.add_finite_element("MIX_BCVALUE", MF_ZERO);
1192  fluxes_name);
1194  fluxes_name);
1196  fluxes_name);
1198  values_name);
1199 
1201  if (it->getName().compare(0, 4, "HEAD") != 0)
1202  continue;
1204  bcValueMap[it->getMeshsetId()]->eNts, MBTRI, "MIX_BCVALUE");
1205  }
1206 
1207  // Define element to apply essential boundary conditions.
1208  CHKERR mField.add_finite_element("MIX_BCFLUX", MF_ZERO);
1210  fluxes_name);
1212  fluxes_name);
1214  fluxes_name);
1216  values_name);
1217 
1219  if (it->getName().compare(0, 4, "FLUX") != 0)
1220  continue;
1222  bcFluxMap[it->getMeshsetId()]->eNts, MBTRI, "MIX_BCFLUX");
1223  }
1224 
1226  }
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string &name_row)=0
set field col which finite element use
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string &name, const bool recursive=true)=0
add entities to finite element
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
MaterialsDoubleMap dMatMap
materials database
BcMap bcValueMap
Store boundary condition for head capillary pressure.
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
#define _IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet having a particular BC meshset in a moFEM field.
#define CHKERR
Inline error check.
Definition: definitions.h:604
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string &name_filed)=0
set finite element field data
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string &name_row)=0
set field row which finite element use
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415

◆ buildProblem()

MoFEMErrorCode MixTransport::UnsaturatedFlowElement::buildProblem ( Range  zero_flux_ents,
BitRefLevel  ref_level = BitRefLevel().set(0) 
)

Build problem.

Parameters
ref_levelmesh refinement on which mesh problem you like to built.
Returns
error code
Examples
UnsaturatedFlow.hpp.

Definition at line 1233 of file UnsaturatedFlow.hpp.

1234  {
1236 
1237  // Build fields
1239  // Build finite elements
1241  CHKERR mField.build_finite_elements("MIX_BCFLUX");
1242  CHKERR mField.build_finite_elements("MIX_BCVALUE");
1243  // Build adjacencies of degrees of freedom and elements
1244  CHKERR mField.build_adjacencies(ref_level);
1245 
1246  // create DM instance
1247  CHKERR DMCreate(PETSC_COMM_WORLD, &dM);
1248  // setting that DM is type of DMMOFEM, i.e. MOFEM implementation manages DM
1249  CHKERR DMSetType(dM, "DMMOFEM");
1250  // mesh is portioned, each process keeps only part of problem
1251  CHKERR DMMoFEMSetIsPartitioned(dM, PETSC_TRUE);
1252  // creates problem in DM
1253  CHKERR DMMoFEMCreateMoFEM(dM, &mField, "MIX", ref_level);
1254  // discretised problem creates square matrix (that makes some optimizations)
1255  CHKERR DMMoFEMSetIsPartitioned(dM, PETSC_TRUE);
1256  // set DM options from command line
1257  CHKERR DMSetFromOptions(dM);
1258  // add finite elements
1259  CHKERR DMMoFEMAddElement(dM, "MIX");
1260  CHKERR DMMoFEMAddElement(dM, "MIX_BCFLUX");
1261  CHKERR DMMoFEMAddElement(dM, "MIX_BCVALUE");
1262  // constructor data structures
1263  CHKERR DMSetUp(dM);
1264 
1265  // remove zero flux dofs
1266  CHKERR mField.getInterface<ProblemsManager>()->removeDofsOnEntities(
1267  "MIX", "FLUXES", zero_flux_ents);
1268 
1269  PetscSection section;
1270  CHKERR mField.getInterface<ISManager>()->sectionCreate("MIX", &section);
1271  CHKERR DMSetDefaultSection(dM, section);
1272  CHKERR DMSetDefaultGlobalSection(dM, section);
1273  // CHKERR PetscSectionView(section,PETSC_VIEWER_STDOUT_WORLD);
1274  CHKERR PetscSectionDestroy(&section);
1275 
1277  }
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
PetscErrorCode DMMoFEMAddElement(DM dm, const char fe_name[])
add element to dm
Definition: DMMMoFEM.cpp:425
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
MoFEMErrorCode getInterface(const MOFEMuuid &uuid, IFACE *&iface) const
Get interface by uuid and return reference to pointer of interface.
PetscErrorCode DMMoFEMCreateMoFEM(DM dm, MoFEM::Interface *m_field_ptr, const char problem_name[], const MoFEM::BitRefLevel bit_level, const MoFEM::BitRefLevel bit_mask=MoFEM::BitRefLevel().set())
Must be called by user to set MoFEM data structures.
Definition: DMMMoFEM.cpp:105
#define CHKERR
Inline error check.
Definition: definitions.h:604
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
DM dM
Discrete manager for unsaturated flow problem.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415
PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned)
Definition: DMMMoFEM.cpp:982
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elementsBuild finite element data structures. Have to be run before problem and adjacenc...

◆ calculateEssentialBc()

MoFEMErrorCode MixTransport::UnsaturatedFlowElement::calculateEssentialBc ( )

Calculate boundary conditions for fluxes.

Returns
Error code
Examples
UnsaturatedFlow.hpp.

Definition at line 1602 of file UnsaturatedFlow.hpp.

1602  {
1604  // clear vectors
1605  CHKERR VecZeroEntries(D0);
1606  CHKERR VecGhostUpdateBegin(D0, INSERT_VALUES, SCATTER_FORWARD);
1607  CHKERR VecGhostUpdateEnd(D0, INSERT_VALUES, SCATTER_FORWARD);
1608  // clear essential bc indices, it could have dofs from other mesh refinement
1609  bcIndices.clear();
1610  // set operator to calculate essential boundary conditions
1611  CHKERR DMoFEMLoopFiniteElements(dM, "MIX_BCFLUX", feFaceBc);
1612  CHKERR VecGhostUpdateBegin(D0, INSERT_VALUES, SCATTER_REVERSE);
1613  CHKERR VecGhostUpdateEnd(D0, INSERT_VALUES, SCATTER_REVERSE);
1614  CHKERR VecAssemblyBegin(D0);
1615  CHKERR VecAssemblyEnd(D0);
1616  double norm2D0;
1617  CHKERR VecNorm(D0, NORM_2, &norm2D0);
1618  // CHKERR VecView(D0,PETSC_VIEWER_STDOUT_WORLD);
1619  PetscPrintf(PETSC_COMM_WORLD, "norm2D0 = %6.4e\n", norm2D0);
1621  }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
boost::shared_ptr< ForcesAndSourcesCore > feFaceBc
Elemnet to calculate essential bc.
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method)
Executes FEMethod for finite elements in DM.
Definition: DMMMoFEM.cpp:507
#define CHKERR
Inline error check.
Definition: definitions.h:604
DM dM
Discrete manager for unsaturated flow problem.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415

◆ calculateInitialPc()

MoFEMErrorCode MixTransport::UnsaturatedFlowElement::calculateInitialPc ( )

Calculate inital pressure head distribution.

Returns
Error code
Examples
UnsaturatedFlow.hpp.

Definition at line 1627 of file UnsaturatedFlow.hpp.

1627  {
1629  // clear vectors
1630  CHKERR VecZeroEntries(D1);
1631  CHKERR VecGhostUpdateBegin(D1, INSERT_VALUES, SCATTER_FORWARD);
1632  CHKERR VecGhostUpdateEnd(D1, INSERT_VALUES, SCATTER_FORWARD);
1633  // Calculate initial pressure head on each element
1635  // Assemble vector
1636  CHKERR VecGhostUpdateBegin(D1, INSERT_VALUES, SCATTER_REVERSE);
1637  CHKERR VecGhostUpdateEnd(D1, INSERT_VALUES, SCATTER_REVERSE);
1638  CHKERR VecAssemblyBegin(D1);
1639  CHKERR VecAssemblyEnd(D1);
1640  // Calculate norm
1641  double norm2D1;
1642  CHKERR VecNorm(D1, NORM_2, &norm2D1);
1643  // CHKERR VecView(D0,PETSC_VIEWER_STDOUT_WORLD);
1644  PetscPrintf(PETSC_COMM_WORLD, "norm2D1 = %6.4e\n", norm2D1);
1646  }
Vec D1
Vector with inital head capillary pressure.
boost::shared_ptr< ForcesAndSourcesCore > feVolInitialPc
Calculate inital boundary conditions.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method)
Executes FEMethod for finite elements in DM.
Definition: DMMMoFEM.cpp:507
#define CHKERR
Inline error check.
Definition: definitions.h:604
DM dM
Discrete manager for unsaturated flow problem.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415

◆ createMatrices()

MoFEMErrorCode MixTransport::UnsaturatedFlowElement::createMatrices ( )

Create vectors and matrices.

Returns
Error code
Examples
UnsaturatedFlow.hpp.

Definition at line 1568 of file UnsaturatedFlow.hpp.

1568  {
1570  CHKERR DMCreateMatrix(dM, &Aij);
1571  CHKERR DMCreateGlobalVector(dM, &D0);
1572  CHKERR VecDuplicate(D0, &D1);
1573  CHKERR VecDuplicate(D0, &D);
1574  CHKERR VecDuplicate(D0, &F);
1575  int ghosts[] = {0};
1576  int nb_locals = mField.get_comm_rank() == 0 ? 1 : 0;
1577  int nb_ghosts = mField.get_comm_rank() > 0 ? 1 : 0;
1578  CHKERR VecCreateGhost(PETSC_COMM_WORLD, nb_locals, 1, nb_ghosts, ghosts,
1579  &ghostFlux);
1581  }
Vec D1
Vector with inital head capillary pressure.
Vec ghostFlux
Ghost Vector of integrated fluxes.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
virtual int get_comm_rank() const =0
#define CHKERR
Inline error check.
Definition: definitions.h:604
DM dM
Discrete manager for unsaturated flow problem.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415

◆ destroyMatrices()

MoFEMErrorCode MixTransport::UnsaturatedFlowElement::destroyMatrices ( )

Delete matrices and vector when no longer needed.

Returns
error code
Examples
UnsaturatedFlow.hpp.

Definition at line 1587 of file UnsaturatedFlow.hpp.

1587  {
1589  CHKERR MatDestroy(&Aij);
1590  CHKERR VecDestroy(&D);
1591  CHKERR VecDestroy(&D0);
1592  CHKERR VecDestroy(&D1);
1593  CHKERR VecDestroy(&F);
1594  CHKERR VecDestroy(&ghostFlux);
1596  }
Vec D1
Vector with inital head capillary pressure.
Vec ghostFlux
Ghost Vector of integrated fluxes.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
#define CHKERR
Inline error check.
Definition: definitions.h:604
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415

◆ getBcOnFluxes()

MoFEMErrorCode MixTransport::UnsaturatedFlowElement::getBcOnFluxes ( const EntityHandle  ent,
const double  x,
const double  y,
const double  z,
double flux 
)
virtual

essential (Neumann) boundary condition (set fluxes)

Parameters
enthandle to finite element entity
xcoord
ycoord
zcoord
fluxreference to flux which is set by function
Returns
[description]

Reimplemented from MixTransport::MixTransportElement.

Examples
UnsaturatedFlow.hpp.

Definition at line 223 of file UnsaturatedFlow.hpp.

224  {
226  int block_id = -1;
227  if (lastEvalBcFluxEnt == ent) {
228  block_id = lastEvalBcBlockFluxId;
229  } else {
230  for (BcMap::iterator it = bcFluxMap.begin(); it != bcFluxMap.end();
231  it++) {
232  if (it->second->eNts.find(ent) != it->second->eNts.end()) {
233  block_id = it->first;
234  }
235  }
236  lastEvalBcFluxEnt = ent;
237  lastEvalBcBlockFluxId = block_id;
238  }
239  if (block_id >= 0) {
240  if (bcFluxMap.at(block_id)->hookFun) {
241  flux = bcFluxMap.at(block_id)->hookFun(x, y, z);
242  } else {
243  flux = bcFluxMap.at(block_id)->fixValue;
244  }
245  } else {
246  flux = 0;
247  }
249  }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415

◆ getBcOnValues()

MoFEMErrorCode MixTransport::UnsaturatedFlowElement::getBcOnValues ( const EntityHandle  ent,
const int  gg,
const double  x,
const double  y,
const double  z,
double value 
)
virtual

Get value on boundary.

Parameters
ententity handle
ggnumber of integration point
xx-coordinate
yy-coordinate
zz-coordinate
valuereturned value
Returns
error code

Reimplemented from MixTransport::MixTransportElement.

Examples
UnsaturatedFlow.hpp.

Definition at line 181 of file UnsaturatedFlow.hpp.

183  {
185  int block_id = -1;
186  if (lastEvalBcValEnt == ent) {
187  block_id = lastEvalBcBlockValId;
188  } else {
189  for (BcMap::iterator it = bcValueMap.begin(); it != bcValueMap.end();
190  it++) {
191  if (it->second->eNts.find(ent) != it->second->eNts.end()) {
192  block_id = it->first;
193  }
194  }
195  lastEvalBcValEnt = ent;
196  lastEvalBcBlockValId = block_id;
197  }
198  if (block_id >= 0) {
199  if (bcValueMap.at(block_id)->hookFun) {
200  value = bcValueMap.at(block_id)->hookFun(x, y, z);
201  } else {
202  value = bcValueMap.at(block_id)->fixValue;
203  }
204  } else {
205  value = 0;
206  }
208  }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
BcMap bcValueMap
Store boundary condition for head capillary pressure.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415

◆ getMaterial()

virtual MoFEMErrorCode MixTransport::UnsaturatedFlowElement::getMaterial ( const EntityHandle  ent,
int &  block_id 
) const
virtual

For given element handle get material block Id.

Parameters
entfinite element entity handle
block_idreference to returned block id
Returns
error code
Examples
UnsaturatedFlow.hpp.

Definition at line 140 of file UnsaturatedFlow.hpp.

141  {
143  for (MaterialsDoubleMap::const_iterator mit = dMatMap.begin();
144  mit != dMatMap.end(); mit++) {
145  if (mit->second->tEts.find(ent) != mit->second->tEts.end()) {
146  block_id = mit->first;
148  }
149  }
151  "Element not found, no material data");
153  }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:516
MaterialsDoubleMap dMatMap
materials database
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415
virtual MPI_Comm & get_comm() const =0

◆ setFiniteElements()

MoFEMErrorCode MixTransport::UnsaturatedFlowElement::setFiniteElements ( ForcesAndSourcesCore::RuleHookFun  vol_rule = VolRule(),
ForcesAndSourcesCore::RuleHookFun  face_rule = FaceRule() 
)

Create finite element instances.

Parameters
vol_ruleintegration rule for volume element
face_ruleintegration rule for boundary element
Returns
error code
Examples
UnsaturatedFlow.hpp.

Definition at line 1439 of file UnsaturatedFlow.hpp.

1440  {
1442 
1443  // create finite element instances
1444  feFaceBc = boost::shared_ptr<ForcesAndSourcesCore>(
1446  feFaceRhs = boost::shared_ptr<ForcesAndSourcesCore>(
1448  feVolInitialPc = boost::shared_ptr<ForcesAndSourcesCore>(
1450  feVolLhs = boost::shared_ptr<ForcesAndSourcesCore>(
1452  feVolRhs = boost::shared_ptr<ForcesAndSourcesCore>(
1454  // set integration rule to elements
1455  feFaceBc->getRuleHook = face_rule;
1456  feFaceRhs->getRuleHook = face_rule;
1457  feVolInitialPc->getRuleHook = vol_rule;
1458  feVolLhs->getRuleHook = vol_rule;
1459  feVolRhs->getRuleHook = vol_rule;
1460  // set function hook for finite element postprocessing stage
1461  feVolRhs->preProcessHook = preProcessVol(*this, feVolRhs);
1462  feVolLhs->preProcessHook = preProcessVol(*this, feVolLhs);
1463  feVolRhs->postProcessHook = postProcessVol(*this, feVolRhs);
1464  feVolLhs->postProcessHook = postProcessVol(*this, feVolLhs);
1465 
1466  // create method for setting history for fluxes on boundary
1467  scaleMethodFlux = boost::shared_ptr<MethodForForceScaling>(
1468  new TimeForceScale("-flux_history", false));
1469 
1470  // create method for setting history for presure heads on boundary
1471  scaleMethodValue = boost::shared_ptr<MethodForForceScaling>(
1472  new TimeForceScale("-head_history", false));
1473 
1474  // Set operator to calculate essential boundary conditions
1475  feFaceBc->getOpPtrVector().push_back(
1476  new OpEvaluateBcOnFluxes(*this, "FLUXES"));
1477 
1478  // Set operator to calculate initial capillary pressure
1479  feVolInitialPc->getOpPtrVector().push_back(
1480  new OpEvaluateInitiallHead(*this, "VALUES"));
1481 
1482  // set residual face from Neumann terms, i.e. applied pressure
1483  feFaceRhs->getOpPtrVector().push_back(
1484  new OpRhsBcOnValues(*this, "FLUXES", scaleMethodValue));
1485  // set residual finite element operators
1486  headRateAtGaussPts = boost::make_shared<VectorDouble>();
1487  // resAtGaussPts = boost::make_shared<VectorDouble>();
1488  feVolRhs->getOpPtrVector().push_back(new OpCalculateScalarFieldValues(
1489  string("VALUES") + "_t", headRateAtGaussPts, MBTET));
1490  feVolRhs->getOpPtrVector().push_back(
1491  new OpValuesAtGaussPts(*this, "VALUES"));
1492  feVolRhs->getOpPtrVector().push_back(
1493  new OpFluxDivergenceAtGaussPts(*this, "FLUXES"));
1494  feVolRhs->getOpPtrVector().push_back(new OpResidualFlux(*this, "FLUXES"));
1495  feVolRhs->getOpPtrVector().push_back(new OpResidualMass(*this, "VALUES"));
1496  feVolRhs->getOpPtrVector().back().opType =
1497  ForcesAndSourcesCore::UserDataOperator::OPROW;
1498  // set tangent matrix finite element operators
1499  feVolLhs->getOpPtrVector().push_back(new OpCalculateScalarFieldValues(
1500  string("VALUES") + "_t", headRateAtGaussPts, MBTET));
1501  // feVolLhs->getOpPtrVector().push_back(
1502  // new
1503  // OpCalculateScalarFieldValues(string("FLUXES")+"_residual",resAtGaussPts,MBTET)
1504  // );
1505  feVolLhs->getOpPtrVector().push_back(
1506  new OpValuesAtGaussPts(*this, "VALUES"));
1507  feVolLhs->getOpPtrVector().push_back(
1508  new OpFluxDivergenceAtGaussPts(*this, "FLUXES"));
1509  feVolLhs->getOpPtrVector().push_back(
1510  new OpTauDotSigma_HdivHdiv(*this, "FLUXES"));
1511  feVolLhs->getOpPtrVector().push_back(new OpVU_L2L2(*this, "VALUES"));
1512  feVolLhs->getOpPtrVector().push_back(
1513  new OpVDivSigma_L2Hdiv(*this, "VALUES", "FLUXES"));
1514  feVolLhs->getOpPtrVector().push_back(
1515  new OpDivTauU_HdivL2(*this, "FLUXES", "VALUES"));
1516 
1517  // Adding finite elements to DM, time solver will ask for it to assemble
1518  // tangent matrices and residuals
1519  boost::shared_ptr<FEMethod> null;
1520  CHKERR DMMoFEMTSSetIFunction(dM, "MIX_BCVALUE", feFaceRhs, null, null);
1521  CHKERR DMMoFEMTSSetIFunction(dM, "MIX", feVolRhs, null, null);
1522  CHKERR DMMoFEMTSSetIJacobian(dM, "MIX", feVolLhs, null, null);
1523 
1524  // setting up post-processing
1525  boost::shared_ptr<PostProcVolumeOnRefinedMesh> post_process(
1527  CHKERR post_process->generateReferenceElementMesh();
1528  CHKERR post_process->addFieldValuesPostProc("VALUES");
1529  CHKERR post_process->addFieldValuesPostProc("VALUES_t");
1530  // CHKERR post_process->addFieldValuesPostProc("FLUXES_residual");
1531  CHKERR post_process->addFieldValuesPostProc("FLUXES");
1532  post_process->getOpPtrVector().push_back(
1533  new OpValuesAtGaussPts(*this, "VALUES"));
1534  post_process->getOpPtrVector().push_back(
1535  new OpPostProcMaterial(*this, post_process->postProcMesh,
1536  post_process->mapGaussPts, "VALUES"));
1537 
1538  // Integrate fluxes on boundary
1539  boost::shared_ptr<ForcesAndSourcesCore> flux_integrate;
1540  flux_integrate = boost::shared_ptr<ForcesAndSourcesCore>(
1542  flux_integrate->getOpPtrVector().push_back(
1543  new OpIntegrateFluxes(*this, "FLUXES"));
1544  int frequency = 1;
1545  CHKERR PetscOptionsBegin(PETSC_COMM_WORLD, "", "Monitor post proc", "none");
1546  CHKERR PetscOptionsInt(
1547  "-how_often_output",
1548  "frequency how often results are dumped on hard disk", "", frequency,
1549  &frequency, NULL);
1550  ierr = PetscOptionsEnd();
1551  CHKERRG(ierr);
1552 
1553  tsMonitor = boost::shared_ptr<FEMethod>(
1554  new MonitorPostProc(*this, post_process, flux_integrate, frequency));
1555  TsCtx *ts_ctx;
1556  DMMoFEMGetTsCtx(dM, &ts_ctx);
1557  ts_ctx->get_postProcess_to_do_Monitor().push_back(tsMonitor);
1559  }
boost::shared_ptr< ForcesAndSourcesCore > feVolInitialPc
Calculate inital boundary conditions.
Force scale operator for reading two columns.
VolumeElementForcesAndSourcesCoreSwitch< 0 > VolumeElementForcesAndSourcesCore
Volume finite element default.
boost::shared_ptr< ForcesAndSourcesCore > feFaceRhs
Face element apply natural bc.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:552
boost::shared_ptr< ForcesAndSourcesCore > feFaceBc
Elemnet to calculate essential bc.
boost::shared_ptr< FEMethod > tsMonitor
PetscErrorCode DMMoFEMTSSetIJacobian(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS Jacobian evaluation function
Definition: DMMMoFEM.cpp:772
boost::shared_ptr< VectorDouble > headRateAtGaussPts
Vector keeps head rate.
boost::shared_ptr< ForcesAndSourcesCore > feVolRhs
Assemble residual vector.
Post processing.
boost::shared_ptr< MethodForForceScaling > scaleMethodFlux
Method scaling fluxes.
#define CHKERR
Inline error check.
Definition: definitions.h:604
boost::shared_ptr< MethodForForceScaling > scaleMethodValue
Method scaling values.
PetscErrorCode DMMoFEMGetTsCtx(DM dm, MoFEM::TsCtx **ts_ctx)
get MoFEM::TsCtx data structure
Definition: DMMMoFEM.cpp:1001
CHKERR DMMoFEMTSSetIFunction(dm, simple_interface->getDomainFEName(), null, fe_set_coords, null)
boost::shared_ptr< ForcesAndSourcesCore > feVolLhs
Assemble tangent matrix.
DM dM
Discrete manager for unsaturated flow problem.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415
FaceElementForcesAndSourcesCoreSwitch< 0 > FaceElementForcesAndSourcesCore
Face finite element default.

◆ solveProblem()

MoFEMErrorCode MixTransport::UnsaturatedFlowElement::solveProblem ( bool  set_initial_pc = true)

solve problem

Returns
error code
Examples
UnsaturatedFlow.hpp.

Definition at line 1652 of file UnsaturatedFlow.hpp.

1652  {
1654  if (set_initial_pc) {
1655  // Set initial head
1656  CHKERR DMoFEMMeshToLocalVector(dM, D1, INSERT_VALUES, SCATTER_REVERSE);
1657  }
1658 
1659  // Initiate vector from data on the mesh
1660  CHKERR DMoFEMMeshToLocalVector(dM, D, INSERT_VALUES, SCATTER_FORWARD);
1661 
1662  // Create time solver
1663  TS ts;
1664  CHKERR TSCreate(PETSC_COMM_WORLD, &ts);
1665  // Use backward Euler method
1666  CHKERR TSSetType(ts, TSBEULER);
1667  // Set final time
1668  double ftime = 1;
1669  CHKERR TSSetDuration(ts, PETSC_DEFAULT, ftime);
1670  // Setup solver from commabd line
1671  CHKERR TSSetFromOptions(ts);
1672  // Set DM to TS
1673  CHKERR TSSetDM(ts, dM);
1674 #if PETSC_VERSION_GE(3, 7, 0)
1675  CHKERR TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER);
1676 #endif
1677  // Set-up monitor
1678  TsCtx *ts_ctx;
1679  DMMoFEMGetTsCtx(dM, &ts_ctx);
1680  CHKERR TSMonitorSet(ts, TsMonitorSet, ts_ctx, PETSC_NULL);
1681 
1682  // This add SNES monitor, to show error by fields. It is dirty trick
1683  // to add monitor, so code is hidden from doxygen
1684  CHKERR TSSetSolution(ts, D);
1685  CHKERR TSSetUp(ts);
1686  SNES snes;
1687  CHKERR TSGetSNES(ts, &snes);
1688 
1689 #if PETSC_VERSION_GE(3, 7, 0)
1690  {
1691  PetscViewerAndFormat *vf;
1692  CHKERR PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,
1693  PETSC_VIEWER_DEFAULT, &vf);
1694  CHKERR SNESMonitorSet(
1695  snes,
1696  (MoFEMErrorCode(*)(SNES, PetscInt, PetscReal,
1697  void *))SNESMonitorFields,
1698  vf, (MoFEMErrorCode(*)(void **))PetscViewerAndFormatDestroy);
1699  }
1700 #else
1701  {
1702  CHKERR SNESMonitorSet(snes,
1703  (MoFEMErrorCode(*)(SNES, PetscInt, PetscReal,
1704  void *))SNESMonitorFields,
1705  0, 0);
1706  }
1707 #endif
1708 
1709  CHKERR TSSolve(ts, D);
1710 
1711  // Get statistic form TS and print it
1712  CHKERR TSGetTime(ts, &ftime);
1713  PetscInt steps, snesfails, rejects, nonlinits, linits;
1714  CHKERR TSGetTimeStepNumber(ts, &steps);
1715  CHKERR TSGetSNESFailures(ts, &snesfails);
1716  CHKERR TSGetStepRejections(ts, &rejects);
1717  CHKERR TSGetSNESIterations(ts, &nonlinits);
1718  CHKERR TSGetKSPIterations(ts, &linits);
1719  PetscPrintf(PETSC_COMM_WORLD,
1720  "steps %D (%D rejected, %D SNES fails), ftime %g, nonlinits "
1721  "%D, linits %D\n",
1722  steps, rejects, snesfails, ftime, nonlinits, linits);
1723 
1725  }
Vec D1
Vector with inital head capillary pressure.
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition: DMMMoFEM.cpp:445
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
PetscErrorCode TsMonitorSet(TS ts, PetscInt step, PetscReal t, Vec u, void *ctx)
Set monitor for TS solver.
Definition: TsCtx.cpp:209
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:67
#define CHKERR
Inline error check.
Definition: definitions.h:604
PetscErrorCode DMMoFEMGetTsCtx(DM dm, MoFEM::TsCtx **ts_ctx)
get MoFEM::TsCtx data structure
Definition: DMMMoFEM.cpp:1001
DM dM
Discrete manager for unsaturated flow problem.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415

Member Data Documentation

◆ bcFluxMap

BcMap MixTransport::UnsaturatedFlowElement::bcFluxMap
Examples
UnsaturatedFlow.hpp.

Definition at line 210 of file UnsaturatedFlow.hpp.

◆ bcValueMap

BcMap MixTransport::UnsaturatedFlowElement::bcValueMap

Store boundary condition for head capillary pressure.

Examples
UnsaturatedFlow.hpp.

Definition at line 166 of file UnsaturatedFlow.hpp.

◆ bcVecIds

std::vector<int> MixTransport::UnsaturatedFlowElement::bcVecIds
Examples
UnsaturatedFlow.hpp.

Definition at line 1316 of file UnsaturatedFlow.hpp.

◆ bcVecVals

VectorDouble MixTransport::UnsaturatedFlowElement::bcVecVals
Examples
UnsaturatedFlow.hpp.

Definition at line 1317 of file UnsaturatedFlow.hpp.

◆ D1

Vec MixTransport::UnsaturatedFlowElement::D1

Vector with inital head capillary pressure.

Examples
UnsaturatedFlow.hpp.

Definition at line 1561 of file UnsaturatedFlow.hpp.

◆ dM

DM MixTransport::UnsaturatedFlowElement::dM

Discrete manager for unsaturated flow problem.

Examples
UnsaturatedFlow.hpp.

Definition at line 117 of file UnsaturatedFlow.hpp.

◆ dMatMap

MaterialsDoubleMap MixTransport::UnsaturatedFlowElement::dMatMap

materials database

Examples
UnsaturatedFlow.hpp.

Definition at line 132 of file UnsaturatedFlow.hpp.

◆ feFaceBc

boost::shared_ptr<ForcesAndSourcesCore> MixTransport::UnsaturatedFlowElement::feFaceBc

Elemnet to calculate essential bc.

Examples
UnsaturatedFlow.hpp.

Definition at line 1280 of file UnsaturatedFlow.hpp.

◆ feFaceRhs

boost::shared_ptr<ForcesAndSourcesCore> MixTransport::UnsaturatedFlowElement::feFaceRhs

Face element apply natural bc.

Examples
UnsaturatedFlow.hpp.

Definition at line 1282 of file UnsaturatedFlow.hpp.

◆ feVolInitialPc

boost::shared_ptr<ForcesAndSourcesCore> MixTransport::UnsaturatedFlowElement::feVolInitialPc

Calculate inital boundary conditions.

Examples
UnsaturatedFlow.hpp.

Definition at line 1284 of file UnsaturatedFlow.hpp.

◆ feVolLhs

boost::shared_ptr<ForcesAndSourcesCore> MixTransport::UnsaturatedFlowElement::feVolLhs

Assemble tangent matrix.

Examples
UnsaturatedFlow.hpp.

Definition at line 1287 of file UnsaturatedFlow.hpp.

◆ feVolRhs

boost::shared_ptr<ForcesAndSourcesCore> MixTransport::UnsaturatedFlowElement::feVolRhs

Assemble residual vector.

Examples
UnsaturatedFlow.hpp.

Definition at line 1286 of file UnsaturatedFlow.hpp.

◆ ghostFlux

Vec MixTransport::UnsaturatedFlowElement::ghostFlux

Ghost Vector of integrated fluxes.

Examples
UnsaturatedFlow.hpp.

Definition at line 1562 of file UnsaturatedFlow.hpp.

◆ headRateAtGaussPts

boost::shared_ptr<VectorDouble> MixTransport::UnsaturatedFlowElement::headRateAtGaussPts

Vector keeps head rate.

Examples
UnsaturatedFlow.hpp.

Definition at line 1296 of file UnsaturatedFlow.hpp.

◆ lastEvalBcBlockFluxId

int MixTransport::UnsaturatedFlowElement::lastEvalBcBlockFluxId
Examples
UnsaturatedFlow.hpp.

Definition at line 212 of file UnsaturatedFlow.hpp.

◆ lastEvalBcBlockValId

int MixTransport::UnsaturatedFlowElement::lastEvalBcBlockValId
Examples
UnsaturatedFlow.hpp.

Definition at line 169 of file UnsaturatedFlow.hpp.

◆ lastEvalBcFluxEnt

EntityHandle MixTransport::UnsaturatedFlowElement::lastEvalBcFluxEnt
Examples
UnsaturatedFlow.hpp.

Definition at line 211 of file UnsaturatedFlow.hpp.

◆ lastEvalBcValEnt

EntityHandle MixTransport::UnsaturatedFlowElement::lastEvalBcValEnt
Examples
UnsaturatedFlow.hpp.

Definition at line 168 of file UnsaturatedFlow.hpp.

◆ scaleMethodFlux

boost::shared_ptr<MethodForForceScaling> MixTransport::UnsaturatedFlowElement::scaleMethodFlux

Method scaling fluxes.

Examples
UnsaturatedFlow.hpp.

Definition at line 1289 of file UnsaturatedFlow.hpp.

◆ scaleMethodValue

boost::shared_ptr<MethodForForceScaling> MixTransport::UnsaturatedFlowElement::scaleMethodValue

Method scaling values.

Examples
UnsaturatedFlow.hpp.

Definition at line 1291 of file UnsaturatedFlow.hpp.

◆ tsMonitor

boost::shared_ptr<FEMethod> MixTransport::UnsaturatedFlowElement::tsMonitor

Element used by TS monitor to postprocess results at time step

Examples
UnsaturatedFlow.hpp.

Definition at line 1292 of file UnsaturatedFlow.hpp.

◆ vecValsOnBc

VectorDouble MixTransport::UnsaturatedFlowElement::vecValsOnBc
Examples
UnsaturatedFlow.hpp.

Definition at line 1317 of file UnsaturatedFlow.hpp.


The documentation for this struct was generated from the following file: