v0.14.0
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 SeepageOps.hpp
17  */
18 
19 namespace SeepageOps {
20 
21 template <int DIM>
23  : public AssemblyDomainEleOp { // changed opfaceele to AssemblyDomainEleOp
24 public:
25  OpDomainRhsHydrostaticStress(std::string field_name1,
26  boost::shared_ptr<VectorDouble> h_ptr,
27  double specific_weight_water = 9.81)
28  : AssemblyDomainEleOp(field_name1, field_name1, DomainEleOp::OPROW),
29  hPtr(h_ptr), specificWeightWater(specific_weight_water) {}
30 
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  }
78 
79 private:
80  // VectorDouble locRhs;
82  boost::shared_ptr<VectorDouble> hPtr;
83 };
84 
85 } // namespace SeepageOps
SeepageOps::OpDomainRhsHydrostaticStress::OpDomainRhsHydrostaticStress
OpDomainRhsHydrostaticStress(std::string field_name1, boost::shared_ptr< VectorDouble > h_ptr, double specific_weight_water=9.81)
Definition: SeepageOps.hpp:25
SeepageOps::OpDomainRhsHydrostaticStress::iNtegrate
MoFEMErrorCode iNtegrate(DataForcesAndSourcesCore::EntData &data)
Definition: SeepageOps.hpp:31
SeepageOps
Definition: SeepageOps.hpp:19
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
SeepageOps::OpDomainRhsHydrostaticStress::specificWeightWater
double specificWeightWater
Definition: SeepageOps.hpp:81
MoFEM::OpBaseImpl::locF
VectorDouble locF
local entity vector
Definition: FormsIntegrators.hpp:251
MoFEM::OpBaseImpl
Definition: FormsIntegrators.hpp:178
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
EntData
EntitiesFieldData::EntData EntData
Definition: child_and_parent.cpp:37
FTensor::Index< 'i', DIM >
DomainEleOp
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
SeepageOps::OpDomainRhsHydrostaticStress
Definition: SeepageOps.hpp:22