v0.14.0
build_composite_problem.cpp
Go to the documentation of this file.
1 /** \file build_composite_problems.cpp
2  \example build_composite_problems.cpp
3  \brief Atom test for building composite problems
4 
5 */
6 
7 
8 
9 #include <MoFEM.hpp>
10 using namespace MoFEM;
11 
12 static char help[] = "...\n\n";
13 
14 int main(int argc, char *argv[]) {
15 
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 
33  // read mesh and create moab and mofem data structures
34 
35  moab::Core mb_instance;
36  moab::Interface &moab = mb_instance;
37  const char *option;
38  option = "";
39  CHKERR moab.load_file(mesh_file_name, 0, option);
40 
41 
42  MoFEM::Core core(moab);
43  MoFEM::Interface &m_field = core;
44 
45  EntityHandle root_set = moab.get_root_set();
46  Range tets;
47  CHKERR moab.get_entities_by_type(root_set, MBTET, tets, false);
48 
49  ProblemsManager *prb_mng_ptr = m_field.getInterface<ProblemsManager>();
50  CommInterface *comm_interface_ptr = m_field.getInterface<CommInterface>();
51  CHKERR comm_interface_ptr->partitionMesh(
52  tets, 3, 2, m_field.get_comm_size(), NULL, NULL, NULL);
53 
54  EntityHandle part_set;
55  CHKERR moab.create_meshset(MESHSET_SET, part_set);
56 
57  ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
58  if (pcomm == NULL)
59  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
60  "Communicator should be allocated");
61 
62  Tag part_tag = pcomm->part_tag();
63  Range proc_ents;
64  Range tagged_sets;
65  CHKERR m_field.get_moab().get_entities_by_type_and_tag(
66  0, MBENTITYSET, &part_tag, NULL, 1, tagged_sets,
67  moab::Interface::UNION);
68  for (Range::iterator mit = tagged_sets.begin(); mit != tagged_sets.end();
69  mit++) {
70  int part;
71  CHKERR moab.tag_get_data(part_tag, &*mit, 1, &part);
72  if (part == m_field.get_comm_rank()) {
73  // pcomm->partition_sets().insert(*mit);
74  CHKERR moab.get_entities_by_type(*mit, MBTET, proc_ents, true);
75  CHKERR moab.add_entities(part_set, proc_ents);
76  }
77  }
78 
79  Skinner skin(&m_field.get_moab());
80  Range tets_skin;
81  CHKERR skin.find_skin(0, tets, false, tets_skin);
82  Range proc_ents_skin[4];
83  proc_ents_skin[3] = proc_ents;
84  CHKERR skin.find_skin(0, proc_ents, false, proc_ents_skin[2]);
85  proc_ents_skin[2] = subtract(proc_ents_skin[2], tets_skin);
86  CHKERR moab.get_adjacencies(proc_ents_skin[2], 1, false, proc_ents_skin[1],
87  moab::Interface::UNION);
88  CHKERR moab.get_connectivity(proc_ents_skin[1], proc_ents_skin[0], true);
89  for (int dd = 0; dd != 3; dd++) {
90  CHKERR moab.add_entities(part_set, proc_ents_skin[dd]);
91  }
92 
93  if (0) {
94  std::ostringstream file_skin;
95  file_skin << "out_skin_" << m_field.get_comm_rank() << ".vtk";
96  EntityHandle meshset_skin;
97  CHKERR moab.create_meshset(MESHSET_SET, meshset_skin);
98  CHKERR moab.add_entities(meshset_skin, proc_ents_skin[2]);
99  CHKERR moab.add_entities(meshset_skin, proc_ents_skin[1]);
100  CHKERR moab.add_entities(meshset_skin, proc_ents_skin[0]);
101  CHKERR moab.write_file(file_skin.str().c_str(), "VTK", "", &meshset_skin,
102  1);
103  }
104 
105  CHKERR pcomm->resolve_shared_ents(0, proc_ents, 3, -1, proc_ents_skin);
106  Range owned_tets = proc_ents;
107 
108  if (0) {
109  std::ostringstream file_owned;
110  file_owned << "out_owned_" << m_field.get_comm_rank() << ".vtk";
111  EntityHandle meshset_owned;
112  CHKERR moab.create_meshset(MESHSET_SET, meshset_owned);
113  CHKERR moab.add_entities(meshset_owned, owned_tets);
114  CHKERR moab.write_file(file_owned.str().c_str(), "VTK", "",
115  &meshset_owned, 1);
116  }
117  Range shared_ents;
118  // Get entities shared with all other processors
119  CHKERR pcomm->get_shared_entities(-1, shared_ents);
120 
121  if (0) {
122  std::ostringstream file_shared_owned;
123  file_shared_owned << "out_shared_owned_" << m_field.get_comm_rank()
124  << ".vtk";
125  EntityHandle meshset_shared_owned;
126  CHKERR moab.create_meshset(MESHSET_SET, meshset_shared_owned);
127  CHKERR moab.add_entities(meshset_shared_owned, shared_ents);
128  CHKERR moab.write_file(file_shared_owned.str().c_str(), "VTK", "",
129  &meshset_shared_owned, 1);
130  }
131 
132  // set entities bit level
133  BitRefLevel bit_level0;
134  bit_level0.set(0);
135  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
136  part_set, 3, bit_level0);
137 
138  // Fields
139  CHKERR m_field.add_field("F1", H1, AINSWORTH_LEGENDRE_BASE, 1);
140  CHKERR m_field.add_field("F2", H1, AINSWORTH_LEGENDRE_BASE, 1);
141  // add entities to field
142  CHKERR m_field.add_ents_to_field_by_type(part_set, MBTET, "F1");
143  CHKERR m_field.add_ents_to_field_by_type(part_set, MBTET, "F2");
144  int order = 4;
145  CHKERR m_field.set_field_order(part_set, MBTET, "F1", order);
146  CHKERR m_field.set_field_order(part_set, MBTRI, "F1", order);
147  CHKERR m_field.set_field_order(part_set, MBEDGE, "F1", order);
148  CHKERR m_field.set_field_order(part_set, MBVERTEX, "F1", 1);
149  CHKERR m_field.set_field_order(part_set, MBTET, "F2", order);
150  CHKERR m_field.set_field_order(part_set, MBTRI, "F2", order);
151  CHKERR m_field.set_field_order(part_set, MBEDGE, "F2", order);
152  CHKERR m_field.set_field_order(part_set, MBVERTEX, "F2", 1);
153  CHKERR m_field.build_fields();
154 
155  // Elements
156  CHKERR m_field.add_finite_element("E1");
157  CHKERR m_field.add_finite_element("E2");
158  CHKERR m_field.modify_finite_element_add_field_row("E1", "F1");
159  CHKERR m_field.modify_finite_element_add_field_col("E1", "F1");
160  CHKERR m_field.modify_finite_element_add_field_row("E2", "F2");
161  CHKERR m_field.modify_finite_element_add_field_col("E2", "F2");
162  CHKERR m_field.add_ents_to_finite_element_by_type(part_set, MBTET, "E1");
163  CHKERR m_field.add_ents_to_finite_element_by_type(part_set, MBTET, "E2");
164  CHKERR m_field.build_finite_elements();
165  CHKERR m_field.build_adjacencies(bit_level0);
166 
167  // Problems
168  CHKERR m_field.add_problem("P1");
169  CHKERR m_field.add_problem("P2");
170  // set refinement level for problem
171  CHKERR m_field.modify_problem_ref_level_add_bit("P1", bit_level0);
172  CHKERR m_field.modify_problem_ref_level_add_bit("P2", bit_level0);
173  CHKERR m_field.modify_problem_add_finite_element("P1", "E1");
174  CHKERR m_field.modify_problem_add_finite_element("P2", "E2");
175 
176  // Build problems
177  CHKERR prb_mng_ptr->buildProblemOnDistributedMesh("P1", true, 1);
178  CHKERR prb_mng_ptr->buildProblemOnDistributedMesh("P2", true, 1);
179  CHKERR prb_mng_ptr->partitionFiniteElements("P1", true, 0,
180  m_field.get_comm_size(), 1);
181  CHKERR prb_mng_ptr->partitionGhostDofs("P1", 1);
182  CHKERR prb_mng_ptr->partitionFiniteElements("P2", true, 0,
183  m_field.get_comm_size(), 1);
184  CHKERR prb_mng_ptr->partitionGhostDofs("P2", 1);
185 
186  if (0) {
189  ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>("P1", m);
190  MatView(m, PETSC_VIEWER_DRAW_WORLD);
191  std::string wait;
192  std::cin >> wait;
193  }
194 
196  ->checkMPIAIJWithArraysMatrixFillIn<PetscGlobalIdx_mi_tag>("P1", -1, -1,
197  1);
199  ->checkMPIAIJWithArraysMatrixFillIn<PetscGlobalIdx_mi_tag>("P2", -1, -1,
200  1);
201 
202  // register new dm type, i.e. mofem
203  DMType dm_name = "MOFEM";
204  CHKERR DMRegister_MoFEM(dm_name);
205  // create dm instance
206  auto dm = createDM(m_field.get_comm(), dm_name);
207 
208  CHKERR DMMoFEMCreateMoFEM(dm, &m_field, "COMP", bit_level0);
209  CHKERR DMSetFromOptions(dm);
210  CHKERR DMMoFEMSetIsPartitioned(dm, PETSC_TRUE);
211  CHKERR DMMoFEMAddElement(dm, "E1");
212  CHKERR DMMoFEMAddElement(dm, "E2");
215  CHKERR DMSetUp(dm);
216 
217  if (0) {
220  ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>("COMP", m);
221  MatView(m, PETSC_VIEWER_DRAW_WORLD);
222  std::string wait;
223  std::cin >> wait;
224  }
225 
227  ->checkMPIAIJWithArraysMatrixFillIn<PetscGlobalIdx_mi_tag>("COMP", -1,
228  -1, 1);
229 
230  }
231  CATCH_ERRORS;
232 
233  // Finish work cleaning memory, getting statistics, etc.
235 
236  return 0;
237 }
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
MYPCOMM_INDEX
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:215
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
EntityHandle
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
help
static char help[]
Definition: build_composite_problem.cpp:12
MoFEM::CoreInterface::get_comm
virtual MPI_Comm & get_comm() const =0
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::CoreInterface::get_comm_rank
virtual int get_comm_rank() const =0
MoFEM.hpp
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:501
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:1975
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
MoFEM::DMMoFEMAddRowCompositeProblem
PetscErrorCode DMMoFEMAddRowCompositeProblem(DM dm, const char prb_name[])
Add problem to composite DM on row.
Definition: DMMoFEM.cpp:371
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
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::get_moab
virtual moab::Interface & get_moab()=0
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
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::createDM
auto createDM(MPI_Comm comm, const std::string dm_type_name)
Creates smart DM object.
Definition: PetscSmartObj.hpp:137
MoFEM::CoreInterface::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
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.
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:47
MoFEM::CoreInterface::get_comm_size
virtual int get_comm_size() const =0
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:118
MoFEM::ProblemsManager::partitionFiniteElements
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
Definition: ProblemsManager.cpp:2167
MoFEM::CoreInterface::modify_problem_ref_level_add_bit
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:23
MoFEM::MatrixManager
Matrix manager is used to build and partition problems.
Definition: MatrixManager.hpp:21
Range
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
FTensor::dd
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
MoFEM::CommInterface
Managing BitRefLevels.
Definition: CommInterface.hpp:21
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1094
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
main
int main(int argc, char *argv[])
Definition: build_composite_problem.cpp:14
MoFEM::ProblemsManager::buildProblemOnDistributedMesh
MoFEMErrorCode buildProblemOnDistributedMesh(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures, assuming that mesh is distributed (collective)
Definition: ProblemsManager.cpp:290
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::CoreInterface::modify_problem_add_finite_element
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
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::SmartPetscObj< Mat >
MoFEM::ProblemsManager::partitionGhostDofs
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
Definition: ProblemsManager.cpp:2339
MoFEM::CoreInterface::add_problem
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
MoFEM::DMMoFEMSetIsPartitioned
PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned)
Definition: DMMoFEM.cpp:1127