|
| v0.14.0
|
Go to the documentation of this file.
11 using namespace MoFEM;
13 static char help[] =
"...\n\n";
33 const int nb_gauss_pts = data.
getN().size1();
35 for (
int gg = 0; gg != nb_gauss_pts; gg++) {
36 const double val = getArea() * getGaussPts()(2, gg);
37 for (
int bb = 0; bb != nb_dofs; bb++) {
38 dIv += val * (t_diff_base_fun(0, 0) + t_diff_base_fun(1, 1) +
39 t_diff_base_fun(2, 2));
60 const int nb_gauss_pts = data.
getN().size1();
64 for (
int gg = 0; gg != nb_gauss_pts; gg++) {
65 const double val = getGaussPts()(1, gg);
66 for (
int bb = 0; bb != nb_dofs; bb++) {
67 fLux += val * t_base_fun(0);
76 int main(
int argc,
char *argv[]) {
85 PetscBool flg_file = PETSC_TRUE;
89 if (flg_file != PETSC_TRUE)
91 "*** ERROR -my_file (MESH FILE NEEDED)");
106 enum bases { AINSWORTH, DEMKOWICZ, LASBASETOP };
107 const char *list_bases[] = {
"ainsworth",
"demkowicz"};
109 PetscInt choice_base_value = AINSWORTH;
111 LASBASETOP, &choice_base_value, &flg);
112 if (flg != PETSC_TRUE)
115 if (choice_base_value == AINSWORTH)
117 else if (choice_base_value == DEMKOWICZ)
156 auto set_edge_elements_entities_on_mesh_skin = [&]() {
159 CHKERR moab.get_entities_by_dimension(0, 2, faces,
false);
162 CHKERR skin.find_skin(0, faces,
false, faces_skin);
167 CHKERR set_edge_elements_entities_on_mesh_skin();
186 auto rule = [&](
int,
int,
int p) {
return 2 * p; };
188 auto calculate_divergence = [&]() {
198 auto calculate_flux = [&]() {
207 const double div = calculate_divergence();
208 const double flux = calculate_flux();
210 PetscPrintf(PETSC_COMM_WORLD,
"Div = %4.3e Flux = %3.4e Error = %4.3e\n",
211 div, flux, div + flux);
213 constexpr
double tol = 1e-8;
214 if (std::abs(div + flux) >
tol)
216 "Test failed (div != flux) %3.4e != %3.4e", div, flux);
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Data on single entity (This is passed as argument to DataOperator::doWork)
virtual MoFEMErrorCode loop_finite_elements(const std::string problem_name, const std::string &fe_name, FEMethod &method, boost::shared_ptr< NumeredEntFiniteElement_multiIndex > fe_ptr=nullptr, MoFEMTypes bh=MF_EXIST, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr(), int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
int main(int argc, char *argv[])
Problem manager is used to build and partition problems.
RuleHookFun getRuleHook
Hook to get rule.
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
virtual MoFEMErrorCode add_ents_to_field_by_type(const Range &ents, const EntityType type, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
Deprecated interface functions.
DeprecatedCoreInterface Interface
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2DiffN(FieldApproximationBase base)
Get derivatives of base functions for Hdiv space.
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string &name, const bool recursive=true)=0
add entities to finite element
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
#define CHKERR
Inline error check.
friend class UserDataOperator
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual moab::Interface & get_moab()=0
implementation of Data Operators for Forces and Sources
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string name_row)=0
set field col which finite element use
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
default operator for TRI element
const VectorInt & getIndices() const
Get global indices of dofs on entity.
friend class UserDataOperator
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
FTensor::Index< 'i', SPACE_DIM > i
MoFEMErrorCode partitionFiniteElements(const std::string name, bool part_from_moab=false, int low_proc=-1, int hi_proc=-1, int verb=VERBOSE)
partition finite elements
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
Add operators pushing bases from local to physical configuration.
default operator for EDGE element
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
#define CATCH_ERRORS
Catch errors.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
OpDivergence(double &div)
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
@ HCURL
field with continuous tangents
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
FieldApproximationBase
approximation base
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
virtual MoFEMErrorCode modify_problem_add_finite_element(const std::string name_problem, const std::string &fe_name)=0
add finite element to problem, this add entities assigned to finite element to a particular problem
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
@ MOFEM_ATOM_TEST_INVALID
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
FTensor::Index< 'i', 3 > i
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
virtual MoFEMErrorCode add_field(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_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...