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

Public Attributes

DM dM
 Discrete manager for unsaturated flow problem.
 
MaterialsDoubleMap dMatMap
 materials database
 
BcMap bcValueMap
 Store boundary condition for head capillary pressure.
 
EntityHandle lastEvalBcValEnt
 
int lastEvalBcBlockValId
 
BcMap bcFluxMap
 
EntityHandle lastEvalBcFluxEnt
 
int lastEvalBcBlockFluxId
 
boost::shared_ptr< ForcesAndSourcesCorefeFaceBc
 Elemnet to calculate essential bc.
 
boost::shared_ptr< ForcesAndSourcesCorefeFaceRhs
 Face element apply natural bc.
 
boost::shared_ptr< ForcesAndSourcesCorefeVolInitialPc
 Calculate inital boundary conditions.
 
boost::shared_ptr< ForcesAndSourcesCorefeVolRhs
 Assemble residual vector.
 
boost::shared_ptr< ForcesAndSourcesCorefeVolLhs
 Assemble tangent matrix.
 
boost::shared_ptr< MethodForForceScalingscaleMethodFlux
 Method scaling fluxes.
 
boost::shared_ptr< MethodForForceScalingscaleMethodValue
 Method scaling values.
 
boost::shared_ptr< FEMethodtsMonitor
 
boost::shared_ptr< VectorDouble > headRateAtGaussPts
 Vector keeps head rate.
 
std::vector< int > bcVecIds
 
VectorDouble bcVecVals
 
VectorDouble vecValsOnBc
 
Vec D1
 Vector with inital head capillary pressure.
 
Vec ghostFlux
 Ghost Vector of integrated fluxes.
 
- Public Attributes inherited from MixTransport::MixTransportElement
MoFEM::InterfacemField
 
MyVolumeFE feVol
 Instance of volume element.
 
MyTriFE feTri
 Instance of surface element.
 
VectorDouble valuesAtGaussPts
 values at integration points on element
 
MatrixDouble valuesGradientAtGaussPts
 gradients at integration points on element
 
VectorDouble divergenceAtGaussPts
 divergence at integration points on element
 
MatrixDouble fluxesAtGaussPts
 fluxes at integration points on element
 
set< int > bcIndices
 
std::map< int, BlockDatasetOfBlocks
 maps block set id with appropriate BlockData
 
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
mofem/tutorials/cor-0to1_unsaturated_transport/unsaturated_transport.cpp.

Definition at line 102 of file UnsaturatedFlow.hpp.

Member Typedef Documentation

◆ BcMap

typedef map<int, boost::shared_ptr<BcData> > MixTransport::UnsaturatedFlowElement::BcMap

◆ MaterialsDoubleMap

Constructor & Destructor Documentation

◆ UnsaturatedFlowElement()

MixTransport::UnsaturatedFlowElement::UnsaturatedFlowElement ( MoFEM::Interface m_field)
inline

Definition at line 106 of file UnsaturatedFlow.hpp.

107 : MixTransportElement(m_field), dM(PETSC_NULLPTR), lastEvalBcValEnt(0),
MixTransportElement(MoFEM::Interface &m_field)
construction of this data structure
DM dM
Discrete manager for unsaturated flow problem.

◆ ~UnsaturatedFlowElement()

MixTransport::UnsaturatedFlowElement::~UnsaturatedFlowElement ( )
inline
Examples
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 111 of file UnsaturatedFlow.hpp.

111 {
112 if (dM != PETSC_NULLPTR) {
113 CHKERR DMDestroy(&dM);
114 CHKERRABORT(PETSC_COMM_WORLD, ierr);
115 }
116 }
static PetscErrorCode ierr
#define CHKERR
Inline error check.

Member Function Documentation

◆ addFields()

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

add fields

Examples
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 1136 of file UnsaturatedFlow.hpp.

1137 {
1139 // Fields
1142 CHKERR mField.add_field(values + "_t", L2, AINSWORTH_LEGENDRE_BASE, 1);
1143 // CHKERR mField.add_field(fluxes+"_residual",L2,AINSWORTH_LEGENDRE_BASE,1);
1144
1145 // meshset consisting all entities in mesh
1146 EntityHandle root_set = mField.get_moab().get_root_set();
1147 // add entities to field
1148
1150 if (it->getName().compare(0, 4, "SOIL") != 0)
1151 continue;
1152 CHKERR mField.add_ents_to_field_by_type(dMatMap[it->getMeshsetId()]->tEts,
1153 MBTET, fluxes);
1154 CHKERR mField.add_ents_to_field_by_type(dMatMap[it->getMeshsetId()]->tEts,
1155 MBTET, values);
1156 CHKERR mField.add_ents_to_field_by_type(dMatMap[it->getMeshsetId()]->tEts,
1157 MBTET, values + "_t");
1158 // CHKERR mField.add_ents_to_field_by_type(
1159 // dMatMap[it->getMeshsetId()]->tEts,MBTET,fluxes+"_residual"
1160 // );
1161 }
1162
1163 CHKERR mField.set_field_order(root_set, MBTET, fluxes, order + 1);
1164 CHKERR mField.set_field_order(root_set, MBTRI, fluxes, order + 1);
1165 CHKERR mField.set_field_order(root_set, MBTET, values, order);
1166 CHKERR mField.set_field_order(root_set, MBTET, values + "_t", order);
1167 // CHKERR mField.set_field_order(root_set,MBTET,fluxes+"_residual",order);
1169 }
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
@ DEMKOWICZ_JACOBI_BASE
Definition definitions.h:66
@ L2
field with C-1 continuity
Definition definitions.h:88
@ HDIV
field with continuous normal traction
Definition definitions.h:87
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ BLOCKSET
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
constexpr int order
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.
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.
#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.
MaterialsDoubleMap dMatMap
materials database
virtual moab::Interface & get_moab()=0
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.

◆ addFiniteElements()

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

add finite elements

Examples
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 1172 of file UnsaturatedFlow.hpp.

1173 {
1175
1176 // Define element "MIX". Note that this element will work with fluxes_name
1177 // and values_name. This reflect bilinear form for the problem
1186 values_name + "_t");
1187 // CHKERR
1188 // mField.modify_finite_element_add_field_data("MIX",fluxes_name+"_residual");
1189
1191 if (it->getName().compare(0, 4, "SOIL") != 0)
1192 continue;
1194 dMatMap[it->getMeshsetId()]->tEts, MBTET, "MIX");
1195 }
1196
1197 // Define element to integrate natural boundary conditions, i.e. set values.
1198 CHKERR mField.add_finite_element("MIX_BCVALUE", MF_ZERO);
1200 fluxes_name);
1202 fluxes_name);
1204 fluxes_name);
1206 values_name);
1207
1209 if (it->getName().compare(0, 4, "HEAD") != 0)
1210 continue;
1212 bcValueMap[it->getMeshsetId()]->eNts, MBTRI, "MIX_BCVALUE");
1213 }
1214
1215 // Define element to apply essential boundary conditions.
1216 CHKERR mField.add_finite_element("MIX_BCFLUX", MF_ZERO);
1218 fluxes_name);
1220 fluxes_name);
1222 fluxes_name);
1224 values_name);
1225
1227 if (it->getName().compare(0, 4, "FLUX") != 0)
1228 continue;
1230 bcFluxMap[it->getMeshsetId()]->eNts, MBTRI, "MIX_BCFLUX");
1231 }
1232
1234 }
@ MF_ZERO
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
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
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
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
BcMap bcValueMap
Store boundary condition for head capillary pressure.

◆ buildProblem()

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

Build problem.

Parameters
ref_levelmesh refinement on which mesh problem you like to built.
Returns
error code
Examples
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 1241 of file UnsaturatedFlow.hpp.

1242 {
1244
1245 // Build fields
1247 // Build finite elements
1249 CHKERR mField.build_finite_elements("MIX_BCFLUX");
1250 CHKERR mField.build_finite_elements("MIX_BCVALUE");
1251 // Build adjacencies of degrees of freedom and elements
1252 CHKERR mField.build_adjacencies(ref_level);
1253
1254 // create DM instance
1255 CHKERR DMCreate(PETSC_COMM_WORLD, &dM);
1256 // setting that DM is type of DMMOFEM, i.e. MOFEM implementation manages DM
1257 CHKERR DMSetType(dM, "DMMOFEM");
1258 // mesh is portioned, each process keeps only part of problem
1259 CHKERR DMMoFEMSetIsPartitioned(dM, PETSC_TRUE);
1260 // creates problem in DM
1261 CHKERR DMMoFEMCreateMoFEM(dM, &mField, "MIX", ref_level);
1262 // discretised problem creates square matrix (that makes some optimizations)
1263 CHKERR DMMoFEMSetIsPartitioned(dM, PETSC_TRUE);
1264 // set DM options from command line
1265 CHKERR DMSetFromOptions(dM);
1266 // add finite elements
1267 CHKERR DMMoFEMAddElement(dM, "MIX");
1268 CHKERR DMMoFEMAddElement(dM, "MIX_BCFLUX");
1269 CHKERR DMMoFEMAddElement(dM, "MIX_BCVALUE");
1270 // constructor data structures
1271 CHKERR DMSetUp(dM);
1272
1273 // remove zero flux dofs
1274 CHKERR mField.getInterface<ProblemsManager>()->removeDofsOnEntities(
1275 "MIX", "FLUXES", zero_flux_ents);
1276
1277 PetscSection section;
1278 CHKERR mField.getInterface<ISManager>()->sectionCreate("MIX", &section);
1279 CHKERR DMSetSection(dM, section);
1280 // CHKERR PetscSectionView(section,PETSC_VIEWER_STDOUT_WORLD);
1281 CHKERR PetscSectionDestroy(&section);
1282
1284 }
PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned)
Definition DMMoFEM.cpp:1113
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Definition DMMoFEM.cpp:488
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 DMMoFEM.cpp:114
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.

◆ calculateEssentialBc()

MoFEMErrorCode MixTransport::UnsaturatedFlowElement::calculateEssentialBc ( )
inline

Calculate boundary conditions for fluxes.

Returns
Error code

Definition at line 1648 of file UnsaturatedFlow.hpp.

1648 {
1650 // clear vectors
1651 CHKERR VecZeroEntries(D0);
1652 CHKERR VecGhostUpdateBegin(D0, INSERT_VALUES, SCATTER_FORWARD);
1653 CHKERR VecGhostUpdateEnd(D0, INSERT_VALUES, SCATTER_FORWARD);
1654 // clear essential bc indices, it could have dofs from other mesh refinement
1655 bcIndices.clear();
1656 // set operator to calculate essential boundary conditions
1658 CHKERR VecAssemblyBegin(D0);
1659 CHKERR VecAssemblyEnd(D0);
1660 CHKERR VecGhostUpdateBegin(D0, INSERT_VALUES, SCATTER_FORWARD);
1661 CHKERR VecGhostUpdateEnd(D0, INSERT_VALUES, SCATTER_FORWARD);
1662 double norm2D0;
1663 CHKERR VecNorm(D0, NORM_2, &norm2D0);
1664 // CHKERR VecView(D0,PETSC_VIEWER_STDOUT_WORLD);
1665 PetscPrintf(PETSC_COMM_WORLD, "norm2D0 = %6.4e\n", norm2D0);
1667 }
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition DMMoFEM.cpp:576
boost::shared_ptr< ForcesAndSourcesCore > feFaceBc
Elemnet to calculate essential bc.

◆ calculateInitialPc()

MoFEMErrorCode MixTransport::UnsaturatedFlowElement::calculateInitialPc ( )
inline

Calculate inital pressure head distribution.

Returns
Error code

Definition at line 1673 of file UnsaturatedFlow.hpp.

1673 {
1675 // clear vectors
1676 CHKERR VecZeroEntries(D1);
1677 CHKERR VecGhostUpdateBegin(D1, INSERT_VALUES, SCATTER_FORWARD);
1678 CHKERR VecGhostUpdateEnd(D1, INSERT_VALUES, SCATTER_FORWARD);
1679 // Calculate initial pressure head on each element
1681 // Assemble vector
1682 CHKERR VecAssemblyBegin(D1);
1683 CHKERR VecAssemblyEnd(D1);
1684 CHKERR VecGhostUpdateBegin(D1, INSERT_VALUES, SCATTER_FORWARD);
1685 CHKERR VecGhostUpdateEnd(D1, INSERT_VALUES, SCATTER_FORWARD);
1686 // Calculate norm
1687 double norm2D1;
1688 CHKERR VecNorm(D1, NORM_2, &norm2D1);
1689 // CHKERR VecView(D0,PETSC_VIEWER_STDOUT_WORLD);
1690 PetscPrintf(PETSC_COMM_WORLD, "norm2D1 = %6.4e\n", norm2D1);
1692 }
Vec D1
Vector with inital head capillary pressure.
boost::shared_ptr< ForcesAndSourcesCore > feVolInitialPc
Calculate inital boundary conditions.

◆ createMatrices()

MoFEMErrorCode MixTransport::UnsaturatedFlowElement::createMatrices ( )
inline

Create vectors and matrices.

Returns
Error code

Definition at line 1614 of file UnsaturatedFlow.hpp.

1614 {
1616 CHKERR DMCreateMatrix(dM, &Aij);
1617 CHKERR DMCreateGlobalVector(dM, &D0);
1618 CHKERR VecDuplicate(D0, &D1);
1619 CHKERR VecDuplicate(D0, &D);
1620 CHKERR VecDuplicate(D0, &F);
1621 int ghosts[] = {0};
1622 int nb_locals = mField.get_comm_rank() == 0 ? 1 : 0;
1623 int nb_ghosts = mField.get_comm_rank() > 0 ? 1 : 0;
1624 CHKERR VecCreateGhost(PETSC_COMM_WORLD, nb_locals, 1, nb_ghosts, ghosts,
1625 &ghostFlux);
1627 }
Vec ghostFlux
Ghost Vector of integrated fluxes.
virtual int get_comm_rank() const =0

◆ destroyMatrices()

MoFEMErrorCode MixTransport::UnsaturatedFlowElement::destroyMatrices ( )
inline

Delete matrices and vector when no longer needed.

Returns
error code

Definition at line 1633 of file UnsaturatedFlow.hpp.

1633 {
1635 CHKERR MatDestroy(&Aij);
1636 CHKERR VecDestroy(&D);
1637 CHKERR VecDestroy(&D0);
1638 CHKERR VecDestroy(&D1);
1639 CHKERR VecDestroy(&F);
1640 CHKERR VecDestroy(&ghostFlux);
1642 }

◆ getBcOnFluxes()

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

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
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 210 of file UnsaturatedFlow.hpp.

211 {
213 int block_id = -1;
214 if (lastEvalBcFluxEnt == ent) {
215 block_id = lastEvalBcBlockFluxId;
216 } else {
217 for (BcMap::iterator it = bcFluxMap.begin(); it != bcFluxMap.end();
218 it++) {
219 if (it->second->eNts.find(ent) != it->second->eNts.end()) {
220 block_id = it->first;
221 }
222 }
223 lastEvalBcFluxEnt = ent;
224 lastEvalBcBlockFluxId = block_id;
225 }
226 if (block_id >= 0) {
227 if (bcFluxMap.at(block_id)->hookFun) {
228 flux = bcFluxMap.at(block_id)->hookFun(x, y, z);
229 } else {
230 flux = bcFluxMap.at(block_id)->fixValue;
231 }
232 } else {
233 flux = 0;
234 }
236 }

◆ getBcOnValues()

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

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
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 168 of file UnsaturatedFlow.hpp.

170 {
172 int block_id = -1;
173 if (lastEvalBcValEnt == ent) {
174 block_id = lastEvalBcBlockValId;
175 } else {
176 for (BcMap::iterator it = bcValueMap.begin(); it != bcValueMap.end();
177 it++) {
178 if (it->second->eNts.find(ent) != it->second->eNts.end()) {
179 block_id = it->first;
180 }
181 }
182 lastEvalBcValEnt = ent;
183 lastEvalBcBlockValId = block_id;
184 }
185 if (block_id >= 0) {
186 if (bcValueMap.at(block_id)->hookFun) {
187 value = bcValueMap.at(block_id)->hookFun(x, y, z);
188 } else {
189 value = bcValueMap.at(block_id)->fixValue;
190 }
191 } else {
192 value = 0;
193 }
195 }

◆ getMaterial()

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

For given element handle get material block Id.

Parameters
entfinite element entity handle
block_idreference to returned block id
Returns
error code
Examples
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 127 of file UnsaturatedFlow.hpp.

128 {
130 for (MaterialsDoubleMap::const_iterator mit = dMatMap.begin();
131 mit != dMatMap.end(); mit++) {
132 if (mit->second->tEts.find(ent) != mit->second->tEts.end()) {
133 block_id = mit->first;
135 }
136 }
138 "Element not found, no material data");
140 }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
virtual MPI_Comm & get_comm() const =0

◆ setFiniteElements()

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

Create finite element instances.

Parameters
vol_ruleintegration rule for volume element
face_ruleintegration rule for boundary element
Returns
error code
Examples
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 1445 of file UnsaturatedFlow.hpp.

1446 {
1448
1449 // create finite element instances
1450 feFaceBc = boost::shared_ptr<ForcesAndSourcesCore>(
1452 feFaceRhs = boost::shared_ptr<ForcesAndSourcesCore>(
1454 feVolInitialPc = boost::shared_ptr<ForcesAndSourcesCore>(
1455 new VolumeElementForcesAndSourcesCore(mField));
1456 feVolLhs = boost::shared_ptr<ForcesAndSourcesCore>(
1457 new VolumeElementForcesAndSourcesCore(mField));
1458 feVolRhs = boost::shared_ptr<ForcesAndSourcesCore>(
1459 new VolumeElementForcesAndSourcesCore(mField));
1460 // set integration rule to elements
1461 feFaceBc->getRuleHook = face_rule;
1462 feFaceRhs->getRuleHook = face_rule;
1463 feVolInitialPc->getRuleHook = vol_rule;
1464 feVolLhs->getRuleHook = vol_rule;
1465 feVolRhs->getRuleHook = vol_rule;
1466 // set function hook for finite element postprocessing stage
1467 feVolRhs->preProcessHook = preProcessVol(*this, feVolRhs);
1468 feVolLhs->preProcessHook = preProcessVol(*this, feVolLhs);
1469 feVolRhs->postProcessHook = postProcessVol(*this, feVolRhs);
1470 feVolLhs->postProcessHook = postProcessVol(*this, feVolLhs);
1471
1472 // Set Piola transform
1473 CHKERR AddHOOps<2, 3, 3>::add(feFaceBc->getOpPtrVector(), {HDIV});
1474 CHKERR AddHOOps<2, 3, 3>::add(feFaceRhs->getOpPtrVector(), {HDIV});
1475
1476 // create method for setting history for fluxes on boundary
1477 scaleMethodFlux = boost::shared_ptr<MethodForForceScaling>(
1478 new TimeForceScale("-flux_history", false));
1479
1480 // create method for setting history for presure heads on boundary
1481 scaleMethodValue = boost::shared_ptr<MethodForForceScaling>(
1482 new TimeForceScale("-head_history", false));
1483
1484 // Set operator to calculate essential boundary conditions
1485 feFaceBc->getOpPtrVector().push_back(
1486 new OpEvaluateBcOnFluxes(*this, "FLUXES"));
1487
1488 // Set operator to calculate initial capillary pressure
1489 CHKERR AddHOOps<3, 3, 3>::add(feVolInitialPc->getOpPtrVector(), {HDIV, L2});
1490 feVolInitialPc->getOpPtrVector().push_back(
1491 new OpEvaluateInitiallHead(*this, "VALUES"));
1492
1493 // set residual face from Neumann terms, i.e. applied pressure
1494 feFaceRhs->getOpPtrVector().push_back(
1495 new OpRhsBcOnValues(*this, "FLUXES", scaleMethodValue));
1496 // set residual finite element operators
1497 headRateAtGaussPts = boost::make_shared<VectorDouble>();
1498
1499
1500 // resAtGaussPts = boost::make_shared<VectorDouble>();
1501 CHKERR AddHOOps<3, 3, 3>::add(feVolRhs->getOpPtrVector(), {HDIV, L2});
1502 feVolRhs->getOpPtrVector().push_back(new OpCalculateScalarFieldValues(
1503 string("VALUES") + "_t", headRateAtGaussPts, MBTET));
1504 feVolRhs->getOpPtrVector().push_back(
1505 new OpValuesAtGaussPts(*this, "VALUES"));
1506 feVolRhs->getOpPtrVector().push_back(
1507 new OpFluxDivergenceAtGaussPts(*this, "FLUXES"));
1508 feVolRhs->getOpPtrVector().push_back(new OpResidualFlux(*this, "FLUXES"));
1509 feVolRhs->getOpPtrVector().push_back(new OpResidualMass(*this, "VALUES"));
1510 feVolRhs->getOpPtrVector().back().opType =
1511 ForcesAndSourcesCore::UserDataOperator::OPROW;
1512 // set tangent matrix finite element operators
1513 CHKERR AddHOOps<3, 3, 3>::add(feVolLhs->getOpPtrVector(), {HDIV, L2});
1514 feVolLhs->getOpPtrVector().push_back(new OpCalculateScalarFieldValues(
1515 string("VALUES") + "_t", headRateAtGaussPts, MBTET));
1516 // feVolLhs->getOpPtrVector().push_back(
1517 // new
1518 // OpCalculateScalarFieldValues(string("FLUXES")+"_residual",resAtGaussPts,MBTET)
1519 // );
1520 feVolLhs->getOpPtrVector().push_back(
1521 new OpValuesAtGaussPts(*this, "VALUES"));
1522 feVolLhs->getOpPtrVector().push_back(
1523 new OpFluxDivergenceAtGaussPts(*this, "FLUXES"));
1524 feVolLhs->getOpPtrVector().push_back(
1525 new OpTauDotSigma_HdivHdiv(*this, "FLUXES"));
1526 feVolLhs->getOpPtrVector().push_back(new OpVU_L2L2(*this, "VALUES"));
1527 feVolLhs->getOpPtrVector().push_back(
1528 new OpVDivSigma_L2Hdiv(*this, "VALUES", "FLUXES"));
1529 feVolLhs->getOpPtrVector().push_back(
1530 new OpDivTauU_HdivL2(*this, "FLUXES", "VALUES"));
1531
1532 // Adding finite elements to DM, time solver will ask for it to assemble
1533 // tangent matrices and residuals
1534 boost::shared_ptr<FEMethod> null;
1535 CHKERR DMMoFEMTSSetIFunction(dM, "MIX_BCVALUE", feFaceRhs, null, null);
1536 CHKERR DMMoFEMTSSetIFunction(dM, "MIX", feVolRhs, null, null);
1537 CHKERR DMMoFEMTSSetIJacobian(dM, "MIX", feVolLhs, null, null);
1538
1539 // setting up post-processing
1540
1541
1542 auto get_post_process_ele = [&]() {
1543 auto post_process = boost::make_shared<PostProcEle>(mField);
1544
1545 CHKERR AddHOOps<3, 3, 3>::add(post_process->getOpPtrVector(), {HDIV, L2});
1546
1547 auto values_ptr = boost::make_shared<VectorDouble>();
1548 auto grad_ptr = boost::make_shared<MatrixDouble>();
1549 auto flux_ptr = boost::make_shared<MatrixDouble>();
1550
1551 post_process->getOpPtrVector().push_back(
1552 new OpCalculateScalarFieldValues("VALUES", values_ptr));
1553 post_process->getOpPtrVector().push_back(
1554 new OpCalculateScalarFieldGradient<3>("VALUES", grad_ptr));
1555 post_process->getOpPtrVector().push_back(
1556 new OpCalculateHVecVectorField<3>("FLUXES", flux_ptr));
1557
1558 using OpPPMap = OpPostProcMapInMoab<3, 3>;
1559
1560 post_process->getOpPtrVector().push_back(
1561
1562 new OpPPMap(post_process->getPostProcMesh(),
1563 post_process->getMapGaussPts(),
1564
1565 {{"VALUES", values_ptr}},
1566
1567 {{"GRAD", grad_ptr}, {"FLUXES", flux_ptr}},
1568
1569 {}, {}
1570
1571 ));
1572
1573 return post_process;
1574 };
1575
1576 auto post_process = get_post_process_ele();
1577
1578 post_process->getOpPtrVector().push_back(
1579 new OpValuesAtGaussPts(*this, "VALUES"));
1580 post_process->getOpPtrVector().push_back(
1581 new OpPostProcMaterial(*this, post_process->getPostProcMesh(),
1582 post_process->getMapGaussPts(), "VALUES"));
1583
1584 // Integrate fluxes on boundary
1585 boost::shared_ptr<ForcesAndSourcesCore> flux_integrate;
1586 flux_integrate = boost::shared_ptr<ForcesAndSourcesCore>(
1588 CHKERR AddHOOps<2, 3, 3>::add(flux_integrate->getOpPtrVector(), {HDIV});
1589 flux_integrate->getOpPtrVector().push_back(
1590 new OpIntegrateFluxes(*this, "FLUXES"));
1591 int frequency = 1;
1592 PetscOptionsBegin(PETSC_COMM_WORLD, "", "Monitor post proc", "none");
1593 CHKERR PetscOptionsInt(
1594 "-how_often_output",
1595 "frequency how often results are dumped on hard disk", "", frequency,
1596 &frequency, NULL);
1597 PetscOptionsEnd();
1598
1599 tsMonitor = boost::shared_ptr<FEMethod>(
1600 new MonitorPostProc(*this, post_process, flux_integrate, frequency));
1601 TsCtx *ts_ctx;
1605 }
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 DMMoFEM.cpp:790
PetscErrorCode DMMoFEMGetTsCtx(DM dm, MoFEM::TsCtx **ts_ctx)
get MoFEM::TsCtx data structure
Definition DMMoFEM.cpp:1132
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 DMMoFEM.cpp:843
MoFEM::TsCtx * ts_ctx
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
boost::shared_ptr< ForcesAndSourcesCore > feVolLhs
Assemble tangent matrix.
boost::shared_ptr< VectorDouble > headRateAtGaussPts
Vector keeps head rate.
boost::shared_ptr< MethodForForceScaling > scaleMethodFlux
Method scaling fluxes.
boost::shared_ptr< ForcesAndSourcesCore > feFaceRhs
Face element apply natural bc.
boost::shared_ptr< ForcesAndSourcesCore > feVolRhs
Assemble residual vector.
boost::shared_ptr< FEMethod > tsMonitor
boost::shared_ptr< MethodForForceScaling > scaleMethodValue
Method scaling values.
BasicMethodsSequence & getPostProcessMonitor()
Get the postProcess to do Monitor object.
Definition TsCtx.hpp:148
Force scale operator for reading two columns.

◆ solveProblem()

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

solve problem

Returns
error code

Definition at line 1698 of file UnsaturatedFlow.hpp.

1698 {
1700 if (set_initial_pc) {
1701 // Set initial head
1702 CHKERR DMoFEMMeshToLocalVector(dM, D1, INSERT_VALUES, SCATTER_REVERSE);
1703 }
1704
1705 // Initiate vector from data on the mesh
1706 CHKERR DMoFEMMeshToLocalVector(dM, D, INSERT_VALUES, SCATTER_FORWARD);
1707
1708 // Create time solver
1709 TS ts;
1710 CHKERR TSCreate(PETSC_COMM_WORLD, &ts);
1711 // Use backward Euler method
1712 CHKERR TSSetType(ts, TSBEULER);
1713 // Set final time
1714 double ftime = 1;
1715 CHKERR TSSetMaxTime(ts, ftime);
1716 // Setup solver from commabd line
1717 CHKERR TSSetFromOptions(ts);
1718 // Set DM to TS
1719 CHKERR TSSetDM(ts, dM);
1720#if PETSC_VERSION_GE(3, 7, 0)
1721 CHKERR TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER);
1722#endif
1723 // Set-up monitor
1724 TsCtx *ts_ctx;
1726 CHKERR TSMonitorSet(ts, TsMonitorSet, ts_ctx, PETSC_NULLPTR);
1727
1728 // This add SNES monitor, to show error by fields. It is dirty trick
1729 // to add monitor, so code is hidden from doxygen
1730 CHKERR TSSetSolution(ts, D);
1731 CHKERR TSSetUp(ts);
1732 SNES snes;
1733 CHKERR TSGetSNES(ts, &snes);
1734
1735#if PETSC_VERSION_GE(3, 7, 0)
1736 {
1737 PetscViewerAndFormat *vf;
1738 CHKERR PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,
1739 PETSC_VIEWER_DEFAULT, &vf);
1740 CHKERR SNESMonitorSet(
1741 snes,
1742 (MoFEMErrorCode(*)(SNES, PetscInt, PetscReal,
1743 void *))SNESMonitorFields,
1744 vf, (MoFEMErrorCode(*)(void **))PetscViewerAndFormatDestroy);
1745 }
1746#else
1747 {
1748 CHKERR SNESMonitorSet(snes,
1749 (MoFEMErrorCode(*)(SNES, PetscInt, PetscReal,
1750 void *))SNESMonitorFields,
1751 0, 0);
1752 }
1753#endif
1754
1755 CHKERR TSSolve(ts, D);
1756
1757 // Get statistic form TS and print it
1758 CHKERR TSGetTime(ts, &ftime);
1759 PetscInt steps, snesfails, rejects, nonlinits, linits;
1760 CHKERR TSGetStepNumber(ts, &steps);
1761 CHKERR TSGetSNESFailures(ts, &snesfails);
1762 CHKERR TSGetStepRejections(ts, &rejects);
1763 CHKERR TSGetSNESIterations(ts, &nonlinits);
1764 CHKERR TSGetKSPIterations(ts, &linits);
1765 PetscPrintf(PETSC_COMM_WORLD,
1766 "steps %d (%d rejected, %d SNES fails), ftime %g, nonlinits "
1767 "%d, linits %d\n",
1768 steps, rejects, snesfails, ftime, nonlinits, linits);
1769
1771 }
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition DMMoFEM.cpp:514
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.

Member Data Documentation

◆ bcFluxMap

BcMap MixTransport::UnsaturatedFlowElement::bcFluxMap

◆ bcValueMap

BcMap MixTransport::UnsaturatedFlowElement::bcValueMap

Store boundary condition for head capillary pressure.

Examples
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 153 of file UnsaturatedFlow.hpp.

◆ bcVecIds

std::vector<int> MixTransport::UnsaturatedFlowElement::bcVecIds

◆ bcVecVals

VectorDouble MixTransport::UnsaturatedFlowElement::bcVecVals

◆ D1

Vec MixTransport::UnsaturatedFlowElement::D1

Vector with inital head capillary pressure.

Examples
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 1607 of file UnsaturatedFlow.hpp.

◆ dM

DM MixTransport::UnsaturatedFlowElement::dM

Discrete manager for unsaturated flow problem.

Examples
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 104 of file UnsaturatedFlow.hpp.

◆ dMatMap

MaterialsDoubleMap MixTransport::UnsaturatedFlowElement::dMatMap

◆ feFaceBc

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

Elemnet to calculate essential bc.

Examples
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 1287 of file UnsaturatedFlow.hpp.

◆ feFaceRhs

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

Face element apply natural bc.

Examples
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 1289 of file UnsaturatedFlow.hpp.

◆ feVolInitialPc

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

Calculate inital boundary conditions.

Examples
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 1291 of file UnsaturatedFlow.hpp.

◆ feVolLhs

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

Assemble tangent matrix.

Examples
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 1294 of file UnsaturatedFlow.hpp.

◆ feVolRhs

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

Assemble residual vector.

Examples
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 1293 of file UnsaturatedFlow.hpp.

◆ ghostFlux

Vec MixTransport::UnsaturatedFlowElement::ghostFlux

Ghost Vector of integrated fluxes.

Examples
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 1608 of file UnsaturatedFlow.hpp.

◆ headRateAtGaussPts

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

Vector keeps head rate.

Examples
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 1303 of file UnsaturatedFlow.hpp.

◆ lastEvalBcBlockFluxId

int MixTransport::UnsaturatedFlowElement::lastEvalBcBlockFluxId

◆ lastEvalBcBlockValId

int MixTransport::UnsaturatedFlowElement::lastEvalBcBlockValId

◆ lastEvalBcFluxEnt

EntityHandle MixTransport::UnsaturatedFlowElement::lastEvalBcFluxEnt

◆ lastEvalBcValEnt

EntityHandle MixTransport::UnsaturatedFlowElement::lastEvalBcValEnt

◆ scaleMethodFlux

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

Method scaling fluxes.

Examples
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 1296 of file UnsaturatedFlow.hpp.

◆ scaleMethodValue

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

Method scaling values.

Examples
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 1298 of file UnsaturatedFlow.hpp.

◆ tsMonitor

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

Element used by TS monitor to postprocess results at time step

Examples
mofem/tutorials/cor-0to1_unsaturated_transport/src/UnsaturatedFlow.hpp.

Definition at line 1299 of file UnsaturatedFlow.hpp.

◆ vecValsOnBc

VectorDouble MixTransport::UnsaturatedFlowElement::vecValsOnBc

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