v0.14.0
Loading...
Searching...
No Matches
poisson_2d_nonhomogeneous.hpp
Go to the documentation of this file.
1#ifndef __POISSON2DNONHOMOGENEOUS_HPP__
2#define __POISSON2DNONHOMOGENEOUS_HPP__
3#include <MoFEM.hpp>
4using namespace MoFEM;
5
8
11
13
15
17
18typedef boost::function<double(const double, const double, const double)>
20//! [OpDomainLhs]
21struct OpDomainLhs : public OpFaceEle {
22public:
23 OpDomainLhs(std::string row_field_name, std::string col_field_name)
24 : OpFaceEle(row_field_name, col_field_name, OpFaceEle::OPROWCOL) {
25 sYmm = true;
26 }
27
28 MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
29 EntityType col_type, EntData &row_data,
30 EntData &col_data) {
32
33 const int nb_row_dofs = row_data.getIndices().size();
34 const int nb_col_dofs = col_data.getIndices().size();
35
36 if (nb_row_dofs && nb_col_dofs) {
37
38 locLhs.resize(nb_row_dofs, nb_col_dofs, false);
39 locLhs.clear();
40
41 // get element area
42 const double area = getMeasure();
43
44 // get number of integration points
45 const int nb_integration_points = getGaussPts().size2();
46 // get integration weights
48
49 // get derivatives of base functions on row
50 auto t_row_diff_base = row_data.getFTensor1DiffN<2>();
51
52 // START THE LOOP OVER INTEGRATION POINTS TO CALCULATE LOCAL MATRIX
53 for (int gg = 0; gg != nb_integration_points; gg++) {
54
55 const double a = t_w * area;
56
57 for (int rr = 0; rr != nb_row_dofs; ++rr) {
58 // get derivatives of base functions on column
59 auto t_col_diff_base = col_data.getFTensor1DiffN<2>(gg, 0);
60
61 for (int cc = 0; cc != nb_col_dofs; cc++) {
62
63 locLhs(rr, cc) += t_row_diff_base(i) * t_col_diff_base(i) * a;
64
65 // move to the derivatives of the next base functions on column
66 ++t_col_diff_base;
67 }
68
69 // move to the derivatives of the next base functions on row
70 ++t_row_diff_base;
71 }
72
73 // move to the weight of the next integration point
74 ++t_w;
75 }
76
77 // FILL VALUES OF LOCAL MATRIX ENTRIES TO THE GLOBAL MATRIX
78
79 // Fill value to local stiffness matrix ignoring boundary DOFs
81 getKSPB(), row_data, col_data, &locLhs(0, 0), ADD_VALUES);
82
83 // Fill values of symmetric local stiffness matrix
84 if (row_side != col_side || row_type != col_type) {
85 transLocLhs.resize(nb_col_dofs, nb_row_dofs, false);
86 noalias(transLocLhs) = trans(locLhs);
88 getKSPB(), col_data, row_data, &transLocLhs(0, 0), ADD_VALUES);
89 }
90 }
91
93 }
94
95private:
96 boost::shared_ptr<std::vector<unsigned char>> boundaryMarker;
97 MatrixDouble locLhs, transLocLhs;
98};
99//! [OpDomainLhs]
100
101//! [OpDomainRhs]
102struct OpDomainRhs : public OpFaceEle {
103public:
104 OpDomainRhs(std::string field_name, ScalarFunc source_term_function)
106 sourceTermFunc(source_term_function) {}
107
108 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
110
111 const int nb_dofs = data.getIndices().size();
112
113 if (nb_dofs) {
114
115 locRhs.resize(nb_dofs, false);
116 locRhs.clear();
117
118 // get element area
119 const double area = getMeasure();
120
121 // get number of integration points
122 const int nb_integration_points = getGaussPts().size2();
123 // get integration weights
124 auto t_w = getFTensor0IntegrationWeight();
125 // get coordinates of the integration point
126 auto t_coords = getFTensor1CoordsAtGaussPts();
127
128 // get base functions
129 auto t_base = data.getFTensor0N();
130
131 // START THE LOOP OVER INTEGRATION POINTS TO CALCULATE LOCAL VECTOR
132 for (int gg = 0; gg != nb_integration_points; gg++) {
133
134 const double a = t_w * area;
135 double body_source =
136 sourceTermFunc(t_coords(0), t_coords(1), t_coords(2));
137
138 for (int rr = 0; rr != nb_dofs; rr++) {
139
140 locRhs[rr] += t_base * body_source * a;
141
142 // move to the next base function
143 ++t_base;
144 }
145
146 // move to the weight of the next integration point
147 ++t_w;
148 // move to the coordinates of the next integration point
149 ++t_coords;
150 }
151
152 // FILL VALUES OF THE GLOBAL VECTOR ENTRIES FROM THE LOCAL ONES
154 getKSPf(), data, &*locRhs.begin(), ADD_VALUES);
155 }
156
158 }
159
160private:
162 VectorDouble locRhs;
163};
164//! [OpDomainRhs]
165
166//! [OpBoundaryLhs]
167struct OpBoundaryLhs : public OpEdgeEle {
168public:
169 OpBoundaryLhs(std::string row_field_name, std::string col_field_name)
170 : OpEdgeEle(row_field_name, col_field_name, OpEdgeEle::OPROWCOL) {
171 sYmm = true;
172 }
173
174 MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
175 EntityType col_type, EntData &row_data,
176 EntData &col_data) {
178
179 const int nb_row_dofs = row_data.getIndices().size();
180 const int nb_col_dofs = col_data.getIndices().size();
181
182 if (nb_row_dofs && nb_col_dofs) {
183
184 locBoundaryLhs.resize(nb_row_dofs, nb_col_dofs, false);
185 locBoundaryLhs.clear();
186
187 // get (boundary) element length
188 const double edge = getMeasure();
189
190 // get number of integration points
191 const int nb_integration_points = getGaussPts().size2();
192 // get integration weights
193 auto t_w = getFTensor0IntegrationWeight();
194
195 // get base functions on row
196 auto t_row_base = row_data.getFTensor0N();
197
198 // START THE LOOP OVER INTEGRATION POINTS TO CALCULATE LOCAL MATRIX
199 for (int gg = 0; gg != nb_integration_points; gg++) {
200 const double a = t_w * edge;
201
202 for (int rr = 0; rr != nb_row_dofs; ++rr) {
203 // get base functions on column
204 auto t_col_base = col_data.getFTensor0N(gg, 0);
205
206 for (int cc = 0; cc != nb_col_dofs; cc++) {
207
208 locBoundaryLhs(rr, cc) += t_row_base * t_col_base * a;
209
210 // move to the next base functions on column
211 ++t_col_base;
212 }
213 // move to the next base functions on row
214 ++t_row_base;
215 }
216
217 // move to the weight of the next integration point
218 ++t_w;
219 }
220
221 // FILL VALUES OF LOCAL MATRIX ENTRIES TO THE GLOBAL MATRIX
223 getKSPB(), row_data, col_data, &locBoundaryLhs(0, 0), ADD_VALUES);
224 if (row_side != col_side || row_type != col_type) {
225 transLocBoundaryLhs.resize(nb_col_dofs, nb_row_dofs, false);
226 noalias(transLocBoundaryLhs) = trans(locBoundaryLhs);
228 getKSPB(), col_data, row_data, &transLocBoundaryLhs(0, 0),
229 ADD_VALUES);
230 }
231 }
232
234 }
235
236private:
238};
239//! [OpBoundaryLhs]
240
241//! [OpBoundaryRhs]
242struct OpBoundaryRhs : public OpEdgeEle {
243public:
244 OpBoundaryRhs(std::string field_name, ScalarFunc boundary_function)
246 boundaryFunc(boundary_function) {}
247
248 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
250
251 const int nb_dofs = data.getIndices().size();
252
253 if (nb_dofs) {
254
255 locBoundaryRhs.resize(nb_dofs, false);
256 locBoundaryRhs.clear();
257
258 // get (boundary) element length
259 const double edge = getMeasure();
260
261 // get number of integration points
262 const int nb_integration_points = getGaussPts().size2();
263 // get integration weights
264 auto t_w = getFTensor0IntegrationWeight();
265 // get coordinates at integration point
266 auto t_coords = getFTensor1CoordsAtGaussPts();
267
268 // get base function
269 auto t_base = data.getFTensor0N();
270
271 // START THE LOOP OVER INTEGRATION POINTS TO CALCULATE LOCAL VECTOR
272 for (int gg = 0; gg != nb_integration_points; gg++) {
273
274 const double a = t_w * edge;
275 double boundary_term =
276 boundaryFunc(t_coords(0), t_coords(1), t_coords(2));
277
278 for (int rr = 0; rr != nb_dofs; rr++) {
279
280 locBoundaryRhs[rr] += t_base * boundary_term * a;
281
282 // move to the next base function
283 ++t_base;
284 }
285
286 // move to the weight of the next integration point
287 ++t_w;
288 // move to the coordinates of the next integration point
289 ++t_coords;
290 }
291
292 // FILL VALUES OF LOCAL VECTOR ENTRIES TO THE GLOBAL VECTOR
294 getKSPf(), data, &*locBoundaryRhs.begin(), ADD_VALUES);
295 }
296
298 }
299
300private:
302 VectorDouble locBoundaryRhs;
303};
304//! [OpBoundaryRhs]
305}; // namespace Poisson2DNonhomogeneousOperators
306
307#endif //__POISSON2DNONHOMOGENEOUS_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 dofs 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.
OpBoundaryLhs(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.
OpBoundaryRhs(std::string field_name, ScalarFunc boundary_function)
OpDomainLhs(std::string row_field_name, std::string col_field_name)
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.
OpDomainRhs(std::string field_name, ScalarFunc source_term_function)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.