v0.9.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/basic_finite_elements/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 (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< intbcVecIds
 
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< intbcIndices
 
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 114 of file UnsaturatedFlow.hpp.

Member Typedef Documentation

◆ BcMap

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

Definition at line 164 of file UnsaturatedFlow.hpp.

◆ MaterialsDoubleMap

Examples
UnsaturatedFlow.hpp.

Definition at line 130 of file UnsaturatedFlow.hpp.

Constructor & Destructor Documentation

◆ UnsaturatedFlowElement()

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

Definition at line 118 of file UnsaturatedFlow.hpp.

119  : MixTransportElement(m_field), dM(PETSC_NULL), lastEvalBcValEnt(0),
121  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 123 of file UnsaturatedFlow.hpp.

123  {
124  if (dM != PETSC_NULL) {
125  CHKERR DMDestroy(&dM);
126  CHKERRABORT(PETSC_COMM_WORLD, ierr);
127  }
128  }
#define CHKERR
Inline error check.
Definition: definitions.h:596
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 1123 of file UnsaturatedFlow.hpp.

1124  {
1126  // Fields
1129  CHKERR mField.add_field(values + "_t", L2, AINSWORTH_LEGENDRE_BASE, 1);
1130  // CHKERR mField.add_field(fluxes+"_residual",L2,AINSWORTH_LEGENDRE_BASE,1);
1131 
1132  // meshset consisting all entities in mesh
1133  EntityHandle root_set = mField.get_moab().get_root_set();
1134  // add entities to field
1135 
1137  if (it->getName().compare(0, 4, "SOIL") != 0)
1138  continue;
1139  CHKERR mField.add_ents_to_field_by_type(dMatMap[it->getMeshsetId()]->tEts,
1140  MBTET, fluxes);
1141  CHKERR mField.add_ents_to_field_by_type(dMatMap[it->getMeshsetId()]->tEts,
1142  MBTET, values);
1143  CHKERR mField.add_ents_to_field_by_type(dMatMap[it->getMeshsetId()]->tEts,
1144  MBTET, values + "_t");
1145  // CHKERR mField.add_ents_to_field_by_type(
1146  // dMatMap[it->getMeshsetId()]->tEts,MBTET,fluxes+"_residual"
1147  // );
1148  }
1149 
1150  CHKERR mField.set_field_order(root_set, MBTET, fluxes, order + 1);
1151  CHKERR mField.set_field_order(root_set, MBTRI, fluxes, order + 1);
1152  CHKERR mField.set_field_order(root_set, MBTET, values, order);
1153  CHKERR mField.set_field_order(root_set, MBTET, values + "_t", order);
1154  // CHKERR mField.set_field_order(root_set,MBTET,fluxes+"_residual",order);
1156  }
field with continuous normal traction
Definition: definitions.h:173
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:477
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
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:146
#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:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
field with C-1 continuity
Definition: definitions.h:174

◆ addFiniteElements()

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

add finite elements

Examples
UnsaturatedFlow.hpp.

Definition at line 1159 of file UnsaturatedFlow.hpp.

1160  {
1162 
1163  // Define element "MIX". Note that this element will work with fluxes_name
1164  // and values_name. This reflect bilinear form for the problem
1173  values_name + "_t");
1174  // CHKERR
1175  // mField.modify_finite_element_add_field_data("MIX",fluxes_name+"_residual");
1176 
1178  if (it->getName().compare(0, 4, "SOIL") != 0)
1179  continue;
1181  dMatMap[it->getMeshsetId()]->tEts, MBTET, "MIX");
1182  }
1183 
1184  // Define element to integrate natural boundary conditions, i.e. set values.
1185  CHKERR mField.add_finite_element("MIX_BCVALUE", MF_ZERO);
1187  fluxes_name);
1189  fluxes_name);
1191  fluxes_name);
1193  values_name);
1194 
1196  if (it->getName().compare(0, 4, "HEAD") != 0)
1197  continue;
1199  bcValueMap[it->getMeshsetId()]->eNts, MBTRI, "MIX_BCVALUE");
1200  }
1201 
1202  // Define element to apply essential boundary conditions.
1203  CHKERR mField.add_finite_element("MIX_BCFLUX", MF_ZERO);
1205  fluxes_name);
1207  fluxes_name);
1209  fluxes_name);
1211  values_name);
1212 
1214  if (it->getName().compare(0, 4, "FLUX") != 0)
1215  continue;
1217  bcFluxMap[it->getMeshsetId()]->eNts, MBTRI, "MIX_BCFLUX");
1218  }
1219 
1221  }
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:477
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:596
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:407

◆ buildProblem()

MoFEMErrorCode MixTransport::UnsaturatedFlowElement::buildProblem ( 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 1228 of file UnsaturatedFlow.hpp.

1228  {
1230 
1231  // Build fields
1233  // Build finite elements
1235  CHKERR mField.build_finite_elements("MIX_BCFLUX");
1236  CHKERR mField.build_finite_elements("MIX_BCVALUE");
1237  // Build adjacencies of degrees of freedom and elements
1238  CHKERR mField.build_adjacencies(ref_level);
1239 
1240  // create DM instance
1241  CHKERR DMCreate(PETSC_COMM_WORLD, &dM);
1242  // setting that DM is type of DMMOFEM, i.e. MOFEM implementation manages DM
1243  CHKERR DMSetType(dM, "DMMOFEM");
1244  // mesh is portioned, each process keeps only part of problem
1245  CHKERR DMMoFEMSetIsPartitioned(dM, PETSC_TRUE);
1246  // creates problem in DM
1247  CHKERR DMMoFEMCreateMoFEM(dM, &mField, "MIX", ref_level);
1248  // discretised problem creates square matrix (that makes some optimizations)
1249  CHKERR DMMoFEMSetIsPartitioned(dM, PETSC_TRUE);
1250  // set DM options from command line
1251  CHKERR DMSetFromOptions(dM);
1252  // add finite elements
1253  CHKERR DMMoFEMAddElement(dM, "MIX");
1254  CHKERR DMMoFEMAddElement(dM, "MIX_BCFLUX");
1255  CHKERR DMMoFEMAddElement(dM, "MIX_BCVALUE");
1256  // constructor data structures
1257  CHKERR DMSetUp(dM);
1258 
1259  PetscSection section;
1260  CHKERR mField.getInterface<ISManager>()->sectionCreate("MIX", &section);
1261  CHKERR DMSetDefaultSection(dM, section);
1262  CHKERR DMSetDefaultGlobalSection(dM, section);
1263  // CHKERR PetscSectionView(section,PETSC_VIEWER_STDOUT_WORLD);
1264  CHKERR PetscSectionDestroy(&section);
1265 
1267  }
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:414
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
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:109
#define CHKERR
Inline error check.
Definition: definitions.h:596
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:407
PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned)
Definition: DMMMoFEM.cpp:971
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 1598 of file UnsaturatedFlow.hpp.

1598  {
1600  // clear vectors
1601  CHKERR VecZeroEntries(D0);
1602  CHKERR VecGhostUpdateBegin(D0, INSERT_VALUES, SCATTER_FORWARD);
1603  CHKERR VecGhostUpdateEnd(D0, INSERT_VALUES, SCATTER_FORWARD);
1604  // clear essential bc indices, it could have dofs from other mesh refinement
1605  bcIndices.clear();
1606  // set operator to calculate essential boundary conditions
1607  CHKERR DMoFEMLoopFiniteElements(dM, "MIX_BCFLUX", feFaceBc);
1608  CHKERR VecGhostUpdateBegin(D0, INSERT_VALUES, SCATTER_REVERSE);
1609  CHKERR VecGhostUpdateEnd(D0, INSERT_VALUES, SCATTER_REVERSE);
1610  CHKERR VecAssemblyBegin(D0);
1611  CHKERR VecAssemblyEnd(D0);
1612  double norm2D0;
1613  CHKERR VecNorm(D0, NORM_2, &norm2D0);
1614  // CHKERR VecView(D0,PETSC_VIEWER_STDOUT_WORLD);
1615  PetscPrintf(PETSC_COMM_WORLD, "norm2D0 = %6.4e\n", norm2D0);
1617  }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
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:496
#define CHKERR
Inline error check.
Definition: definitions.h:596
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:407

◆ calculateInitialPc()

MoFEMErrorCode MixTransport::UnsaturatedFlowElement::calculateInitialPc ( )

Calculate inital pressure head distribution.

Returns
Error code
Examples
UnsaturatedFlow.hpp.

Definition at line 1623 of file UnsaturatedFlow.hpp.

1623  {
1625  // clear vectors
1626  CHKERR VecZeroEntries(D1);
1627  CHKERR VecGhostUpdateBegin(D1, INSERT_VALUES, SCATTER_FORWARD);
1628  CHKERR VecGhostUpdateEnd(D1, INSERT_VALUES, SCATTER_FORWARD);
1629  // Calculate initial pressure head on each element
1631  // Assemble vector
1632  CHKERR VecGhostUpdateBegin(D1, INSERT_VALUES, SCATTER_REVERSE);
1633  CHKERR VecGhostUpdateEnd(D1, INSERT_VALUES, SCATTER_REVERSE);
1634  CHKERR VecAssemblyBegin(D1);
1635  CHKERR VecAssemblyEnd(D1);
1636  // Calculate norm
1637  double norm2D1;
1638  CHKERR VecNorm(D1, NORM_2, &norm2D1);
1639  // CHKERR VecView(D0,PETSC_VIEWER_STDOUT_WORLD);
1640  PetscPrintf(PETSC_COMM_WORLD, "norm2D1 = %6.4e\n", norm2D1);
1642  }
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:477
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method)
Executes FEMethod for finite elements in DM.
Definition: DMMMoFEM.cpp:496
#define CHKERR
Inline error check.
Definition: definitions.h:596
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:407

◆ createMatrices()

MoFEMErrorCode MixTransport::UnsaturatedFlowElement::createMatrices ( )

Create vectors and matrices.

Returns
Error code
Examples
UnsaturatedFlow.hpp.

Definition at line 1564 of file UnsaturatedFlow.hpp.

1564  {
1566  CHKERR DMCreateMatrix(dM, &Aij);
1567  CHKERR DMCreateGlobalVector(dM, &D0);
1568  CHKERR VecDuplicate(D0, &D1);
1569  CHKERR VecDuplicate(D0, &D);
1570  CHKERR VecDuplicate(D0, &F);
1571  int ghosts[] = {0};
1572  int nb_locals = mField.get_comm_rank() == 0 ? 1 : 0;
1573  int nb_ghosts = mField.get_comm_rank() > 0 ? 1 : 0;
1574  CHKERR VecCreateGhost(PETSC_COMM_WORLD, nb_locals, 1, nb_ghosts, ghosts,
1575  &ghostFlux);
1577  }
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:477
virtual int get_comm_rank() const =0
#define CHKERR
Inline error check.
Definition: definitions.h:596
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:407

◆ destroyMatrices()

MoFEMErrorCode MixTransport::UnsaturatedFlowElement::destroyMatrices ( )

Delete matrices and vector when no longer needed.

Returns
error code
Examples
UnsaturatedFlow.hpp.

Definition at line 1583 of file UnsaturatedFlow.hpp.

1583  {
1585  CHKERR MatDestroy(&Aij);
1586  CHKERR VecDestroy(&D);
1587  CHKERR VecDestroy(&D0);
1588  CHKERR VecDestroy(&D1);
1589  CHKERR VecDestroy(&F);
1590  CHKERR VecDestroy(&ghostFlux);
1592  }
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:477
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407

◆ 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 222 of file UnsaturatedFlow.hpp.

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

◆ 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 180 of file UnsaturatedFlow.hpp.

182  {
184  int block_id = -1;
185  if (lastEvalBcValEnt == ent) {
186  block_id = lastEvalBcBlockValId;
187  } else {
188  for (BcMap::iterator it = bcValueMap.begin(); it != bcValueMap.end();
189  it++) {
190  if (it->second->eNts.find(ent) != it->second->eNts.end()) {
191  block_id = it->first;
192  }
193  }
194  lastEvalBcValEnt = ent;
195  lastEvalBcBlockValId = block_id;
196  }
197  if (block_id >= 0) {
198  if (bcValueMap.at(block_id)->hookFun) {
199  value = bcValueMap.at(block_id)->hookFun(x, y, z);
200  } else {
201  value = bcValueMap.at(block_id)->fixValue;
202  }
203  } else {
204  value = 0;
205  }
207  }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
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:407

◆ 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 139 of file UnsaturatedFlow.hpp.

140  {
142  for (MaterialsDoubleMap::const_iterator mit = dMatMap.begin();
143  mit != dMatMap.end(); mit++) {
144  if (mit->second->tEts.find(ent) != mit->second->tEts.end()) {
145  block_id = mit->first;
147  }
148  }
150  "Element not found, no material data");
152  }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
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:407
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 1435 of file UnsaturatedFlow.hpp.

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

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

Member Data Documentation

◆ bcFluxMap

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

Definition at line 209 of file UnsaturatedFlow.hpp.

◆ bcValueMap

BcMap MixTransport::UnsaturatedFlowElement::bcValueMap

Store boundary condition for head capillary pressure.

Examples
UnsaturatedFlow.hpp.

Definition at line 165 of file UnsaturatedFlow.hpp.

◆ bcVecIds

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

Definition at line 1306 of file UnsaturatedFlow.hpp.

◆ bcVecVals

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

Definition at line 1307 of file UnsaturatedFlow.hpp.

◆ D1

Vec MixTransport::UnsaturatedFlowElement::D1

Vector with inital head capillary pressure.

Examples
UnsaturatedFlow.hpp.

Definition at line 1557 of file UnsaturatedFlow.hpp.

◆ dM

DM MixTransport::UnsaturatedFlowElement::dM

Discrete manager for unsaturated flow problem.

Examples
UnsaturatedFlow.hpp.

Definition at line 116 of file UnsaturatedFlow.hpp.

◆ dMatMap

MaterialsDoubleMap MixTransport::UnsaturatedFlowElement::dMatMap

materials database

Examples
UnsaturatedFlow.hpp.

Definition at line 131 of file UnsaturatedFlow.hpp.

◆ feFaceBc

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

Elemnet to calculate essential bc.

Examples
UnsaturatedFlow.hpp.

Definition at line 1270 of file UnsaturatedFlow.hpp.

◆ feFaceRhs

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

Face element apply natural bc.

Examples
UnsaturatedFlow.hpp.

Definition at line 1272 of file UnsaturatedFlow.hpp.

◆ feVolInitialPc

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

Calculate inital boundary conditions.

Examples
UnsaturatedFlow.hpp.

Definition at line 1274 of file UnsaturatedFlow.hpp.

◆ feVolLhs

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

Assemble tangent matrix.

Examples
UnsaturatedFlow.hpp.

Definition at line 1277 of file UnsaturatedFlow.hpp.

◆ feVolRhs

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

Assemble residual vector.

Examples
UnsaturatedFlow.hpp.

Definition at line 1276 of file UnsaturatedFlow.hpp.

◆ ghostFlux

Vec MixTransport::UnsaturatedFlowElement::ghostFlux

Ghost Vector of integrated fluxes.

Examples
UnsaturatedFlow.hpp.

Definition at line 1558 of file UnsaturatedFlow.hpp.

◆ headRateAtGaussPts

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

Vector keeps head rate.

Examples
UnsaturatedFlow.hpp.

Definition at line 1286 of file UnsaturatedFlow.hpp.

◆ lastEvalBcBlockFluxId

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

Definition at line 211 of file UnsaturatedFlow.hpp.

◆ lastEvalBcBlockValId

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

Definition at line 168 of file UnsaturatedFlow.hpp.

◆ lastEvalBcFluxEnt

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

Definition at line 210 of file UnsaturatedFlow.hpp.

◆ lastEvalBcValEnt

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

Definition at line 167 of file UnsaturatedFlow.hpp.

◆ scaleMethodFlux

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

Method scaling fluxes.

Examples
UnsaturatedFlow.hpp.

Definition at line 1279 of file UnsaturatedFlow.hpp.

◆ scaleMethodValue

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

Method scaling values.

Examples
UnsaturatedFlow.hpp.

Definition at line 1281 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 1282 of file UnsaturatedFlow.hpp.

◆ vecValsOnBc

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

Definition at line 1307 of file UnsaturatedFlow.hpp.


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