v0.14.0
Loading...
Searching...
No Matches
hdiv_check_approx_in_3d.cpp
Go to the documentation of this file.
1/**
2 * \file hdiv_check_approx_in_3d
3 * \example hdiv_check_approx_in_3d.cpp
4 *
5 * Approximate vectorial polynomial in 3D and check derivatives
6 *
7 */
8
9#include <MoFEM.hpp>
10
11using namespace MoFEM;
12
13static char help[] = "...\n\n";
14
15constexpr int BASE_DIM = 3;
16constexpr int SPACE_DIM = 3;
17
19using EleOp = Ele::UserDataOperator;
20
21/**
22 * @brief OPerator to integrate mass matrix for least square approximation
23 *
24 */
27
28/**
29 * @brief Operator to integrate the right hand side matrix for the problem
30 *
31 */
33 GAUSS>::OpSource<BASE_DIM, SPACE_DIM>;
34
35constexpr double a0 = 0.0;
36constexpr double a1 = 2.0;
37constexpr double a2 = -15.0 * a0;
38constexpr double a3 = -20.0 / 6 * a1;
39constexpr double a4 = 15.0 * a0;
40constexpr double a5 = a1;
41constexpr double a6 = -a0;
42struct ApproxFunctions {
43
44 static FTensor::Tensor1<double, BASE_DIM> fUn(const double x, const double y,
45 double z) {
47 // x
48 x + y*y + x*x*x,
49 // y
50 y + z*x + y*y*y,
51 // z
52 z + x*x + z*z*z);
53 }
54
56 diffFun(const double x, const double y, const double z) {
58 // x,x
59 1. + 3 * x * x,
60 // x,y
61 2. * y,
62 // x, z
63 0.,
64 // y,x
65 z,
66 // y,y
67 1. + 3 * y * y,
68 // y,z
69 x,
70 // z
71 2 * x, 0., 1. + 3 * z * z);
72 }
73
74}; // namespace ApproxFunctions
75
76struct OpCheckValsDiffVals : public EleOp {
77 boost::shared_ptr<MatrixDouble> ptrVals;
78 boost::shared_ptr<MatrixDouble> ptrGrad;
79
80 OpCheckValsDiffVals(boost::shared_ptr<MatrixDouble> ptr_vals,
81 boost::shared_ptr<MatrixDouble> ptr_grad)
82 : EleOp(NOSPACE, OPSPACE), ptrVals(ptr_vals), ptrGrad(ptr_grad) {}
83
87
88 MoFEMErrorCode doWork(int side, EntityType type,
91
92 const double eps = 1e-6;
93 const int nb_gauss_pts = getGaussPts().size2();
94
95 auto t_vals_from_op = getFTensor1FromMat<SPACE_DIM>(*ptrVals);
97
98 for (int gg = 0; gg != nb_gauss_pts; gg++) {
99 const double x = getCoordsAtGaussPts()(gg, 0);
100 const double y = getCoordsAtGaussPts()(gg, 1);
101 const double z = getCoordsAtGaussPts()(gg, 2);
102
103 // Check approximation
105 delta_val(i) = t_vals_from_op(i) - ApproxFunctions::fUn(x, y, z)(i);
106 double err_val = sqrt(delta_val(i) * delta_val(i));
107 MOFEM_LOG("SELF", Sev::verbose) << "Approximation error: " << err_val;
108
109 if (err_val > eps)
110 SETERRQ1(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID, "Wrong value %4.3e",
111 err_val);
112
113
115 delta_diff_val(i, j) =
116 t_grad_from_op(i, j) - ApproxFunctions::diffFun(x, y, z)(i, j);
117 double err_diff_val = sqrt(delta_diff_val(i, j) * delta_diff_val(i, j));
118 if (err_diff_val > eps)
119 SETERRQ1(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
120 "Wrong derivative of value %4.3e", err_diff_val);
121
122 ++t_vals_from_op;
123 ++t_grad_from_op;
124 }
126 }
127};
128
129int main(int argc, char *argv[]) {
130
131 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
132
133 try {
134
135 DMType dm_name = "DMMOFEM";
136 CHKERR DMRegister_MoFEM(dm_name);
137
138 moab::Core mb_instance;
139 moab::Interface &moab = mb_instance;
140
141 // Create MoFEM instance
142 MoFEM::Core core(moab);
143 MoFEM::Interface &m_field = core;
144
145 Simple *simple = m_field.getInterface<Simple>();
146 PipelineManager *pipeline_mng = m_field.getInterface<PipelineManager>();
147 CHKERR simple->getOptions();
148 CHKERR simple->loadFile();
149
150 // Declare elements
151 enum bases { AINSWORTH, DEMKOWICZ, LASBASETOP };
152 const char *list_bases[] = {"ainsworth", "demkowicz"};
153 PetscBool flg;
154 PetscInt choice_base_value = AINSWORTH;
155 CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list_bases,
156 LASBASETOP, &choice_base_value, &flg);
157
158 if (flg != PETSC_TRUE)
159 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "base not set");
161 if (choice_base_value == AINSWORTH)
163 else if (choice_base_value == DEMKOWICZ)
165
166 AinsworthOrderHooks::broken_nbfacetri_edge_hdiv = [](int p) { return p; };
167 AinsworthOrderHooks::broken_nbfacetri_face_hdiv = [](int p) { return p; };
171 return p;
172 };
173
174 CHKERR simple->addDomainField("FIELD1", HDIV, base, 1);
175 constexpr int order = 4;
176 CHKERR simple->setFieldOrder("FIELD1", order);
177 CHKERR simple->setUp();
178 auto dm = simple->getDM();
179
180 MatrixDouble vals, diff_vals;
181
182 auto assemble_matrices_and_vectors = [&]() {
184
186 pipeline_mng->getOpDomainRhsPipeline(), {HDIV});
187 pipeline_mng->getOpDomainRhsPipeline().push_back(
188 new OpDomainSource("FIELD1", ApproxFunctions::fUn));
189
191 pipeline_mng->getOpDomainLhsPipeline(), {HDIV});
192 pipeline_mng->getOpDomainLhsPipeline().push_back(
193
194 new OpDomainMass("FIELD1", "FIELD1",
195 [](double x, double, double) { return 1; })
196
197 );
198
199 auto integration_rule = [](int, int, int p_data) {
200 return 2 * p_data + 1;
201 };
202 CHKERR pipeline_mng->setDomainRhsIntegrationRule(integration_rule);
203 CHKERR pipeline_mng->setDomainLhsIntegrationRule(integration_rule);
204
206 };
207
208 auto solve_problem = [&] {
210 auto solver = pipeline_mng->createKSP();
211 CHKERR KSPSetFromOptions(solver);
212 CHKERR KSPSetUp(solver);
213
214 auto D = createDMVector(dm);
215 auto F = vectorDuplicate(D);
216
217 CHKERR KSPSolve(solver, F, D);
218 CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
219 CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
220 CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
222 };
223
224 auto check_solution = [&] {
226
227 pipeline_mng->getOpDomainLhsPipeline().clear();
228 pipeline_mng->getOpDomainRhsPipeline().clear();
229
230 auto ptr_values = boost::make_shared<MatrixDouble>();
231 auto ptr_grad = boost::make_shared<MatrixDouble>();
232
233 // Change H-curl to H-div in 2D, and apply Piola transform
235 pipeline_mng->getOpDomainRhsPipeline(), {HDIV});
236
237 // Calculate field values at integration points
238 pipeline_mng->getOpDomainRhsPipeline().push_back(
239 new OpCalculateHVecVectorField<BASE_DIM>("FIELD1", ptr_values));
240 // Gradient
241 pipeline_mng->getOpDomainRhsPipeline().push_back(
243 ptr_grad));
244 pipeline_mng->getOpDomainRhsPipeline().push_back(
245 new OpCheckValsDiffVals(ptr_values, ptr_grad));
246
247 CHKERR pipeline_mng->loopFiniteElements();
248
250 };
251
252 CHKERR assemble_matrices_and_vectors();
253 CHKERR solve_problem();
254 CHKERR check_solution();
255 }
257
259}
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69
int main()
static const double eps
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, FIELD_DIM > OpDomainSource
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, FIELD_DIM > OpDomainMass
#define CATCH_ERRORS
Catch errors.
FieldApproximationBase
approximation base
Definition definitions.h:58
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base nme:nme847.
Definition definitions.h:60
@ DEMKOWICZ_JACOBI_BASE
Definition definitions.h:66
@ NOSPACE
Definition definitions.h:83
@ HDIV
field with continuous normal traction
Definition definitions.h:87
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_IMPOSSIBLE_CASE
Definition definitions.h:35
@ MOFEM_ATOM_TEST_INVALID
Definition definitions.h:40
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
constexpr int order
@ F
auto integration_rule
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition DMMoFEM.cpp:523
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition DMMoFEM.cpp:43
auto createDMVector(DM dm)
Get smart vector from DM.
Definition DMMoFEM.hpp:1099
#define MOFEM_LOG(channel, severity)
Log.
constexpr double a5
constexpr double a2
static char help[]
constexpr double a4
constexpr int SPACE_DIM
FormsIntegrators< EleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< BASE_DIM, SPACE_DIM > OpDomainSource
Operator to integrate the right hand side matrix for the problem.
constexpr double a0
constexpr int BASE_DIM
constexpr double a3
MoFEM::PipelineManager::ElementsAndOpsByDim< SPACE_DIM >::DomainEle Ele
FormsIntegrators< EleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< BASE_DIM, SPACE_DIM > OpDomainMass
OPerator to integrate mass matrix for least square approximation.
constexpr double a1
constexpr double a6
double D
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
FTensor::Tensor1< FTensor::PackPtr< T *, S >, Tensor_Dim > getFTensor1FromMat(ublas::matrix< T, L, A > &data)
Get tensor rank 1 (vector) form data matrix.
FTensor::Tensor2< FTensor::PackPtr< double *, 1 >, Tensor_Dim1, Tensor_Dim2 > getFTensor2FromMat(MatrixDouble &data)
Get tensor rank 2 (matrix) form data matrix.
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
constexpr int SPACE_DIM
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpMass
[Only used with Hooke equation (linear material model)]
Definition seepage.cpp:56
static FTensor::Tensor2< double, BASE_DIM, SPACE_DIM > diffFun(const double x, const double y, const double z)
static FTensor::Tensor1< double, BASE_DIM > fUn(const double x, const double y, double z)
Add operators pushing bases from local to physical configuration.
static boost::function< int(int)> broken_nbvolumetet_edge_hdiv
Definition Hdiv.hpp:27
static boost::function< int(int)> broken_nbvolumetet_face_hdiv
Definition Hdiv.hpp:28
static boost::function< int(int)> broken_nbfacetri_face_hdiv
Definition Hdiv.hpp:26
static boost::function< int(int)> broken_nbvolumetet_volume_hdiv
Definition Hdiv.hpp:29
static boost::function< int(int)> broken_nbfacetri_edge_hdiv
Definition Hdiv.hpp:25
Core (interface) class.
Definition Core.hpp:82
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition Core.cpp:72
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:112
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
MatrixDouble & getCoordsAtGaussPts()
Gauss points and weight, matrix (nb. of points x 3)
@ OPSPACE
operator do Work is execute on space data
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
Get vector field for H-div approximation.
Calculate gradient of vector field.
PipelineManager interface.
Simple interface for fast problem set-up.
Definition Simple.hpp:27
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
FTensor::Index< 'i', 3 > i
FTENSOR_INDEX(SPACE_DIM, i)
boost::shared_ptr< MatrixDouble > ptrVals
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
FTENSOR_INDEX(SPACE_DIM, j)
FTENSOR_INDEX(SPACE_DIM, k)
FTensor::Index< 'j', 2 > j
OpCheckValsDiffVals(boost::shared_ptr< MatrixDouble > ptr_vals, boost::shared_ptr< MatrixDouble > ptr_grad)
boost::shared_ptr< MatrixDouble > ptrGrad
FTensor::Index< 'k', 2 > k