v0.15.4
Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
SeepageOps::OpDomainRhsHydrostaticStress< DIM > Struct Template Reference

#include "tutorials/adv-5_poroelasticity/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 > p_ptr, std::function< double(const double, const double, const double)> multiplier_func)
 
MoFEMErrorCode iNtegrate (DataForcesAndSourcesCore::EntData &data)
 

Private Attributes

boost::shared_ptr< VectorDouble > pPtr
 
std::function< double(const double, const double, const double)> multiplier_func
 

Detailed Description

template<int DIM>
struct SeepageOps::OpDomainRhsHydrostaticStress< DIM >
Examples
mofem/tutorials/adv-5_poroelasticity/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 >  p_ptr,
std::function< double(const double, const double, const double)>  multiplier_func 
)
inline

Definition at line 25 of file SeepageOps.hpp.

28 : AssemblyDomainEleOp(field_name1, field_name1, DomainEleOp::OPROW),
boost::shared_ptr< VectorDouble > pPtr
std::function< double(const double, const double, const double)> multiplier_func
FormsIntegrators< DomainEleOp >::Assembly< A >::OpBase AssemblyDomainEleOp

Member Function Documentation

◆ iNtegrate()

template<int DIM>
MoFEMErrorCode SeepageOps::OpDomainRhsHydrostaticStress< DIM >::iNtegrate ( DataForcesAndSourcesCore::EntData &  data)
inline
Examples
mofem/tutorials/adv-5_poroelasticity/src/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
38 auto &nf = AssemblyDomainEleOp::locF; //access to local vector
39 // get element measure
40 const double measure = getMeasure();
41
42 // get number of integration points
43 const int nb_integration_points = getGaussPts().size2();
44 // get integration weights
45 auto t_w = getFTensor0IntegrationWeight();
46
47 // get base function
48 auto t_base_diff = data.getFTensor1DiffN<DIM>();
49
50 //constexpr double g_acceleration = 9.81;
51
52 FTensor::Index<'i', DIM> i;
53
54 auto t_p = getFTensor0FromVec(*pPtr);
55 for (int gg = 0; gg != nb_integration_points; gg++) {
56
57 auto t_nf = getFTensor1FromPtr<DIM>(&nf[0]);
58
59 const double biot = multiplier_func(0.0, 0.0, 0.0);
60 const double a = t_w * measure * t_p * biot;
61
62 for (int rr = 0; rr != nb_dofs / DIM; rr++) { //loop over dofs (rr = alpha)
63 t_nf(i) -= t_base_diff(i) * a; //loop over dimensions
64
65 // move to the next base function
66 ++t_base_diff; // moves the pointer to the next shape function
67 ++t_nf;
68 }
69
70 // move to the weight of the next integration point
71 ++t_w;
72 ++t_p;
73 }
74
75 }
76
78 }
constexpr double a
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
FTensor::Index< 'i', SPACE_DIM > i

Member Data Documentation

◆ multiplier_func

template<int DIM>
std::function<double(const double, const double, const double)> SeepageOps::OpDomainRhsHydrostaticStress< DIM >::multiplier_func
private

◆ pPtr

template<int DIM>
boost::shared_ptr<VectorDouble> SeepageOps::OpDomainRhsHydrostaticStress< DIM >::pPtr
private

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