v0.13.2
Loading...
Searching...
No Matches
Functions | Variables
mofem_part.cpp File Reference
#include <MoFEM.hpp>

Go to the source code of this file.

Functions

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

Variables

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

Function Documentation

◆ main()

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

Definition at line 14 of file mofem_part.cpp.

14 {
15
16 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
17
18 try {
19
20 moab::Core mb_instance;
21 moab::Interface &moab = mb_instance;
22 ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
23 if (pcomm == NULL)
24 pcomm = new ParallelComm(&moab, PETSC_COMM_WORLD);
25
26 // global variables
27 char mesh_file_name[255];
28 char mesh_out_file[255] = "out.h5m";
29 PetscBool flg_file = PETSC_FALSE;
30 PetscBool flg_n_part = PETSC_FALSE;
31 PetscBool flg_part = PETSC_FALSE;
32 PetscBool only_tags = PETSC_FALSE;
33 PetscInt n_partas = 1;
34 PetscBool create_lower_dim_ents = PETSC_TRUE;
35 PetscInt dim = 3;
36 PetscInt adj_dim = 2;
37
38 ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "", "none", "none");
39 CHKERRQ(ierr);
40
41 CHKERR PetscOptionsString("-my_file", "mesh file name", "", "mesh.h5m",
42 mesh_file_name, 255, &flg_file);
43 if (flg_file != PETSC_TRUE)
44 CHKERR PetscOptionsString("-file_name", "mesh file name", "", "mesh.h5m",
45 mesh_file_name, 255, &flg_file);
46 if (flg_file != PETSC_TRUE)
47 SETERRQ(PETSC_COMM_SELF, 1,
48 "*** ERROR -my_file (-file_name) (mesh file needed)");
49
50 const char *option;
51 option = "";
52 CHKERR moab.load_file(mesh_file_name, 0, option);
53
54 CHKERR PetscOptionsString("-output_file", "output mesh file name", "",
55 "out.h5m", mesh_out_file, 255, PETSC_NULL);
56 CHKERR PetscOptionsInt("-my_nparts", "number of parts", "", n_partas,
57 &n_partas, &flg_n_part);
58 CHKERR PetscOptionsInt("-nparts", "number of parts", "", n_partas,
59 &n_partas, &flg_part);
60
61 if (!flg_n_part && !flg_part)
62 SETERRQ(PETSC_COMM_SELF, 1,
63 "*** ERROR partitioning number not given (-nparts)");
64
65 auto get_nb_ents_by_dim = [&](const int dim) {
66 int nb;
67 CHKERR moab.get_number_entities_by_dimension(0, dim, nb);
68 return nb;
69 };
70 for (; dim >= 0; dim--) {
71 if (get_nb_ents_by_dim(dim))
72 break;
73 }
74
75 if (!dim)
76 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
77 "Dimension of entities to partition not found");
78
79 CHKERR PetscOptionsInt("-dim", "entities dim", "", dim, &dim, PETSC_NULL);
80 adj_dim = dim - 1;
81 CHKERR PetscOptionsInt("-adj_dim", "adjacency dim", "", adj_dim, &adj_dim,
82 PETSC_NULL);
83 CHKERR PetscOptionsBool(
84 "-my_create_lower_dim_ents", "if tru create lower dimension entireties",
85 "", create_lower_dim_ents, &create_lower_dim_ents, PETSC_NULL);
86 CHKERR PetscOptionsBool("-block_tags", "only block and meshsests tags", "",
87 only_tags, &only_tags, PETSC_NULL);
88
89 ierr = PetscOptionsEnd();
90 CHKERRQ(ierr);
91
92 // Create MoFEM database
93 MoFEM::Core core(moab);
94 MoFEM::Interface &m_field = core;
95
96 auto meshsets_interface_ptr = m_field.getInterface<MeshsetsManager>();
97 CHKERR meshsets_interface_ptr->setMeshsetFromFile();
98
99 MOFEM_LOG_CHANNEL("WORLD");
100 MOFEM_LOG_TAG("WORLD", "mofem_part")
101 for (CubitMeshSet_multiIndex::iterator cit =
102 meshsets_interface_ptr->getBegin();
103 cit != meshsets_interface_ptr->getEnd(); cit++)
104 MOFEM_LOG("WORLD", Sev::inform) << *cit;
105 MOFEM_LOG_CHANNEL("WORLD");
106
107 {
108 Range ents_dim;
109 CHKERR moab.get_entities_by_dimension(0, dim, ents_dim, false);
110 if (create_lower_dim_ents) {
111 if (dim == 3) {
112 Range faces;
113 CHKERR moab.get_adjacencies(ents_dim, 2, true, faces,
114 moab::Interface::UNION);
115 }
116 if (dim > 2) {
117 Range edges;
118 CHKERR moab.get_adjacencies(ents_dim, 1, true, edges,
119 moab::Interface::UNION);
120 }
121 }
122 CommInterface *comm_interafce_ptr = m_field.getInterface<CommInterface>();
123 CHKERR comm_interafce_ptr->partitionMesh(
124 ents_dim, dim, adj_dim, n_partas, nullptr, nullptr, nullptr, VERBOSE);
125 }
126
127 auto get_tag_list = [&]() {
128 std::vector<Tag> tags_list;
129 auto meshsets_mng = m_field.getInterface<MeshsetsManager>();
130 auto &list = meshsets_mng->getMeshsetsMultindex();
131 for (auto &m : list) {
132 auto meshset = m.getMeshset();
133 std::vector<Tag> tmp_tags_list;
134 CHKERR m_field.get_moab().tag_get_tags_on_entity(meshset,
135 tmp_tags_list);
136 for (auto t : tmp_tags_list) {
137 tags_list.push_back(t);
138 }
139 }
140
141 return tags_list;
142 };
143
144 EntityHandle root_mesh = 0;
145 if (m_field.get_comm_rank() == 0) {
146 if (only_tags) {
147 auto tags_list = get_tag_list();
148 std::sort(tags_list.begin(), tags_list.end());
149 auto new_end = std::unique(tags_list.begin(), tags_list.end());
150 tags_list.resize(std::distance(tags_list.begin(), new_end));
151 tags_list.push_back(pcomm->part_tag());
152 // tags_list.push_back(m_field.get_moab().globalId_tag());
153 CHKERR moab.write_file(mesh_out_file, "MOAB", "", &root_mesh, 1,
154 &*tags_list.begin(), tags_list.size());
155 } else {
156 CHKERR moab.write_file(mesh_out_file, "MOAB", "");
157 }
158 }
159 }
161
163 CHKERRQ(ierr);
164
165 return 0;
166}
@ VERBOSE
Definition: definitions.h:209
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:215
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
#define CHKERR
Inline error check.
Definition: definitions.h:535
const int dim
FTensor::Index< 'm', SPACE_DIM > m
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
SeverityLevel
Severity levels.
Definition: LogManager.hpp:33
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
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.
char mesh_file_name[255]
static char help[]
Definition: mofem_part.cpp:12
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
multi_index_container< CubitMeshSets, indexed_by< hashed_unique< tag< Meshset_mi_tag >, member< CubitMeshSets, EntityHandle, &CubitMeshSets::meshset > >, ordered_non_unique< tag< CubitMeshsetType_mi_tag >, const_mem_fun< CubitMeshSets, unsigned long int, &CubitMeshSets::getBcTypeULong > >, ordered_non_unique< tag< CubitMeshsetMaskedType_mi_tag >, const_mem_fun< CubitMeshSets, unsigned long int, &CubitMeshSets::getMaksedBcTypeULong > >, ordered_non_unique< tag< CubitMeshsets_name >, const_mem_fun< CubitMeshSets, std::string, &CubitMeshSets::getName > >, hashed_unique< tag< Composite_Cubit_msId_And_MeshsetType_mi_tag >, composite_key< CubitMeshSets, const_mem_fun< CubitMeshSets, int, &CubitMeshSets::getMeshsetId >, const_mem_fun< CubitMeshSets, unsigned long int, &CubitMeshSets::getMaksedBcTypeULong > > > > > CubitMeshSet_multiIndex
Stores data about meshsets (see CubitMeshSets) storing data about boundary conditions,...
constexpr double t
plate stiffness
Definition: plate.cpp:59
Managing BitRefLevels.
virtual moab::Interface & get_moab()=0
virtual int get_comm_rank() const =0
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.
Interface for managing meshsets containing materials and boundary conditions.
CubitMeshSet_multiIndex & getMeshsetsMultindex()
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.

Variable Documentation

◆ help

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

Definition at line 12 of file mofem_part.cpp.