v0.14.0
Functions | Variables
dm_mofem.cpp File Reference

Atom test for Data Manager Interface in MoFEM. More...

#include <MoFEM.hpp>

Go to the source code of this file.

Functions

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

Variables

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

Detailed Description

Atom test for Data Manager Interface in MoFEM.

Definition in file dm_mofem.cpp.

Function Documentation

◆ main()

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

Definition at line 13 of file dm_mofem.cpp.

13  {
14 
15  // initialize petsc
16  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
17 
18  try {
19 
20  PetscBool flg = PETSC_TRUE;
21  char mesh_file_name[255];
22 #if PETSC_VERSION_GE(3, 6, 4)
23  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
24  255, &flg);
25 #else
26  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
27  mesh_file_name, 255, &flg);
28 #endif
29  if (flg != PETSC_TRUE) {
30  SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
31  }
32  const char *option;
33  option = "";
34 
35  // register new dm type, i.e. mofem
36  DMType dm_name = "DMMOFEM";
37  CHKERR DMRegister_MoFEM(dm_name);
38 
39  // Create dm instance
40  DM dm;
41  CHKERR DMCreate(PETSC_COMM_WORLD, &dm);
42  CHKERR DMSetType(dm, dm_name);
43 
44  // read mesh and create moab and mofem data structures
45  moab::Core mb_instance;
46  moab::Interface &moab = mb_instance;
47  CHKERR moab.load_file(mesh_file_name, 0, option);
48  MoFEM::Core core(moab);
49  MoFEM::Interface &m_field = core;
50 
51  EntityHandle root_set = moab.get_root_set();
52  // add all entities to database, all of them will be used
53  BitRefLevel bit_level0;
54  bit_level0.set(0);
55  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
56  root_set, 3, bit_level0);
57  // define & build field
58  const int field_rank = 1;
59  CHKERR m_field.add_field("FIELD", H1, AINSWORTH_LEGENDRE_BASE, field_rank);
60  // add entities to field
61  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FIELD");
62  // set app. order
63  int order = 4;
64  CHKERR m_field.set_field_order(root_set, MBTET, "FIELD", order);
65  CHKERR m_field.set_field_order(root_set, MBTRI, "FIELD", order);
66  CHKERR m_field.set_field_order(root_set, MBEDGE, "FIELD", order);
67  CHKERR m_field.set_field_order(root_set, MBVERTEX, "FIELD", 1);
68  // build data structures for fields
69  CHKERR m_field.build_fields();
70 
71  // define & build finite elements
72  CHKERR m_field.add_finite_element("FE");
73  // Define rows/cols and element data
74  CHKERR m_field.modify_finite_element_add_field_row("FE", "FIELD");
75  CHKERR m_field.modify_finite_element_add_field_col("FE", "FIELD");
76  CHKERR m_field.modify_finite_element_add_field_data("FE", "FIELD");
77  // add entities to finite element
78  CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBTET, "FE");
79  // build finite elements
80  CHKERR m_field.build_finite_elements();
81  // build adjacencies
82  CHKERR m_field.build_adjacencies(bit_level0);
83 
84  // set dm data structure which created mofem data structures
85  CHKERR DMMoFEMCreateMoFEM(dm, &m_field, dm_name, bit_level0);
86  CHKERR DMSetFromOptions(dm);
87  CHKERR DMMoFEMAddElement(dm, "FE");
88  CHKERR DMSetUp(dm);
89 
90  Mat m;
91  Vec l, g;
92 
93  CHKERR DMCreateGlobalVector(dm, &g);
94  CHKERR DMCreateLocalVector(dm, &l);
95  CHKERR DMCreateMatrix(dm, &m);
96 
97  // glob loc
98  CHKERR VecSet(g, 1.1);
99  CHKERR VecGhostUpdateBegin(g, INSERT_VALUES, SCATTER_FORWARD);
100  CHKERR VecGhostUpdateEnd(g, INSERT_VALUES, SCATTER_FORWARD);
101 
102  CHKERR DMGlobalToLocalBegin(dm, g, ADD_VALUES, l);
103  CHKERR DMGlobalToLocalEnd(dm, g, ADD_VALUES, l);
104 
105  // loc glob
106  CHKERR DMLocalToGlobalBegin(dm, l, ADD_VALUES, g);
107  CHKERR DMLocalToGlobalEnd(dm, l, ADD_VALUES, g);
108 
109  PetscViewer viewer;
110  CHKERR PetscViewerASCIIOpen(PETSC_COMM_WORLD, "dm_mofem.txt", &viewer);
111  const double chop = 1e-4;
112  CHKERR VecChop(g, chop);
113  VecView(g, viewer);
114  CHKERR PetscViewerDestroy(&viewer);
115 
116  CHKERR VecDestroy(&g);
117  CHKERR VecDestroy(&l);
118  CHKERR MatDestroy(&m);
119 
120  // destroy dm
121  CHKERR DMDestroy(&dm);
122  }
123  CATCH_ERRORS;
124 
125  // finish work cleaning memory, getting statistics, ect.
127 
128  return 0;
129 }

Variable Documentation

◆ help

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

Definition at line 11 of file dm_mofem.cpp.

MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
g
constexpr double g
Definition: shallow_wave.cpp:63
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
EntityHandle
MoFEM::CoreInterface::modify_finite_element_add_field_row
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
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MoFEM::CoreInterface::add_ents_to_field_by_type
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.
MoFEM::DMMoFEMAddElement
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Definition: DMMoFEM.cpp:497
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
MoFEM::CoreInterface::add_ents_to_finite_element_by_type
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
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::add_finite_element
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
MoFEM::CoreInterface::modify_finite_element_add_field_col
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
MoFEM::CoreInterface::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:43
MoFEM::DMMoFEMCreateMoFEM
PetscErrorCode DMMoFEMCreateMoFEM(DM dm, MoFEM::Interface *m_field_ptr, const char problem_name[], const MoFEM::BitRefLevel bit_level, const MoFEM::BitRefLevel bit_mask=MoFEM::BitRefLevel().set())
Must be called by user to set MoFEM data structures.
Definition: DMMoFEM.cpp:114
MoFEM::CoreInterface::modify_finite_element_add_field_data
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:22
help
static char help[]
Definition: dm_mofem.cpp:11
MoFEM::CoreTmp< 0 >::Initialize
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
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::CoreInterface::build_fields
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
MoFEM::CoreInterface::build_adjacencies
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
MoFEM::CoreInterface::set_field_order
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.
MoFEM::CoreInterface::add_field
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.
l
FTensor::Index< 'l', 3 > l
Definition: matrix_function.cpp:21