v0.15.4
Loading...
Searching...
No Matches
SeepageOps.hpp
Go to the documentation of this file.
1/* This file is part of MoFEM.
2 * MoFEM is free software: you can redistribute it and/or modify it under
3 * the terms of the GNU Lesser General Public License as published by the
4 * Free Software Foundation, either version 3 of the License, or (at your
5 * option) any later version.
6 *
7 * MoFEM is distributed in the hope that it will be useful, but WITHOUT
8 * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
9 * FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
10 * License for more details.
11 *
12 * You should have received a copy of the GNU Lesser General Public
13 * License along with MoFEM. If not, see <http://www.gnu.org/licenses/>. */
14
15/** \file SeepageOps.hpp
16 * \example mofem/tutorials/adv-5_poroelasticity/src/SeepageOps.hpp
17 */
18
19namespace SeepageOps {
20
21template <int DIM>
23 : public AssemblyDomainEleOp { // changed opfaceele to AssemblyDomainEleOp
24public:
25 OpDomainRhsHydrostaticStress(std::string field_name1,
26 boost::shared_ptr<VectorDouble> p_ptr,
27 std::function<double(const double, const double, const double)> multiplier_func)
28 : AssemblyDomainEleOp(field_name1, field_name1, DomainEleOp::OPROW),
30
31 MoFEMErrorCode iNtegrate(DataForcesAndSourcesCore::EntData &data) {
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 }
79
80private:
81 boost::shared_ptr<VectorDouble> pPtr;
82 std::function<double(const double, const double, const double)> multiplier_func;
83};
84
85} // namespace SeepageOps
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
MoFEMErrorCode iNtegrate(DataForcesAndSourcesCore::EntData &data)
boost::shared_ptr< VectorDouble > pPtr
std::function< double(const double, const double, const double)> multiplier_func
OpDomainRhsHydrostaticStress(std::string field_name1, boost::shared_ptr< VectorDouble > p_ptr, std::function< double(const double, const double, const double)> multiplier_func)