v0.14.0
Loading...
Searching...
No Matches
dg_projection.cpp
Go to the documentation of this file.
1/**
2 * \example dg_projection.cpp
3 *
4 * Testing DG projection operators
5 *
6 */
7
8#include <MoFEM.hpp>
9
10using namespace MoFEM;
11
12static char help[] = "...\n\n";
13
14constexpr char FIELD_NAME[] = "U";
15constexpr int BASE_DIM = 1;
16constexpr int FIELD_DIM = 1;
17constexpr int SPACE_DIM = 2;
18constexpr int order = 2;
19
20template <int DIM> struct ElementsAndOps {};
21
22template <> struct ElementsAndOps<2> {
25};
26
27using DomainEle = ElementsAndOps<SPACE_DIM>::DomainEle; ///< Finite elenent type
29 DomainEle::UserDataOperator; ///< Finire element operator type
30using EntData = EntitiesFieldData::EntData; ///< Data on entities
31
32/**
33 * @brief Function to approximate
34 *
35 */
36auto fun = [](const double x, const double y, const double z) {
37 return x + y + x * x + y * y;
38};
39
40/**
41 * @brief OPerator to integrate mass matrix for least square approximation
42 *
43 */
46
47/**
48 * @brief Operator to integrate the right hand side matrix for the problem
49 *
50 */
53
54struct AtomTest {
55
56 AtomTest(MoFEM::Interface &m_field) : mField(m_field) {}
57
59
60private:
63
69
70 /**
71 * @brief Collected data use d by operator to evaluate errors for the test
72 *
73 */
74 struct CommonData {
75 boost::shared_ptr<MatrixDouble> invJacPtr;
76 boost::shared_ptr<VectorDouble> approxVals;
77 boost::shared_ptr<MatrixDouble> approxGradVals;
78 boost::shared_ptr<MatrixDouble> approxHessianVals;
80 };
81
82 /**
83 * @brief Operator to evaluate errors
84 *
85 */
86 struct OpError;
87};
88
89/**
90 * @brief Operator to evaluate errors
91 *
92 */
93struct AtomTest::OpError : public DomainEleOp {
94 boost::shared_ptr<CommonData> commonDataPtr;
95 OpError(boost::shared_ptr<MatrixDouble> data_ptr)
96 : DomainEleOp(NOSPACE, OPSPACE), dataPtr(data_ptr) {}
97
98 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
100
101 const int nb_integration_pts = getGaussPts().size2();
102 auto t_val = getFTensor1FromMat<1>(
103 *(dataPtr)); // get function approximation at gauss pts
104 auto t_coords = getFTensor1CoordsAtGaussPts(); // get coordinates of
105 // integration points
106
107 for (int gg = 0; gg != nb_integration_pts; ++gg) {
108
109 // Calculate errors
110
111 double diff = t_val(0) - fun(t_coords(0), t_coords(1), t_coords(2));
112 constexpr double eps = 1e-8;
113 if (std::abs(diff) > eps) {
114 MOFEM_LOG("SELF", Sev::error) << "Wrong function value " << diff;
115 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
116 "Wrong function value");
117 }
118
119 // move data to next integration point
120 ++t_val;
121 ++t_coords;
122 }
123
124 MOFEM_LOG("SELF", Sev::noisy) << "All is OK";
125
127 }
128
129 private:
130 boost::shared_ptr<MatrixDouble> dataPtr;
131};
132
133//! [Run programme]
142}
143//! [Run programme]
144
145//! [Read mesh]
148
152
154}
155//! [Read mesh]
156
157//! [Set up problem]
160 // Add field
165
167}
168//! [Set up problem]
169
170//! [Push operators to pipeline]
173
174 auto rule = [](int, int, int p) -> int { return 2 * p; };
175
177 CHKERR pipeline_mng->setDomainLhsIntegrationRule(rule);
178 CHKERR pipeline_mng->setDomainRhsIntegrationRule(rule);
179
180 auto beta = [](const double, const double, const double) { return 1; };
181 pipeline_mng->getOpDomainLhsPipeline().push_back(
183 pipeline_mng->getOpDomainRhsPipeline().push_back(
185
187}
188//! [Push operators to pipeline]
189
190//! [Solve]
194
195 MOFEM_LOG("WORLD", Sev::inform) << "Solve problem";
196
197 auto solver = pipeline_mng->createKSP();
198 CHKERR KSPSetFromOptions(solver);
199 CHKERR KSPSetUp(solver);
200
201 auto dm = simpleInterface->getDM();
202 auto D = createDMVector(dm);
203 auto F = vectorDuplicate(D);
204
205 CHKERR KSPSolve(solver, F, D);
206 CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
207 CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
208 CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
210}
211
212//! [Check results]
216 auto pipeline_mng = mField.getInterface<PipelineManager>();
217 pipeline_mng->getDomainLhsFE().reset();
218 pipeline_mng->getDomainRhsFE().reset();
219 pipeline_mng->getOpDomainRhsPipeline().clear();
220
221 auto rule = [](int, int, int p) -> int { return 2 * p + 1; };
222 CHKERR pipeline_mng->setDomainRhsIntegrationRule(
223 rule); // set integration rule
224
225 auto entity_data_l2 = boost::make_shared<EntitiesFieldData>(
226 MBENTITYSET); // entity data shared between
227 // physical and post proc
228 // elements
229 auto mass_ptr = boost::make_shared<MatrixDouble>(); // integrated mass matrix
230 // of post proc element
231 auto coeffs_ptr =
232 boost::make_shared<MatrixDouble>(); // vector of coeffs shared between
233 // physical and post proc elements
234 auto data_ptr =
235 boost::make_shared<MatrixDouble>(); // data stored at integration points
236 // of the physical element and
237 // evaluated at integration points of
238 // the post proc element
239
240 auto op_this =
241 new OpLoopThis<DomainEle>(mField, simple->getDomainFEName(), Sev::noisy);
242 pipeline_mng->getOpDomainRhsPipeline().push_back(op_this); // 1
243 pipeline_mng->getOpDomainRhsPipeline().push_back(new OpDGProjectionEvaluation(
244 data_ptr, coeffs_ptr, entity_data_l2, AINSWORTH_LEGENDRE_BASE,
245 L2)); // 5
246 pipeline_mng->getOpDomainRhsPipeline().push_back(
247 new OpError(data_ptr)); // 6
248
249 auto fe_physics_ptr = op_this->getThisFEPtr();
250 fe_physics_ptr->getRuleHook = [](int, int, int p) { return 2 * p; };
251
252 fe_physics_ptr->getOpPtrVector().push_back(new OpDGProjectionMassMatrix(
253 order, mass_ptr, entity_data_l2, AINSWORTH_LEGENDRE_BASE, L2)); // 2
254 fe_physics_ptr->getOpPtrVector().push_back(
256 data_ptr)); // 3
257 fe_physics_ptr->getOpPtrVector().push_back(
258 new OpDGProjectionCoefficients(data_ptr, coeffs_ptr, mass_ptr,
259 entity_data_l2, AINSWORTH_LEGENDRE_BASE,
260 L2)); // 4
261
262 CHKERR pipeline_mng->loopFiniteElements();
263
265}
266//! [Check results]
267
268int main(int argc, char *argv[]) {
269
270 // Initialisation of MoFEM/PETSc and MOAB data structures
271 MoFEM::Core::Initialize(&argc, &argv, NULL, help);
272
273 try {
274
275 //! [Register MoFEM discrete manager in PETSc]
276 DMType dm_name = "DMMOFEM";
277 CHKERR DMRegister_MoFEM(dm_name);
278 //! [Register MoFEM discrete manager in PETSc
279
280 //! [Create MoAB]
281 moab::Core mb_instance; ///< mesh database
282 moab::Interface &moab = mb_instance; ///< mesh database interface
283 //! [Create MoAB]
284
285 //! [Create MoFEM]
286 MoFEM::Core core(moab); ///< finite element database
287 MoFEM::Interface &m_field = core; ///< finite element database insterface
288 //! [Create MoFEM]
289
290 //! [AtomTest]
291 AtomTest ex(m_field);
292 CHKERR ex.runProblem();
293 //! [AtomTest]
294 }
296
298}
static Index< 'p', 3 > p
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition: acoustic.cpp:69
int main()
Definition: adol-c_atom.cpp:46
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.
Definition: definitions.h:372
@ 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 ...
Definition: definitions.h:346
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
#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
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:509
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:47
auto createDMVector(DM dm)
Get smart vector from DM.
Definition: DMMoFEM.hpp:1003
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.
Definition: LogManager.hpp:308
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.
Definition: Exceptions.hpp:56
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)
Operator for linear form, usually to calculate values on right hand side.
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:112
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
@ OPSPACE
operator do Work is execute on space data
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
Get values at integration pts for tensor filed 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 loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
Definition: Simple.cpp:194
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:264
MoFEMErrorCode getOptions()
get options
Definition: Simple.cpp:180
MoFEMErrorCode getDM(DM *dm)
Get DM.
Definition: Simple.cpp:667
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Definition: Simple.cpp:473
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
Definition: Simple.cpp:611
intrusive_ptr for managing petsc objects
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.