v0.14.0
Loading...
Searching...
No Matches
dm_create_subdm.cpp
Go to the documentation of this file.
1/** \file dm_create_subdm.cpp
2 \example dm_create_subdm.cpp
3 \brief Atom test for Data Manager Interface and create sub-problem
4*/
5
6#include <MoFEM.hpp>
7
8using namespace MoFEM;
9
10static char help[] = "...\n\n";
11
12static const bool debug = false;
13
14int main(int argc, char *argv[]) {
15
16 // initialize petsc
17 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
18
19 try {
20
21 PetscBool flg = PETSC_TRUE;
22 char mesh_file_name[255];
23#if PETSC_VERSION_GE(3, 6, 4)
24 CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
25 255, &flg);
26#else
27 CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
28 mesh_file_name, 255, &flg);
29#endif
30 if (flg != PETSC_TRUE) {
31 SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
32 }
33
34 // register new dm type, i.e. mofem
35 DMType dm_name = "DMMOFEM";
36 CHKERR DMRegister_MoFEM(dm_name);
37 DMType dm_name_sub = "DMSUB";
38 CHKERR DMRegister_MoFEM(dm_name_sub);
39
40 // craete dm instance
41 DM dm;
42 CHKERR DMCreate(PETSC_COMM_WORLD, &dm);
43 CHKERR DMSetType(dm, dm_name);
44
45 DM subdm0, subdm1;
46 CHKERR DMCreate(PETSC_COMM_WORLD, &subdm0);
47 CHKERR DMSetType(subdm0, dm_name_sub);
48 CHKERR DMCreate(PETSC_COMM_WORLD, &subdm1);
49 CHKERR DMSetType(subdm1, dm_name_sub);
50
51 // read mesh and create moab and mofem data structures
52 moab::Core mb_instance;
53 moab::Interface &moab = mb_instance;
54 ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
55 auto moab_comm_wrap =
56 boost::make_shared<WrapMPIComm>(PETSC_COMM_WORLD, false);
57 if (pcomm == NULL)
58 pcomm = new ParallelComm(&moab, moab_comm_wrap->get_comm());
59
60 const char *option;
61 option = "PARALLEL=BCAST_DELETE;PARALLEL_RESOLVE_SHARED_ENTS;PARTITION="
62 "PARALLEL_PARTITION;";
63 CHKERR moab.load_file(mesh_file_name, 0, option);
64
65 MoFEM::Core core(moab);
66 MoFEM::Interface &m_field = core;
67
68 EntityHandle root_set = moab.get_root_set();
69 // add all entities to database, all of them will be used
70 BitRefLevel bit_level0;
71 bit_level0.set(0);
72 CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
73 root_set, 3, bit_level0);
74 // define & build field
75 const int field_rank = 1;
76 CHKERR m_field.add_field("FIELD0", H1, AINSWORTH_LEGENDRE_BASE, field_rank);
77 // add entities to field
78 CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FIELD0");
79 // set app. order
80 CHKERR m_field.set_field_order(root_set, MBVERTEX, "FIELD0", 1);
81 CHKERR m_field.add_field("FIELD1", H1, AINSWORTH_LEGENDRE_BASE, field_rank);
82 // add entities to field
83 CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FIELD1");
84 // set app. order
85 CHKERR m_field.set_field_order(root_set, MBVERTEX, "FIELD1", 1);
86 CHKERR m_field.set_field_order(root_set, MBEDGE, "FIELD1", 2);
87
88 // build data structures for fields
89 CHKERR m_field.build_fields();
90
91 // define & build finite elements
92 CHKERR m_field.add_finite_element("FE00");
93 // Define rows/cols and element data
94 CHKERR m_field.modify_finite_element_add_field_row("FE00", "FIELD0");
95 CHKERR m_field.modify_finite_element_add_field_col("FE00", "FIELD0");
96 CHKERR m_field.modify_finite_element_add_field_data("FE00", "FIELD0");
97 // add entities to finite element
98 CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBTET, "FE00");
99
100 // define & build finite elements
101 CHKERR m_field.add_finite_element("FE11");
102 // Define rows/cols and element data
103 CHKERR m_field.modify_finite_element_add_field_row("FE11", "FIELD1");
104 CHKERR m_field.modify_finite_element_add_field_col("FE11", "FIELD1");
105 CHKERR m_field.modify_finite_element_add_field_data("FE11", "FIELD1");
106 // add entities to finite element
107 CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBTET, "FE11");
108
109 // define & build finite elements
110 CHKERR m_field.add_finite_element("FE01");
111 // Define rows/cols and element data
112 CHKERR m_field.modify_finite_element_add_field_row("FE01", "FIELD0");
113 CHKERR m_field.modify_finite_element_add_field_col("FE01", "FIELD1");
114 CHKERR m_field.modify_finite_element_add_field_data("FE01", "FIELD0");
115 CHKERR m_field.modify_finite_element_add_field_data("FE01", "FIELD1");
116 // add entities to finite element
117 CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBTET, "FE01");
118
119 // build finite elemnts
121 // build adjacencies
122 CHKERR m_field.build_adjacencies(bit_level0);
123
124 // set dm data structure which created mofem data structures
125 CHKERR DMMoFEMCreateMoFEM(dm, &m_field, "MAIN_PROBLEM", bit_level0);
126 CHKERR DMSetFromOptions(dm);
127 CHKERR DMMoFEMAddElement(dm, "FE00");
128 CHKERR DMMoFEMAddElement(dm, "FE11");
129 CHKERR DMMoFEMAddElement(dm, "FE01");
130 CHKERR DMSetUp(dm);
132 ->checkMPIAIJWithArraysMatrixFillIn<PetscGlobalIdx_mi_tag>(
133 "MAIN_PROBLEM", -1, -1, 1);
134
135 int nf;
136 char **field_names;
137 IS *fields;
138 CHKERR DMCreateFieldIS(dm, &nf, &field_names, &fields);
139 for (int f = 0; f != nf; f++) {
140 PetscPrintf(PETSC_COMM_WORLD, "%d field %s\n", f, field_names[f]);
141 CHKERR PetscFree(field_names[f]);
142 CHKERR ISDestroy(&(fields[f]));
143 }
144 CHKERR PetscFree(field_names);
145 CHKERR PetscFree(fields);
146
147 CHKERR DMMoFEMCreateSubDM(subdm0, dm, "SUB0");
148 CHKERR DMMoFEMSetSquareProblem(subdm0, PETSC_TRUE);
149 CHKERR DMMoFEMAddElement(subdm0, "FE11");
150
151 auto verts_ptr = boost::make_shared<Range>();
152 CHKERR moab.get_entities_by_type(0, MBVERTEX, *verts_ptr, true);
153
154 CHKERR DMMoFEMAddSubFieldRow(subdm0, "FIELD1", verts_ptr);
155 CHKERR DMMoFEMAddSubFieldCol(subdm0, "FIELD1", verts_ptr);
156 CHKERR DMSetUp(subdm0);
158 ->checkMPIAIJWithArraysMatrixFillIn<PetscGlobalIdx_mi_tag>("SUB0", -1,
159 -1, 1);
160 if (debug) {
161 Mat A;
162 CHKERR DMCreateMatrix(subdm0, &A);
163 MatView(A, PETSC_VIEWER_DRAW_WORLD);
164 std::string wait;
165 std::cin >> wait;
166 CHKERR MatDestroy(&A);
167 }
168
169 CHKERR DMMoFEMCreateSubDM(subdm1, dm, "SUB1");
170 CHKERR DMMoFEMSetSquareProblem(subdm1, PETSC_FALSE);
171 CHKERR DMMoFEMAddElement(subdm1, "FE01");
172 CHKERR DMMoFEMAddSubFieldRow(subdm1, "FIELD0");
173 CHKERR DMMoFEMAddSubFieldCol(subdm1, "FIELD1");
174 CHKERR DMSetUp(subdm1);
176 ->checkMPIAIJWithArraysMatrixFillIn<PetscGlobalIdx_mi_tag>("SUB1", -1,
177 -1, 1);
178
179 if (debug) {
180 Mat B;
181 CHKERR DMCreateMatrix(subdm1, &B);
182 MatView(B, PETSC_VIEWER_DRAW_WORLD);
183 std::string wait;
184 std::cin >> wait;
185 CHKERR MatDestroy(&B);
186 }
187
188 // destry dm
189 CHKERR DMDestroy(&dm);
190 CHKERR DMDestroy(&subdm0);
191 CHKERR DMDestroy(&subdm1);
192 }
194
195 // finish work cleaning memory, getting statistics, ect.
197
198 return 0;
199}
int main()
Definition: adol-c_atom.cpp:46
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
@ H1
continuous field
Definition: definitions.h:85
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:215
#define CHKERR
Inline error check.
Definition: definitions.h:535
static char help[]
static const bool debug
PetscErrorCode DMMoFEMCreateSubDM(DM subdm, DM dm, const char problem_name[])
Must be called by user to set Sub DM MoFEM data structures.
Definition: DMMoFEM.cpp:219
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Definition: DMMoFEM.cpp:483
PetscErrorCode DMMoFEMSetSquareProblem(DM dm, PetscBool square_problem)
set squared problem
Definition: DMMoFEM.cpp:442
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:118
PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[])
Definition: DMMoFEM.cpp:242
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:47
PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[])
Definition: DMMoFEM.cpp:275
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
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 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
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
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
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
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.
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.
char mesh_file_name[255]
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
constexpr AssemblyType A
Managing BitRefLevels.
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
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.
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.
Matrix manager is used to build and partition problems.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.