OGS
NodeReordering.cpp
Go to the documentation of this file.
1
12#include <tclap/CmdLine.h>
13
14#ifdef USE_PETSC
15#include <mpi.h>
16#endif
17
18#include <algorithm>
19#include <array>
20#include <memory>
21#include <vector>
22
23#include "BaseLib/Algorithm.h"
24#include "InfoLib/GitInfo.h"
28#include "MeshLib/Mesh.h"
29#include "MeshLib/Node.h"
30
40void reverseNodeOrder(std::vector<MeshLib::Element*>& elements,
41 bool const forced)
42{
43 std::size_t n_corrected_elements = 0;
44 std::size_t const nElements(elements.size());
45 for (std::size_t i = 0; i < nElements; ++i)
46 {
47 if (!forced && elements[i]->testElementNodeOrder())
48 {
49 continue;
50 }
51 n_corrected_elements++;
52
53 const unsigned nElemNodes(elements[i]->getNumberOfBaseNodes());
54 std::vector<MeshLib::Node*> nodes(elements[i]->getNodes(),
55 elements[i]->getNodes() + nElemNodes);
56
57 switch (elements[i]->getGeomType())
58 {
60 for (std::size_t j = 0; j < 4; ++j)
61 {
62 elements[i]->setNode(j, nodes[(j + 1) % 4]);
63 }
64 break;
66 elements[i]->setNode(0, nodes[1]);
67 elements[i]->setNode(1, nodes[0]);
68 elements[i]->setNode(2, nodes[3]);
69 elements[i]->setNode(3, nodes[2]);
70 break;
72 for (std::size_t j = 0; j < 3; ++j)
73 {
74 elements[i]->setNode(j, nodes[j + 3]);
75 elements[i]->setNode(j + 3, nodes[j]);
76 }
77 break;
79 for (std::size_t j = 0; j < 4; ++j)
80 {
81 elements[i]->setNode(j, nodes[j + 4]);
82 elements[i]->setNode(j + 4, nodes[j]);
83 }
84 break;
85 default:
86 for (std::size_t j = 0; j < nElemNodes; ++j)
87 {
88 elements[i]->setNode(j, nodes[nElemNodes - j - 1]);
89 }
90 }
91 }
92
93 INFO("Corrected {:d} elements.", n_corrected_elements);
94}
95
99void fixVtkInconsistencies(std::vector<MeshLib::Element*>& elements)
100{
101 std::size_t const nElements(elements.size());
102 for (std::size_t i = 0; i < nElements; ++i)
103 {
104 const unsigned nElemNodes(elements[i]->getNumberOfBaseNodes());
105 std::vector<MeshLib::Node*> nodes(elements[i]->getNodes(),
106 elements[i]->getNodes() + nElemNodes);
107
108 for (std::size_t j = 0; j < nElemNodes; ++j)
109 {
110 if (elements[i]->getGeomType() == MeshLib::MeshElemType::PRISM)
111 {
112 for (std::size_t k = 0; k < 3; ++k)
113 {
114 elements[i]->setNode(k, nodes[k + 3]);
115 elements[i]->setNode(k + 3, nodes[k]);
116 }
117 break;
118 }
119 }
120 }
121}
122
125{
126 std::vector<MeshLib::Node*> base_nodes;
127 std::vector<MeshLib::Node*> nonlinear_nodes;
128 for (MeshLib::Element const* e : mesh.getElements())
129 {
130 for (unsigned i = 0; i < e->getNumberOfBaseNodes(); i++)
131 {
132 base_nodes.push_back(const_cast<MeshLib::Node*>(e->getNode(i)));
133 }
134 for (unsigned i = e->getNumberOfBaseNodes(); i < e->getNumberOfNodes();
135 i++)
136 {
137 nonlinear_nodes.push_back(
138 const_cast<MeshLib::Node*>(e->getNode(i)));
139 }
140 }
141
142 BaseLib::makeVectorUnique(base_nodes,
144 BaseLib::makeVectorUnique(nonlinear_nodes,
146
147 std::vector<MeshLib::Node*>& allnodes =
148 const_cast<std::vector<MeshLib::Node*>&>(mesh.getNodes());
149 allnodes.clear();
150
151 allnodes.insert(allnodes.end(), base_nodes.begin(), base_nodes.end());
152 allnodes.insert(allnodes.end(), nonlinear_nodes.begin(),
153 nonlinear_nodes.end());
154
155 mesh.resetNodeIDs();
156}
157
158int main(int argc, char* argv[])
159{
160 TCLAP::CmdLine cmd(
161 "Reorders mesh nodes in elements to make old or incorrectly ordered "
162 "meshes compatible with OGS6.\n"
163 "Three options are available:\n"
164 "Method 0: Reversing order of nodes for all elements.\n"
165 "Method 1: Reversing order of nodes unless it's perceived correct by "
166 "OGS6 standards. This is the default selection.\n"
167 "Method 2: Fixing node ordering issues between VTK and OGS6 (only "
168 "applies to prism-elements)\n"
169 "Method 3: Re-ordering of mesh node vector such that all base nodes "
170 "are sorted before all nonlinear nodes.\n\n"
171 "OpenGeoSys-6 software, version " +
173 ".\n"
174 "Copyright (c) 2012-2024, OpenGeoSys Community "
175 "(http://www.opengeosys.org)",
177
178 std::vector<int> method_ids{0, 1, 2, 3};
179 TCLAP::ValuesConstraint<int> allowed_values(method_ids);
180 TCLAP::ValueArg<int> method_arg("m", "method",
181 "reordering method selection", false, 1,
182 &allowed_values);
183 cmd.add(method_arg);
184 TCLAP::ValueArg<std::string> output_mesh_arg(
185 "o", "output_mesh", "the name of the output mesh file", true, "",
186 "filename");
187 cmd.add(output_mesh_arg);
188 TCLAP::ValueArg<std::string> input_mesh_arg(
189 "i", "input_mesh", "the name of the input mesh file", true, "",
190 "filename");
191 cmd.add(input_mesh_arg);
192 cmd.parse(argc, argv);
193
194#ifdef USE_PETSC
195 MPI_Init(&argc, &argv);
196#endif
197
198 std::unique_ptr<MeshLib::Mesh> mesh(
199 MeshLib::IO::readMeshFromFile(input_mesh_arg.getValue()));
200
201 if (!mesh)
202 {
203#ifdef USE_PETSC
204 MPI_Finalize();
205#endif
206 return EXIT_FAILURE;
207 }
208
209 INFO("Reordering nodes... ");
210 if (!method_arg.isSet() || method_arg.getValue() < 2)
211 {
212 bool const forced = (method_arg.getValue() == 0);
214 const_cast<std::vector<MeshLib::Element*>&>(mesh->getElements()),
215 forced);
216 }
217 else if (method_arg.getValue() == 2)
218 {
220 const_cast<std::vector<MeshLib::Element*>&>(mesh->getElements()));
221 }
222 else if (method_arg.getValue() == 3)
223 {
225 }
226
227 MeshLib::IO::writeMeshToFile(*mesh, output_mesh_arg.getValue());
228
229 INFO("VTU file written.");
230
231#ifdef USE_PETSC
232 MPI_Finalize();
233#endif
234 return EXIT_SUCCESS;
235}
Definition of the Element class.
Git information.
std::vector< std::size_t > getNodes(GeoLib::Point const &pnt, std::vector< MeshLib::Node * > const &nodes, MeshLib::PropertyVector< int > const &mat_ids, std::pair< int, int > const &mat_limits, std::pair< double, double > const &elevation_limits, MeshLib::Mesh const &mesh)
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
Definition of the Mesh class.
int main(int argc, char *argv[])
void fixVtkInconsistencies(std::vector< MeshLib::Element * > &elements)
void reverseNodeOrder(std::vector< MeshLib::Element * > &elements, bool const forced)
Reverses order of nodes. In particular, this fixes issues between OGS5 and OGS6 meshes.
void reorderNonlinearNodes(MeshLib::Mesh &mesh)
Orders the base nodes of each elements before its non-linear nodes.
Definition of the Node class.
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition Mesh.h:106
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:109
void resetNodeIDs()
Resets the IDs of all mesh-nodes to their position in the node vector.
Definition Mesh.cpp:159
void makeVectorUnique(std::vector< T > &v)
Definition Algorithm.h:176
GITINFOLIB_EXPORT const std::string ogs_version
MeshLib::Mesh * readMeshFromFile(const std::string &file_name, bool const compute_element_neighbors)
int writeMeshToFile(const MeshLib::Mesh &mesh, std::filesystem::path const &file_path, std::set< std::string > variable_output_names)
bool idsComparator(T const a, T const b)
Definition Mesh.h:206
Definition of readMeshFromFile function.