v0.15.0
Loading...
Searching...
No Matches
poisson_2d_lagrange_multiplier.hpp
Go to the documentation of this file.
1#ifndef __POISSON2DLAGRANGEMULTIPLIER_HPP__
2#define __POISSON2DLAGRANGEMULTIPLIER_HPP__
3
4#include <stdlib.h>
5#include <MoFEM.hpp>
6
7using namespace MoFEM;
8
11
14
16
18
20
21// const double body_source = 10;
22typedef boost::function<double(const double, const double, const double)>
24
25struct OpDomainLhsK : public OpFaceEle {
26public:
27 OpDomainLhsK(std::string row_field_name, std::string col_field_name,
28 boost::shared_ptr<std::vector<unsigned char>> boundary_marker = nullptr)
29 : OpFaceEle(row_field_name, col_field_name, OpFaceEle::OPROWCOL),
30 boundaryMarker(boundary_marker) {
31 sYmm = true;
32 }
33
34 MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
35 EntityType col_type, EntData &row_data,
36 EntData &col_data) {
38
39 const int nb_row_dofs = row_data.getIndices().size();
40 const int nb_col_dofs = col_data.getIndices().size();
41
42 if (nb_row_dofs && nb_col_dofs) {
43
44 locLhs.resize(nb_row_dofs, nb_col_dofs, false);
45 locLhs.clear();
46
47 // get element area
48 const double area = getMeasure();
49
50 // get number of integration points
51 const int nb_integration_points = getGaussPts().size2();
52 // get integration weights
54
55 // get derivatives of base functions on row
56 auto t_row_diff_base = row_data.getFTensor1DiffN<2>();
57
58 // START THE LOOP OVER INTEGRATION POINTS TO CALCULATE LOCAL MATRIX
59 for (int gg = 0; gg != nb_integration_points; gg++) {
60
61 const double a = t_w * area;
62
63 for (int rr = 0; rr != nb_row_dofs; ++rr) {
64 // get derivatives of base functions on column
65 auto t_col_diff_base = col_data.getFTensor1DiffN<2>(gg, 0);
66
67 for (int cc = 0; cc != nb_col_dofs; cc++) {
68
69 locLhs(rr, cc) += t_row_diff_base(i) * t_col_diff_base(i) * a;
70
71 // move to the derivatives of the next base functions on column
72 ++t_col_diff_base;
73 }
74
75 // move to the derivatives of the next base functions on row
76 ++t_row_diff_base;
77 }
78
79 // move to the weight of the next integration point
80 ++t_w;
81 }
82
83 // FILL VALUES OF LOCAL MATRIX ENTRIES TO THE GLOBAL MATRIX
84
85 CHKERR MatSetValues(getKSPB(), row_data, col_data, &locLhs(0, 0),
86 ADD_VALUES);
87
88 // Fill values of symmetric stiffness matrix in global system of equations
89 if (row_side != col_side || row_type != col_type) {
90 transLocLhs.resize(nb_col_dofs, nb_row_dofs, false);
91 noalias(transLocLhs) = trans(locLhs);
92 CHKERR MatSetValues(getKSPB(), col_data, row_data, &transLocLhs(0, 0),
93 ADD_VALUES);
94 }
95
96 }
97
99 }
100
101private:
102 boost::shared_ptr<std::vector<unsigned char>> boundaryMarker;
104};
105
106struct OpDomainRhsF : public OpFaceEle {
107public:
108 OpDomainRhsF(std::string field_name, ScalarFunc source_term_function,
109 boost::shared_ptr<std::vector<unsigned char>> boundary_marker = nullptr)
111 sourceTermFunc(source_term_function), boundaryMarker(boundary_marker) {}
112
113 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
115
116 const int nb_dofs = data.getIndices().size();
117
118 if (nb_dofs) {
119
120 locRhs.resize(nb_dofs, false);
121 locRhs.clear();
122
123 // get element area
124 const double area = getMeasure();
125
126 // get number of integration points
127 const int nb_integration_points = getGaussPts().size2();
128 // get integration weights
129 auto t_w = getFTensor0IntegrationWeight();
130 // get coordinates of the integration point
131 auto t_coords = getFTensor1CoordsAtGaussPts();
132
133 // get base functions
134 auto t_base = data.getFTensor0N();
135
136 // START THE LOOP OVER INTEGRATION POINTS TO CALCULATE LOCAL VECTOR
137 for (int gg = 0; gg != nb_integration_points; gg++) {
138
139 const double a = t_w * area;
140 double body_source =
141 sourceTermFunc(t_coords(0), t_coords(1), t_coords(2));
142
143 for (int rr = 0; rr != nb_dofs; rr++) {
144
145 locRhs[rr] += t_base * body_source * a;
146
147 // move to the next base function
148 ++t_base;
149 }
150
151 // move to the weight of the next integration point
152 ++t_w;
153 // move to the coordinates of the next integration point
154 ++t_coords;
155 }
156
157 // FILL VALUES OF THE GLOBAL VECTOR ENTRIES FROM THE LOCAL ONES
158
159 CHKERR VecSetValues(getKSPf(), data, &*locRhs.begin(), ADD_VALUES);
160 }
161
163 }
164
165private:
167 boost::shared_ptr<std::vector<unsigned char>> boundaryMarker;
169};
170
171struct OpBoundaryLhsC : public OpEdgeEle {
172public:
173 OpBoundaryLhsC(std::string row_field_name, std::string col_field_name)
174 : OpEdgeEle(row_field_name, col_field_name, OpEdgeEle::OPROWCOL) {
175 sYmm = false;
176 }
177
178 MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
179 EntityType col_type, EntData &row_data,
180 EntData &col_data) {
182
183 const int nb_row_dofs = row_data.getIndices().size();
184 const int nb_col_dofs = col_data.getIndices().size();
185
186 if (nb_row_dofs && nb_col_dofs) {
187
188 locBoundaryLhs.resize(nb_row_dofs, nb_col_dofs, false);
189 locBoundaryLhs.clear();
190
191 // get (boundary) element length
192 const double edge = getMeasure();
193
194 // get number of integration points
195 const int nb_integration_points = getGaussPts().size2();
196 // get integration weights
197 auto t_w = getFTensor0IntegrationWeight();
198
199 // get base functions on row
200 auto t_row_base = row_data.getFTensor0N();
201
202 // START THE LOOP OVER INTEGRATION POINTS TO CALCULATE LOCAL MATRIX
203 for (int gg = 0; gg != nb_integration_points; gg++) {
204 const double a = t_w * edge;
205
206 for (int rr = 0; rr != nb_row_dofs; ++rr) {
207 // get base functions on column
208 auto t_col_base = col_data.getFTensor0N(gg, 0);
209
210 for (int cc = 0; cc != nb_col_dofs; cc++) {
211
212 locBoundaryLhs(rr, cc) += t_row_base * t_col_base * a;
213
214 // move to the next base functions on column
215 ++t_col_base;
216 }
217 // move to the next base functions on row
218 ++t_row_base;
219 }
220
221 // move to the weight of the next integration point
222 ++t_w;
223 }
224
225 // FILL VALUES OF LOCAL MATRIX ENTRIES TO THE GLOBAL MATRIX
226 CHKERR MatSetValues(getKSPB(), row_data, col_data, &locBoundaryLhs(0, 0),
227 ADD_VALUES);
228
229 }
230
232 }
233
234private:
236};
237
238struct OpBoundaryRhsG : public OpEdgeEle {
239public:
240 OpBoundaryRhsG(std::string field_name, ScalarFunc boundary_function)
242 boundaryFunc(boundary_function) {}
243
244 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
246
247 const int nb_dofs = data.getIndices().size();
248
249 if (nb_dofs) {
250
251 locBoundaryRhs.resize(nb_dofs, false);
252 locBoundaryRhs.clear();
253
254 // get (boundary) element length
255 const double edge = getMeasure();
256
257 // get number of integration points
258 const int nb_integration_points = getGaussPts().size2();
259 // get integration weights
260 auto t_w = getFTensor0IntegrationWeight();
261 // get coordinates at integration point
262 auto t_coords = getFTensor1CoordsAtGaussPts();
263
264 // get base function
265 auto t_base = data.getFTensor0N();
266
267 // START THE LOOP OVER INTEGRATION POINTS TO CALCULATE LOCAL VECTOR
268 for (int gg = 0; gg != nb_integration_points; gg++) {
269
270 const double a = t_w * edge;
271 double boundary_term =
272 boundaryFunc(t_coords(0), t_coords(1), t_coords(2));
273
274 for (int rr = 0; rr != nb_dofs; rr++) {
275
276 locBoundaryRhs[rr] += t_base * boundary_term * a;
277
278 // move to the next base function
279 ++t_base;
280 }
281
282 // move to the weight of the next integration point
283 ++t_w;
284 // move to the coordinates of the next integration point
285 ++t_coords;
286 }
287
288 // FILL VALUES OF LOCAL VECTOR ENTRIES TO THE GLOBAL VECTOR
289 CHKERR VecSetValues(getKSPf(), data, &*locBoundaryRhs.begin(),
290 ADD_VALUES);
291 }
292
294 }
295
296private:
299};
300
301}; // namespace Poisson2DLagrangeMultiplierOperators
302
303#endif //__POISSON2DLAGRANGEMULTIPLIER_HPP__
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()
#define CHKERR
Inline error check.
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
boost::function< double(const double, const double, const double)> ScalarFunc
constexpr auto field_name
bool sYmm
If true assume that matrix is symmetric structure.
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
const VectorInt & getIndices() const
Get global indices of degrees of freedom on entity.
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
auto getFTensor0IntegrationWeight()
Get integration weights.
double getMeasure() const
get measure of element
@ OPROW
operator doWork function is executed on FE rows
@ OPROWCOL
operator doWork is executed on FE rows &columns
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
Operator for bi-linear form, usually to calculate values on left hand side.
OpBoundaryLhsC(std::string row_field_name, std::string col_field_name)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpBoundaryRhsG(std::string field_name, ScalarFunc boundary_function)
OpDomainLhsK(std::string row_field_name, std::string col_field_name, boost::shared_ptr< std::vector< unsigned char > > boundary_marker=nullptr)
boost::shared_ptr< std::vector< unsigned char > > boundaryMarker
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
Operator for bi-linear form, usually to calculate values on left hand side.
boost::shared_ptr< std::vector< unsigned char > > boundaryMarker
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpDomainRhsF(std::string field_name, ScalarFunc source_term_function, boost::shared_ptr< std::vector< unsigned char > > boundary_marker=nullptr)