12#include <tclap/CmdLine.h>
40 std::size_t n_corrected_elements = 0;
41 std::size_t
const nElements(elements.size());
42 for (std::size_t i = 0; i < nElements; ++i)
44 if (!forced && elements[i]->testElementNodeOrder())
48 n_corrected_elements++;
50 const unsigned nElemNodes(elements[i]->getNumberOfBaseNodes());
51 std::vector<MeshLib::Node*> nodes(elements[i]->
getNodes(),
52 elements[i]->
getNodes() + nElemNodes);
54 switch (elements[i]->getGeomType())
57 for (std::size_t j = 0; j < 4; ++j)
59 elements[i]->setNode(j, nodes[(j + 1) % 4]);
63 elements[i]->setNode(0, nodes[1]);
64 elements[i]->setNode(1, nodes[0]);
65 elements[i]->setNode(2, nodes[3]);
66 elements[i]->setNode(3, nodes[2]);
69 for (std::size_t j = 0; j < 3; ++j)
71 elements[i]->setNode(j, nodes[j + 3]);
72 elements[i]->setNode(j + 3, nodes[j]);
76 for (std::size_t j = 0; j < 4; ++j)
78 elements[i]->setNode(j, nodes[j + 4]);
79 elements[i]->setNode(j + 4, nodes[j]);
83 for (std::size_t j = 0; j < nElemNodes; ++j)
85 elements[i]->setNode(j, nodes[nElemNodes - j - 1]);
90 INFO(
"Corrected {:d} elements.", n_corrected_elements);
98 std::size_t
const nElements(elements.size());
99 for (std::size_t i = 0; i < nElements; ++i)
101 const unsigned nElemNodes(elements[i]->getNumberOfBaseNodes());
102 std::vector<MeshLib::Node*> nodes(elements[i]->
getNodes(),
103 elements[i]->
getNodes() + nElemNodes);
105 for (std::size_t j = 0; j < nElemNodes; ++j)
109 for (std::size_t k = 0; k < 3; ++k)
111 elements[i]->setNode(k, nodes[k + 3]);
112 elements[i]->setNode(k + 3, nodes[k]);
123 std::vector<MeshLib::Node*> base_nodes;
124 std::vector<MeshLib::Node*> nonlinear_nodes;
127 for (
unsigned i = 0; i < e->getNumberOfBaseNodes(); i++)
129 base_nodes.push_back(
const_cast<MeshLib::Node*
>(e->getNode(i)));
131 for (
unsigned i = e->getNumberOfBaseNodes(); i < e->getNumberOfNodes();
134 nonlinear_nodes.push_back(
144 std::vector<MeshLib::Node*>& allnodes =
145 const_cast<std::vector<MeshLib::Node*>&
>(mesh.
getNodes());
148 allnodes.insert(allnodes.end(), base_nodes.begin(), base_nodes.end());
149 allnodes.insert(allnodes.end(), nonlinear_nodes.begin(),
150 nonlinear_nodes.end());
155int main(
int argc,
char* argv[])
158 "Reorders mesh nodes in elements to make old or incorrectly ordered "
159 "meshes compatible with OGS6.\n"
160 "Three options are available:\n"
161 "Method 0: Reversing order of nodes for all elements.\n"
162 "Method 1: Reversing order of nodes unless it's perceived correct by "
163 "OGS6 standards. This is the default selection.\n"
164 "Method 2: Fixing node ordering issues between VTK and OGS6 (only "
165 "applies to prism-elements)\n"
166 "Method 3: Re-ordering of mesh node vector such that all base nodes "
167 "are sorted before all nonlinear nodes.\n\n"
168 "OpenGeoSys-6 software, version " +
171 "Copyright (c) 2012-2024, OpenGeoSys Community "
172 "(http://www.opengeosys.org)",
175 std::vector<int> method_ids{0, 1, 2, 3};
176 TCLAP::ValuesConstraint<int> allowed_values(method_ids);
177 TCLAP::ValueArg<int> method_arg(
"m",
"method",
178 "reordering method selection",
false, 1,
181 TCLAP::ValueArg<std::string> output_mesh_arg(
182 "o",
"output_mesh",
"the name of the output mesh file",
true,
"",
184 cmd.add(output_mesh_arg);
185 TCLAP::ValueArg<std::string> input_mesh_arg(
186 "i",
"input_mesh",
"the name of the input mesh file",
true,
"",
188 cmd.add(input_mesh_arg);
189 cmd.parse(argc, argv);
193 std::unique_ptr<MeshLib::Mesh> mesh(
201 INFO(
"Reordering nodes... ");
202 if (!method_arg.isSet() || method_arg.getValue() < 2)
204 bool const forced = (method_arg.getValue() == 0);
206 const_cast<std::vector<MeshLib::Element*>&
>(mesh->getElements()),
209 else if (method_arg.getValue() == 2)
212 const_cast<std::vector<MeshLib::Element*>&
>(mesh->getElements()));
214 else if (method_arg.getValue() == 3)
221 INFO(
"VTU file written.");
Definition of the Element class.
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 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.
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
void resetNodeIDs()
Resets the IDs of all mesh-nodes to their position in the node vector.
void makeVectorUnique(std::vector< T > &v)
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 of readMeshFromFile function.