v0.15.0
Loading...
Searching...
No Matches
dg_projection.cpp
Go to the documentation of this file.
1/**
2 * @file dg_projection.cpp
3 * @example atomic_tests/dg_projection.cpp
4 *
5 * Testing DG projection operators
6 *
7 */
8
9#include <MoFEM.hpp>
10
11using namespace MoFEM;
12
13static char help[] = "...\n\n";
14
15constexpr char FIELD_NAME[] = "U";
16constexpr int BASE_DIM = 1;
17constexpr int FIELD_DIM = 1;
18constexpr int SPACE_DIM = 2;
19constexpr int order = 2;
20
21template <int DIM> struct ElementsAndOps {};
22
23template <> struct ElementsAndOps<2> {
26};
27
28using DomainEle = ElementsAndOps<SPACE_DIM>::DomainEle; ///< Finite elenent type
30 DomainEle::UserDataOperator; ///< Finire element operator type
31using EntData = EntitiesFieldData::EntData; ///< Data on entities
32
33/**
34 * @brief Function to approximate
35 *
36 */
37auto fun = [](const double x, const double y, const double z) {
38 return x + y + x * x + y * y;
39};
40
41/**
42 * @brief OPerator to integrate mass matrix for least square approximation
43 *
44 */
47
48/**
49 * @brief Operator to integrate the right hand side matrix for the problem
50 *
51 */
54
55struct AtomTest {
56
57 AtomTest(MoFEM::Interface &m_field) : mField(m_field) {}
58
60
61private:
64
70
71 /**
72 * @brief Collected data use d by operator to evaluate errors for the test
73 *
74 */
75 struct CommonData {
76 boost::shared_ptr<MatrixDouble> invJacPtr;
77 boost::shared_ptr<VectorDouble> approxVals;
78 boost::shared_ptr<MatrixDouble> approxGradVals;
79 boost::shared_ptr<MatrixDouble> approxHessianVals;
81 };
82
83 /**
84 * @brief Operator to evaluate errors
85 *
86 */
87 struct OpError;
88};
89
90/**
91 * @brief Operator to evaluate errors
92 *
93 */
94struct AtomTest::OpError : public DomainEleOp {
95 boost::shared_ptr<CommonData> commonDataPtr;
96 OpError(boost::shared_ptr<MatrixDouble> data_ptr)
97 : DomainEleOp(NOSPACE, OPSPACE), dataPtr(data_ptr) {}
98
99 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
101
102 const int nb_integration_pts = getGaussPts().size2();
103 auto t_val = getFTensor1FromMat<1>(
104 *(dataPtr)); // get function approximation at gauss pts
105 auto t_coords = getFTensor1CoordsAtGaussPts(); // get coordinates of
106 // integration points
107
108 for (int gg = 0; gg != nb_integration_pts; ++gg) {
109
110 // Calculate errors
111
112 double diff = t_val(0) - fun(t_coords(0), t_coords(1), t_coords(2));
113 constexpr double eps = 1e-8;
114 if (std::abs(diff) > eps) {
115 MOFEM_LOG("SELF", Sev::error) << "Wrong function value " << diff;
116 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
117 "Wrong function value");
118 }
119
120 // move data to next integration point
121 ++t_val;
122 ++t_coords;
123 }
124
125 MOFEM_LOG("SELF", Sev::noisy) << "All is OK";
126
128 }
129
130 private:
131 boost::shared_ptr<MatrixDouble> dataPtr;
132};
133
134//! [Run programme]
143}
144//! [Run programme]
145
146//! [Read mesh]
149
153
155}
156//! [Read mesh]
157
158//! [Set up problem]
161 // Add field
166
168}
169//! [Set up problem]
170
171//! [Push operators to pipeline]
174
175 auto rule = [](int, int, int p) -> int { return 2 * p; };
176
178 CHKERR pipeline_mng->setDomainLhsIntegrationRule(rule);
179 CHKERR pipeline_mng->setDomainRhsIntegrationRule(rule);
180
181 auto beta = [](const double, const double, const double) { return 1; };
182 pipeline_mng->getOpDomainLhsPipeline().push_back(
184 pipeline_mng->getOpDomainRhsPipeline().push_back(
186
188}
189//! [Push operators to pipeline]
190
191//! [Solve]
195
196 MOFEM_LOG("WORLD", Sev::inform) << "Solve problem";
197
198 auto solver = pipeline_mng->createKSP();
199 CHKERR KSPSetFromOptions(solver);
200 CHKERR KSPSetUp(solver);
201
202 auto dm = simpleInterface->getDM();
203 auto D = createDMVector(dm);
204 auto F = vectorDuplicate(D);
205
206 CHKERR KSPSolve(solver, F, D);
207 CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
208 CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
209 CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
211}
212
213//! [Check results]
217 auto pipeline_mng = mField.getInterface<PipelineManager>();
218 pipeline_mng->getDomainLhsFE().reset();
219 pipeline_mng->getDomainRhsFE().reset();
220 pipeline_mng->getOpDomainRhsPipeline().clear();
221
222 auto rule = [](int, int, int p) -> int { return 2 * p + 1; };
223 CHKERR pipeline_mng->setDomainRhsIntegrationRule(
224 rule); // set integration rule
225
226 auto entity_data_l2 = boost::make_shared<EntitiesFieldData>(
227 MBENTITYSET); // entity data shared between
228 // physical and post proc
229 // elements
230 auto mass_ptr = boost::make_shared<MatrixDouble>(); // integrated mass matrix
231 // of post proc element
232 auto coeffs_ptr =
233 boost::make_shared<MatrixDouble>(); // vector of coeffs shared between
234 // physical and post proc elements
235 auto data_ptr =
236 boost::make_shared<MatrixDouble>(); // data stored at integration points
237 // of the physical element and
238 // evaluated at integration points of
239 // the post proc element
240
241 auto op_this =
242 new OpLoopThis<DomainEle>(mField, simple->getDomainFEName(), Sev::noisy);
243 pipeline_mng->getOpDomainRhsPipeline().push_back(op_this); // 1
244 pipeline_mng->getOpDomainRhsPipeline().push_back(new OpDGProjectionEvaluation(
245 data_ptr, coeffs_ptr, entity_data_l2, AINSWORTH_LEGENDRE_BASE,
246 L2)); // 5
247 pipeline_mng->getOpDomainRhsPipeline().push_back(
248 new OpError(data_ptr)); // 6
249
250 auto fe_physics_ptr = op_this->getThisFEPtr();
251 fe_physics_ptr->getRuleHook = [](int, int, int p) { return 2 * p; };
252
253 fe_physics_ptr->getOpPtrVector().push_back(new OpDGProjectionMassMatrix(
254 order, mass_ptr, entity_data_l2, AINSWORTH_LEGENDRE_BASE, L2)); // 2
255 fe_physics_ptr->getOpPtrVector().push_back(
257 data_ptr)); // 3
258 fe_physics_ptr->getOpPtrVector().push_back(
259 new OpDGProjectionCoefficients(data_ptr, coeffs_ptr, mass_ptr,
260 entity_data_l2, AINSWORTH_LEGENDRE_BASE,
261 L2)); // 4
262
263 CHKERR pipeline_mng->loopFiniteElements();
264
266}
267//! [Check results]
268
269int main(int argc, char *argv[]) {
270
271 // Initialisation of MoFEM/PETSc and MOAB data structures
272 MoFEM::Core::Initialize(&argc, &argv, NULL, help);
273
274 try {
275
276 //! [Register MoFEM discrete manager in PETSc]
277 DMType dm_name = "DMMOFEM";
278 CHKERR DMRegister_MoFEM(dm_name);
279 //! [Register MoFEM discrete manager in PETSc
280
281 //! [Create MoAB]
282 moab::Core mb_instance; ///< mesh database
283 moab::Interface &moab = mb_instance; ///< mesh database interface
284 //! [Create MoAB]
285
286 //! [Create MoFEM]
287 MoFEM::Core core(moab); ///< finite element database
288 MoFEM::Interface &m_field = core; ///< finite element database insterface
289 //! [Create MoFEM]
290
291 //! [AtomTest]
292 AtomTest ex(m_field);
293 CHKERR ex.runProblem();
294 //! [AtomTest]
295 }
297
299}
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
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, FIELD_DIM > OpDomainMass
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
@ L2
field with C-1 continuity
Definition definitions.h:88
@ H1
continuous field
Definition definitions.h:85
@ NOSPACE
Definition definitions.h:83
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ 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.
static char help[]
constexpr char FIELD_NAME[]
auto fun
Function to approximate.
constexpr int SPACE_DIM
constexpr int FIELD_DIM
constexpr int BASE_DIM
constexpr int order
@ F
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:514
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:1119
boost::ptr_deque< UserDataOperator > & getOpDomainLhsPipeline()
Get the Op Domain Lhs Pipeline object.
SmartPetscObj< KSP > createKSP(SmartPetscObj< DM > dm=nullptr)
Create KSP (linear) solver.
boost::ptr_deque< UserDataOperator > & getOpDomainRhsPipeline()
Get the Op Domain Rhs Pipeline object.
#define MOFEM_LOG(channel, severity)
Log.
constexpr char FIELD_NAME[]
auto fun
Function to approximate.
constexpr int FIELD_DIM
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< BASE_DIM, FIELD_DIM > OpDomainMass
OPerator to integrate mass matrix for least square approximation.
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< BASE_DIM, FIELD_DIM > OpDomainSource
Operator to integrate the right hand side matrix for the problem.
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.
Collected data use d by operator to evaluate errors for the test.
boost::shared_ptr< VectorDouble > approxVals
SmartPetscObj< Vec > L2Vec
boost::shared_ptr< MatrixDouble > invJacPtr
boost::shared_ptr< MatrixDouble > approxHessianVals
boost::shared_ptr< MatrixDouble > approxGradVals
Operator to evaluate errors.
boost::shared_ptr< CommonData > commonDataPtr
OpError(boost::shared_ptr< MatrixDouble > data_ptr)
boost::shared_ptr< MatrixDouble > dataPtr
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Simple * simpleInterface
MoFEMErrorCode checkResults()
[Check results]
MoFEMErrorCode solveSystem()
AtomTest(MoFEM::Interface &m_field)
MoFEMErrorCode setupProblem()
MoFEMErrorCode readMesh()
MoFEMErrorCode assembleSystem()
MoFEM::Interface & mField
MoFEMErrorCode runProblem()
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:118
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
Get values at integration pts for tensor field rank 1, i.e. vector field.
Execute "this" element in the operator.
PipelineManager interface.
MoFEM::FaceElementForcesAndSourcesCore FaceEle
MoFEMErrorCode setDomainRhsIntegrationRule(RuleHookFun rule)
MoFEMErrorCode setDomainLhsIntegrationRule(RuleHookFun rule)
Simple interface for fast problem set-up.
Definition Simple.hpp:27
MoFEMErrorCode addDomainField(const std::string name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_ZERO, int verb=-1)
Add field on domain.
Definition Simple.cpp:261
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
Definition Simple.cpp:191
MoFEMErrorCode getOptions()
get options
Definition Simple.cpp:180
MoFEMErrorCode getDM(DM *dm)
Get DM.
Definition Simple.cpp:800
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Definition Simple.cpp:575
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
Definition Simple.cpp:736
intrusive_ptr for managing petsc objects
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.