v0.15.0
Loading...
Searching...
No Matches
SeepageOps::OpDomainRhsHydrostaticStress< DIM > Struct Template Reference

#include "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)
 

Private Attributes

double specificWeightWater
 
boost::shared_ptr< VectorDouble > hPtr
 

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) {}
boost::shared_ptr< VectorDouble > hPtr
FormsIntegrators< DomainEleOp >::Assembly< A >::OpBase AssemblyDomainEleOp

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 FTensor::Index<'i', DIM> i;
52
53 auto t_h = getFTensor0FromVec(*hPtr);
54 for (int gg = 0; gg != nb_integration_points; gg++) {
55 auto t_nf = getFTensor1FromPtr<DIM>(&nf[0]);
56
57 const double a = t_w * area * specificWeightWater * t_h;
58
59 for (int rr = 0; rr != nb_dofs / DIM; rr++) {
60 t_nf(i) -= t_base_diff(i) * a;
61
62 // move to the next base function
63 ++t_base_diff; // moves the pointer to the next shape function
64 ++t_nf;
65 }
66
67 // move to the weight of the next integration point
68 ++t_w;
69 ++t_h;
70 }
71 }
72
74 }
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

◆ hPtr

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

Definition at line 79 of file SeepageOps.hpp.

◆ specificWeightWater

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

Definition at line 78 of file SeepageOps.hpp.


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