v0.14.0
Public Member Functions | Private Attributes | List of all members
SeepageOps::OpDomainRhsHydrostaticStress< DIM > Struct Template Reference

#include <users_modules/tutorials/adv-5/src/SeepageOps.hpp>

Inheritance diagram for SeepageOps::OpDomainRhsHydrostaticStress< DIM >:
[legend]
Collaboration diagram for SeepageOps::OpDomainRhsHydrostaticStress< DIM >:
[legend]

Public Member Functions

 OpDomainRhsHydrostaticStress (std::string field_name1, boost::shared_ptr< VectorDouble > h_ptr, double specific_weight_water=9.81)
 
MoFEMErrorCode iNtegrate (DataForcesAndSourcesCore::EntData &data)
 
- Public Member Functions inherited from MoFEM::OpBaseImpl< A, EleOp >
 OpBaseImpl (const std::string row_field_name, const std::string col_field_name, const OpType type, boost::shared_ptr< Range > ents_ptr=nullptr)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
 Do calculations for the left hand side. More...
 
MoFEMErrorCode doWork (int row_side, EntityType row_type, EntData &row_data)
 Do calculations for the right hand side. More...
 

Private Attributes

double specificWeightWater
 
boost::shared_ptr< VectorDouble > hPtr
 

Additional Inherited Members

- Public Types inherited from MoFEM::OpBaseImpl< A, EleOp >
using OpType = typename EleOp::OpType
 
using EntData = EntitiesFieldData::EntData
 
using MatSetValuesHook = boost::function< MoFEMErrorCode(ForcesAndSourcesCore::UserDataOperator *op_ptr, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, MatrixDouble &m)>
 
- Public Attributes inherited from MoFEM::OpBaseImpl< A, EleOp >
TimeFun timeScalingFun
 assumes that time variable is set More...
 
FEFun feScalingFun
 assumes that time variable is set More...
 
boost::shared_ptr< RangeentsPtr
 Entities on which element is run. More...
 
- Static Public Attributes inherited from MoFEM::OpBaseImpl< A, EleOp >
static MatSetValuesHook matSetValuesHook
 
- Protected Member Functions inherited from MoFEM::OpBaseImpl< A, EleOp >
template<int DIM>
FTensor::Tensor1< FTensor::PackPtr< double *, DIM >, DIM > getNf ()
 
template<int DIM>
FTensor::Tensor2< FTensor::PackPtr< double *, DIM >, DIM, DIM > getLocMat (const int rr)
 
virtual MoFEMErrorCode iNtegrate (EntData &row_data, EntData &col_data)
 Integrate grad-grad operator. More...
 
virtual MoFEMErrorCode aSsemble (EntData &row_data, EntData &col_data, const bool trans)
 
virtual MoFEMErrorCode iNtegrate (EntData &data)
 Class dedicated to integrate operator. More...
 
virtual MoFEMErrorCode aSsemble (EntData &data)
 
virtual size_t getNbOfBaseFunctions (EntitiesFieldData::EntData &data)
 Get number of base functions. More...
 
- Protected Attributes inherited from MoFEM::OpBaseImpl< A, EleOp >
int nbRows
 number of dofs on rows More...
 
int nbCols
 number if dof on column More...
 
int nbIntegrationPts
 number of integration points More...
 
int nbRowBaseFunctions
 number or row base functions More...
 
int rowSide
 row side number More...
 
int colSide
 column side number More...
 
EntityType rowType
 row type More...
 
EntityType colType
 column type More...
 
bool assembleTranspose
 
bool onlyTranspose
 
MatrixDouble locMat
 local entity block matrix More...
 
MatrixDouble locMatTranspose
 local entity block matrix More...
 
VectorDouble locF
 local entity vector More...
 

Detailed Description

template<int DIM>
struct SeepageOps::OpDomainRhsHydrostaticStress< DIM >

Examples
seepage.cpp.

Definition at line 22 of file SeepageOps.hpp.

Constructor & Destructor Documentation

◆ OpDomainRhsHydrostaticStress()

template<int DIM>
SeepageOps::OpDomainRhsHydrostaticStress< DIM >::OpDomainRhsHydrostaticStress ( std::string  field_name1,
boost::shared_ptr< VectorDouble >  h_ptr,
double  specific_weight_water = 9.81 
)
inline
Examples
SeepageOps.hpp.

Definition at line 25 of file SeepageOps.hpp.

28  : AssemblyDomainEleOp(field_name1, field_name1, DomainEleOp::OPROW),
29  hPtr(h_ptr), specificWeightWater(specific_weight_water) {}

Member Function Documentation

◆ iNtegrate()

template<int DIM>
MoFEMErrorCode SeepageOps::OpDomainRhsHydrostaticStress< DIM >::iNtegrate ( DataForcesAndSourcesCore::EntData data)
inline
Examples
SeepageOps.hpp.

Definition at line 31 of file SeepageOps.hpp.

31  {
33 
34  const int nb_dofs = data.getIndices().size();
35 
36  if (nb_dofs) {
37  // locRhs.resize(nb_dofs, false);
38  // locRhs.clear();
39  auto &nf = AssemblyDomainEleOp::locF;
40  // get element area
41  const double area = getMeasure();
42 
43  // get number of integration points
44  const int nb_integration_points = getGaussPts().size2();
45  // get integration weights
46  auto t_w = getFTensor0IntegrationWeight();
47 
48  // get base function
49  auto t_base_diff = data.getFTensor1DiffN<DIM>();
50 
51  constexpr double g_acceleration = 9.81;
52 
54 
55  auto t_h = getFTensor0FromVec(*hPtr);
56  for (int gg = 0; gg != nb_integration_points; gg++) {
57  auto t_nf = getFTensor1FromPtr<DIM>(&nf[0]);
58 
59  const double a = t_w * area * specificWeightWater * t_h;
60 
61  for (int rr = 0; rr != nb_dofs / DIM; rr++) {
62  t_nf(i) -= t_base_diff(i) * a;
63 
64  // move to the next base function
65  ++t_base_diff; // moves the pointer to the next shape function
66  ++t_nf;
67  }
68 
69  // move to the weight of the next integration point
70  ++t_w;
71  ++t_h;
72  }
73 
74  }
75 
77  }

Member Data Documentation

◆ hPtr

template<int DIM>
boost::shared_ptr<VectorDouble> SeepageOps::OpDomainRhsHydrostaticStress< DIM >::hPtr
private
Examples
SeepageOps.hpp.

Definition at line 82 of file SeepageOps.hpp.

◆ specificWeightWater

template<int DIM>
double SeepageOps::OpDomainRhsHydrostaticStress< DIM >::specificWeightWater
private
Examples
SeepageOps.hpp.

Definition at line 81 of file SeepageOps.hpp.


The documentation for this struct was generated from the following file:
SeepageOps::OpDomainRhsHydrostaticStress::specificWeightWater
double specificWeightWater
Definition: SeepageOps.hpp:81
MoFEM::OpBaseImpl::locF
VectorDouble locF
local entity vector
Definition: FormsIntegrators.hpp:241
a
constexpr double a
Definition: approx_sphere.cpp:30
SeepageOps::OpDomainRhsHydrostaticStress::hPtr
boost::shared_ptr< VectorDouble > hPtr
Definition: SeepageOps.hpp:82
MoFEM::getFTensor0FromVec
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
FTensor::Index< 'i', DIM >
AssemblyDomainEleOp
FormsIntegrators< DomainEleOp >::Assembly< A >::OpBase AssemblyDomainEleOp
Definition: tensor_divergence_operator.cpp:59
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346