v0.14.0
mofem_part.cpp
Go to the documentation of this file.
1 /** \file reading_med.cpp
2 
3  \brief Partition mesh and configuring blocksets
4 
5 */
6 
7 
8 
9 #include <MoFEM.hpp>
10 using namespace MoFEM;
11 
12 static char help[] = "-file_name <filename> : mesh file name\n"
13  "-output_file <filename> : output mesh file name\n"
14  "-nparts <int> : number of parts\n"
15  "-dim <int> : entities dim\n"
16  "-adj_dim <int> : adjacency dim\n"
17  "-my_create_lower_dim_ents <bool> : "
18  "if true create lower dimension entireties\n"
19  "-block_tags <bool> : only block and meshsests tags \n\n";
20 
21 int main(int argc, char *argv[]) {
22 
23  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
24 
25  try {
26 
27  moab::Core mb_instance;
28  moab::Interface &moab = mb_instance;
29  ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
30  if (pcomm == NULL)
31  pcomm = new ParallelComm(&moab, PETSC_COMM_WORLD);
32 
33  // global variables
34  char mesh_file_name[255];
35  char mesh_out_file[255] = "out.h5m";
36  PetscBool flg_file = PETSC_FALSE;
37  PetscBool flg_n_part = PETSC_FALSE;
38  PetscBool flg_part = PETSC_FALSE;
39  PetscBool only_tags = PETSC_FALSE;
40  PetscInt n_partas = 1;
41  PetscBool create_lower_dim_ents = PETSC_TRUE;
42  PetscInt dim = 3;
43  PetscInt adj_dim = 0;
44 
45  ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "", "none", "none");
46  CHKERRQ(ierr);
47 
48  CHKERR PetscOptionsString("-my_file", "mesh file name", "", "mesh.h5m",
49  mesh_file_name, 255, &flg_file);
50  if (flg_file != PETSC_TRUE)
51  CHKERR PetscOptionsString("-file_name", "mesh file name", "", "mesh.h5m",
52  mesh_file_name, 255, &flg_file);
53  if (flg_file != PETSC_TRUE)
54  SETERRQ(PETSC_COMM_SELF, 1,
55  "*** ERROR -my_file (-file_name) (mesh file needed)");
56 
57  const char *option;
58  option = "";
59  CHKERR moab.load_file(mesh_file_name, 0, option);
60 
61  CHKERR PetscOptionsString("-output_file", "output mesh file name", "",
62  "out.h5m", mesh_out_file, 255, PETSC_NULL);
63  CHKERR PetscOptionsInt("-my_nparts", "number of parts", "", n_partas,
64  &n_partas, &flg_n_part);
65  CHKERR PetscOptionsInt("-nparts", "number of parts", "", n_partas,
66  &n_partas, &flg_part);
67 
68  if (!flg_n_part && !flg_part)
69  SETERRQ(PETSC_COMM_SELF, 1,
70  "*** ERROR partitioning number not given (-nparts)");
71 
72  auto get_nb_ents_by_dim = [&](const int dim) {
73  int nb;
74  CHKERR moab.get_number_entities_by_dimension(0, dim, nb);
75  return nb;
76  };
77  for (; dim >= 0; dim--) {
78  if (get_nb_ents_by_dim(dim))
79  break;
80  }
81 
82  if (!dim)
83  SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
84  "Dimension of entities to partition not found");
85 
86  CHKERR PetscOptionsInt("-dim", "entities dim", "", dim, &dim, PETSC_NULL);
87  // adj_dim = dim - 1;
88  CHKERR PetscOptionsInt("-adj_dim", "adjacency dim", "", adj_dim, &adj_dim,
89  PETSC_NULL);
90  CHKERR PetscOptionsBool(
91  "-my_create_lower_dim_ents", "if true create lower dimension entireties",
92  "", create_lower_dim_ents, &create_lower_dim_ents, PETSC_NULL);
93  CHKERR PetscOptionsBool("-block_tags", "only block and meshsests tags", "",
94  only_tags, &only_tags, PETSC_NULL);
95 
96  ierr = PetscOptionsEnd();
97  CHKERRQ(ierr);
98 
99  auto get_dim = [](const Range &ents) -> int {
100  for (int d : {3, 2, 1})
101  if (ents.num_of_dimension(d))
102  return d;
103  return 0;
104  };
105 
106  // Create MoFEM database
107  MoFEM::Core core(moab);
108  MoFEM::Interface &m_field = core;
109 
110  auto meshsets_interface_ptr = m_field.getInterface<MeshsetsManager>();
111  CHKERR meshsets_interface_ptr->setMeshsetFromFile();
112 
113  MOFEM_LOG_CHANNEL("WORLD");
114  MOFEM_LOG_TAG("WORLD", "mofem_part");
115  int min_dim = 3;
116  for (auto cit = meshsets_interface_ptr->getBegin();
117  cit != meshsets_interface_ptr->getEnd(); cit++) {
118  Range ents;
119  CHKERR m_field.get_moab().get_entities_by_handle(cit->getMeshset(), ents,
120  true);
121  min_dim = std::min(min_dim, get_dim(ents));
122  MOFEM_LOG("WORLD", Sev::inform) << *cit;
123  }
124 
125  int g_dim; // global dimension
126  MPI_Allreduce(&min_dim, &g_dim, 1, MPI_INT, MPI_MIN, m_field.get_comm());
127  if (g_dim < adj_dim) {
128  MOFEM_LOG("WORLD", Sev::warning)
129  << "The minimum meshsets dimension is = " << min_dim;
130  }
131  if (adj_dim >= dim) {
132  MOFEM_LOG("WORLD", Sev::warning)
133  << "The -adj_dim >= dim, adj_dim = " << adj_dim << " dim = " << dim;
134  }
135  MOFEM_LOG_CHANNEL("WORLD");
136 
137  {
138  Range ents_dim;
139  CHKERR moab.get_entities_by_dimension(0, dim, ents_dim, false);
140  if (create_lower_dim_ents) {
141  if (dim == 3) {
142  Range faces;
143  CHKERR moab.get_adjacencies(ents_dim, 2, true, faces,
144  moab::Interface::UNION);
145  }
146  if (dim > 2) {
147  Range edges;
148  CHKERR moab.get_adjacencies(ents_dim, 1, true, edges,
149  moab::Interface::UNION);
150  }
151  }
152  CommInterface *comm_interafce_ptr = m_field.getInterface<CommInterface>();
153  CHKERR comm_interafce_ptr->partitionMesh(
154  ents_dim, dim, adj_dim, n_partas, nullptr, nullptr, nullptr, VERBOSE);
155  }
156 
157  auto get_tag_list = [&]() {
158  std::vector<Tag> tags_list;
159  auto meshsets_mng = m_field.getInterface<MeshsetsManager>();
160  auto &list = meshsets_mng->getMeshsetsMultindex();
161  for (auto &m : list) {
162  auto meshset = m.getMeshset();
163  std::vector<Tag> tmp_tags_list;
164  CHKERR m_field.get_moab().tag_get_tags_on_entity(meshset,
165  tmp_tags_list);
166  for (auto t : tmp_tags_list) {
167  tags_list.push_back(t);
168  }
169  }
170 
171  return tags_list;
172  };
173 
174  EntityHandle root_mesh = 0;
175  if (m_field.get_comm_rank() == 0) {
176  if (only_tags) {
177  auto tags_list = get_tag_list();
178  std::sort(tags_list.begin(), tags_list.end());
179  auto new_end = std::unique(tags_list.begin(), tags_list.end());
180  tags_list.resize(std::distance(tags_list.begin(), new_end));
181  tags_list.push_back(pcomm->part_tag());
182  // tags_list.push_back(m_field.get_moab().globalId_tag());
183  CHKERR moab.write_file(mesh_out_file, "MOAB", "", &root_mesh, 1,
184  &*tags_list.begin(), tags_list.size());
185  } else {
186  CHKERR moab.write_file(mesh_out_file, "MOAB", "");
187  }
188  }
189  }
190  CATCH_ERRORS;
191 
193  CHKERRQ(ierr);
194 
195  return 0;
196 }
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MYPCOMM_INDEX
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:228
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
MOFEM_LOG_CHANNEL
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
EntityHandle
MoFEM::CoreInterface::get_comm
virtual MPI_Comm & get_comm() const =0
help
static char help[]
Definition: mofem_part.cpp:12
MoFEM::CoreInterface::get_comm_rank
virtual int get_comm_rank() const =0
MoFEM.hpp
MoFEM::CommInterface::partitionMesh
MoFEMErrorCode partitionMesh(const Range &ents, const int dim, const int adj_dim, const int n_parts, Tag *th_vertex_weights=nullptr, Tag *th_edge_weights=nullptr, Tag *th_part_weights=nullptr, int verb=VERBOSE, const bool debug=false)
Set partition tag to each finite element in the problem.
Definition: CommInterface.cpp:868
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
VERBOSE
@ VERBOSE
Definition: definitions.h:222
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::BcManagerImplTools::get_dim
auto get_dim(const Range &ents)
Definition: BcManager.cpp:10
main
int main(int argc, char *argv[])
Definition: mofem_part.cpp:21
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
t
constexpr double t
plate stiffness
Definition: plate.cpp:58
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:22
Range
MoFEM::MeshsetsManager::getMeshsetsMultindex
CubitMeshSet_multiIndex & getMeshsetsMultindex()
Definition: MeshsetsManager.hpp:229
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
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
MoFEM::CommInterface
Managing BitRefLevels.
Definition: CommInterface.hpp:21
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::MeshsetsManager
Interface for managing meshsets containing materials and boundary conditions.
Definition: MeshsetsManager.hpp:104
sdf_hertz_2d_axisymm_plane.d
float d
Definition: sdf_hertz_2d_axisymm_plane.py:4
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80