v0.15.0
Loading...
Searching...
No Matches
Classes | Typedefs | Functions | Variables
higher_derivatives.cpp File Reference
#include <MoFEM.hpp>

Go to the source code of this file.

Classes

struct  ElementsAndOps< DIM >
 
struct  ElementsAndOps< 2 >
 
struct  AtomTest
 
struct  AtomTest::CommonData
 Collected data use d by operator to evaluate errors for the test. More...
 
struct  AtomTest::OpError< FIELD_DIM >
 Operator to evaluate errors. More...
 

Typedefs

using DomainEle = ElementsAndOps< SPACE_DIM >::DomainEle
 Finite elenent type.
 
using DomainEleOp = DomainEle::UserDataOperator
 Finire element operator type.
 
using EntData = EntitiesFieldData::EntData
 Data on entities.
 
using OpDomainMass = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< BASE_DIM, FIELD_DIM >
 OPerator to integrate mass matrix for least square approximation.
 
using OpDomainSource = FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< BASE_DIM, FIELD_DIM >
 Operator to integrate the right hand side matrix for the problem.
 

Functions

int main (int argc, char *argv[])
 [Check results]
 

Variables

static char help [] = "...\n\n"
 
constexpr char FIELD_NAME [] = "U"
 
constexpr int BASE_DIM = 1
 
constexpr int FIELD_DIM = 1
 
constexpr int SPACE_DIM = 2
 
auto fun
 Function to approximate.
 
auto diff_fun
 Function derivative.
 
auto diff2_fun
 Function second derivative.
 

Typedef Documentation

◆ DomainEle

Finite elenent type.

Definition at line 28 of file higher_derivatives.cpp.

◆ DomainEleOp

Finire element operator type.

Examples
mofem/atom_tests/child_and_parent.cpp.

Definition at line 29 of file higher_derivatives.cpp.

◆ EntData

Data on entities.

Definition at line 31 of file higher_derivatives.cpp.

◆ OpDomainMass

OPerator to integrate mass matrix for least square approximation.

Examples
hanging_node_approx.cpp, higher_derivatives.cpp, and mofem/atom_tests/child_and_parent.cpp.

Definition at line 66 of file higher_derivatives.cpp.

◆ OpDomainSource

using OpDomainSource = FormsIntegrators<DomainEleOp>::Assembly< PETSC>::LinearForm<GAUSS>::OpSource<BASE_DIM, FIELD_DIM>

Operator to integrate the right hand side matrix for the problem.

Examples
hanging_node_approx.cpp, higher_derivatives.cpp, and mofem/atom_tests/child_and_parent.cpp.

Definition at line 73 of file higher_derivatives.cpp.

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)

[Check results]

[Register MoFEM discrete manager in PETSc]

[Register MoFEM discrete manager in PETSc

[Create MoAB]

< mesh database

< mesh database interface

[Create MoAB]

[Create MoFEM]

< finite element database

< finite element database insterface

[Create MoFEM]

[AtomTest]

[AtomTest]

Definition at line 397 of file higher_derivatives.cpp.

397 {
398
399 // Initialisation of MoFEM/PETSc and MOAB data structures
400 MoFEM::Core::Initialize(&argc, &argv, NULL, help);
401
402 try {
403
404 //! [Register MoFEM discrete manager in PETSc]
405 DMType dm_name = "DMMOFEM";
406 CHKERR DMRegister_MoFEM(dm_name);
407 //! [Register MoFEM discrete manager in PETSc
408
409 //! [Create MoAB]
410 moab::Core mb_instance; ///< mesh database
411 moab::Interface &moab = mb_instance; ///< mesh database interface
412 //! [Create MoAB]
413
414 //! [Create MoFEM]
415 MoFEM::Core core(moab); ///< finite element database
416 MoFEM::Interface &m_field = core; ///< finite element database insterface
417 //! [Create MoFEM]
418
419 //! [AtomTest]
420 AtomTest ex(m_field);
421 CHKERR ex.runProblem();
422 //! [AtomTest]
423 }
425
427}
#define CATCH_ERRORS
Catch errors.
#define CHKERR
Inline error check.
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition DMMoFEM.cpp:43
static char help[]
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.

Variable Documentation

◆ BASE_DIM

constexpr int BASE_DIM = 1
constexpr

Definition at line 17 of file higher_derivatives.cpp.

◆ diff2_fun

auto diff2_fun
Initial value:
= [](const double x, const double y, const double z) {
2 + 6 * x + 12 * pow(x, 2), 1.,
1., 2 + 6 * y + 12 * pow(y, 2)};
}

Function second derivative.

Examples
higher_derivatives.cpp.

Definition at line 55 of file higher_derivatives.cpp.

55 {
57 2 + 6 * x + 12 * pow(x, 2), 1.,
58
59 1., 2 + 6 * y + 12 * pow(y, 2)};
60};

◆ diff_fun

auto diff_fun
Initial value:
= [](const double x, const double y, const double z) {
2 * x + y + 3 * pow(x, 2) + 4 * pow(x, 3),
2 * y + x + 3 * pow(y, 2) + 4 * pow(y, 3)};
}

Function derivative.

Examples
higher_derivatives.cpp, prism_polynomial_approximation.cpp, and quad_polynomial_approximation.cpp.

Definition at line 45 of file higher_derivatives.cpp.

45 {
47 2 * x + y + 3 * pow(x, 2) + 4 * pow(x, 3),
48 2 * y + x + 3 * pow(y, 2) + 4 * pow(y, 3)};
49};

◆ FIELD_DIM

constexpr int FIELD_DIM = 1
constexpr

◆ FIELD_NAME

constexpr char FIELD_NAME[] = "U"
constexpr

◆ fun

auto fun
Initial value:
= [](const double x, const double y, const double z) {
return x * x + y * y + x * y + pow(x, 3) + pow(y, 3) + pow(x, 4) + pow(y, 4);
}

Function to approximate.

Examples
higher_derivatives.cpp.

Definition at line 37 of file higher_derivatives.cpp.

37 {
38 return x * x + y * y + x * y + pow(x, 3) + pow(y, 3) + pow(x, 4) + pow(y, 4);
39};

◆ help

char help[] = "...\n\n"
static

Definition at line 14 of file higher_derivatives.cpp.

◆ SPACE_DIM

constexpr int SPACE_DIM = 2
constexpr
Examples
hanging_node_approx.cpp.

Definition at line 19 of file higher_derivatives.cpp.