v0.13.2
Loading...
Searching...
No Matches
poisson_3d_homogeneous.hpp
Go to the documentation of this file.
1/**
2 * @file poisson_3d_homogeneous.hpp
3 * @example poisson_3d_homogeneous.hpp
4 * @brief Operator for 3D poisson problem example
5 * @date 2023-01-31
6 *
7 * @copyright Copyright (c) 2023
8 *
9 */
10
11#ifndef __POISSON3DHOMOGENEOUS_HPP__
12#define __POISSON3DHOMOGENEOUS_HPP__
13
14#include <stdlib.h>
16
18
20
22
24
26
27const double body_source = 1.;
28
30public:
31 OpDomainLhsMatrixK(std::string row_field_name, std::string col_field_name)
32 : OpVolEle(row_field_name, col_field_name, OpVolEle::OPROWCOL) {
33 sYmm = true;
34 }
35
36 MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
37 EntityType col_type, EntData &row_data,
38 EntData &col_data) {
40
41 const int nb_row_dofs = row_data.getIndices().size();
42 const int nb_col_dofs = col_data.getIndices().size();
43
44 if (nb_row_dofs && nb_col_dofs) {
45
46 locLhs.resize(nb_row_dofs, nb_col_dofs, false);
47 locLhs.clear();
48
49 // get element volume
50 const double volume = getMeasure();
51
52 // get number of integration points
53 const int nb_integration_points = getGaussPts().size2();
54 // get integration weights
56
57 // get derivatives of base functions on row
58 auto t_row_diff_base = row_data.getFTensor1DiffN<3>();
59
60 // START THE LOOP OVER INTEGRATION POINTS TO CALCULATE LOCAL MATRIX
61 for (int gg = 0; gg != nb_integration_points; gg++) {
62 const double a = t_w * volume;
63
64 for (int rr = 0; rr != nb_row_dofs; ++rr) {
65 // get derivatives of base functions on column
66 auto t_col_diff_base = col_data.getFTensor1DiffN<3>(gg, 0);
67
68 for (int cc = 0; cc != nb_col_dofs; cc++) {
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 CHKERR MatSetValues(getKSPB(), row_data, col_data, &locLhs(0, 0),
85 ADD_VALUES);
86 if (row_side != col_side || row_type != col_type) {
87 transLocLhs.resize(nb_col_dofs, nb_row_dofs, false);
88 noalias(transLocLhs) = trans(locLhs);
89 CHKERR MatSetValues(getKSPB(), col_data, row_data, &transLocLhs(0, 0),
90 ADD_VALUES);
91 }
92 }
93
95 }
96
97private:
99};
100
102public:
105
106 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
108
109 const int nb_dofs = data.getIndices().size();
110
111 if (nb_dofs) {
112 locRhs.resize(nb_dofs, false);
113 locRhs.clear();
114
115 // get element volume
116 const double volume = getMeasure();
117
118 // get number of integration points
119 const int nb_integration_points = getGaussPts().size2();
120 // get integration weights
121 auto t_w = getFTensor0IntegrationWeight();
122
123 // get base function
124 auto t_base = data.getFTensor0N();
125
126 // START THE LOOP OVER INTEGRATION POINTS TO CALCULATE LOCAL VECTOR
127 for (int gg = 0; gg != nb_integration_points; gg++) {
128 const double a = t_w * volume;
129
130 for (int rr = 0; rr != nb_dofs; rr++) {
131 locRhs[rr] += t_base * body_source * a;
132
133 // move to the next base function
134 ++t_base;
135 }
136
137 // move to the weight of the next integration point
138 ++t_w;
139 }
140
141 // FILL VALUES OF LOCAL VECTOR ENTRIES TO THE GLOBAL VECTOR
142
143 // Ignoring DOFs on boundary (index -1)
144 CHKERR VecSetOption(getKSPf(), VEC_IGNORE_NEGATIVE_INDICES, PETSC_TRUE);
145 CHKERR VecSetValues(getKSPf(), data, &*locRhs.begin(), ADD_VALUES);
146 }
147
149 }
150
151private:
153};
154
155}; // namespace Poisson3DHomogeneousOperators
156
157#endif //__POISSON3DHOMOGENEOUS_HPP__
constexpr double a
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
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.
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 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
OpDomainLhsMatrixK(std::string row_field_name, std::string col_field_name)
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.
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.