OGS
NodeReordering.cpp File Reference

Detailed Description

2013/13/06 KR Initial implementation

Definition in file NodeReordering.cpp.

#include <tclap/CmdLine.h>
#include <algorithm>
#include <array>
#include <memory>
#include <vector>
#include "BaseLib/Algorithm.h"
#include "BaseLib/MPI.h"
#include "InfoLib/GitInfo.h"
#include "MeshLib/Elements/Element.h"
#include "MeshLib/IO/readMeshFromFile.h"
#include "MeshLib/IO/writeMeshToFile.h"
#include "MeshLib/Mesh.h"
#include "MeshLib/Node.h"
Include dependency graph for NodeReordering.cpp:

Go to the source code of this file.

Functions

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 fixVtkInconsistencies (std::vector< MeshLib::Element * > &elements)
 
void reorderNonlinearNodes (MeshLib::Mesh &mesh)
 Orders the base nodes of each elements before its non-linear nodes.
 
int main (int argc, char *argv[])
 

Function Documentation

◆ fixVtkInconsistencies()

void fixVtkInconsistencies ( std::vector< MeshLib::Element * > & elements)

Fixes inconsistencies between VTK's and OGS' node order for prism elements. In particular, this fixes issues between OGS6 meshes with and without InSitu-Lib

Definition at line 96 of file NodeReordering.cpp.

97{
98 std::size_t const nElements(elements.size());
99 for (std::size_t i = 0; i < nElements; ++i)
100 {
101 const unsigned nElemNodes(elements[i]->getNumberOfBaseNodes());
102 std::vector<MeshLib::Node*> nodes(elements[i]->getNodes(),
103 elements[i]->getNodes() + nElemNodes);
104
105 for (std::size_t j = 0; j < nElemNodes; ++j)
106 {
107 if (elements[i]->getGeomType() == MeshLib::MeshElemType::PRISM)
108 {
109 for (std::size_t k = 0; k < 3; ++k)
110 {
111 elements[i]->setNode(k, nodes[k + 3]);
112 elements[i]->setNode(k + 3, nodes[k]);
113 }
114 break;
115 }
116 }
117 }
118}
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)

References getNodes(), and MeshLib::PRISM.

Referenced by main().

◆ main()

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

Definition at line 155 of file NodeReordering.cpp.

156{
157 TCLAP::CmdLine cmd(
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 " +
170 ".\n"
171 "Copyright (c) 2012-2024, OpenGeoSys Community "
172 "(http://www.opengeosys.org)",
174
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,
179 &allowed_values);
180 cmd.add(method_arg);
181 TCLAP::ValueArg<std::string> output_mesh_arg(
182 "o", "output_mesh", "the name of the output mesh file", true, "",
183 "filename");
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, "",
187 "filename");
188 cmd.add(input_mesh_arg);
189 cmd.parse(argc, argv);
190
191 BaseLib::MPI::Setup mpi_setup(argc, argv);
192
193 std::unique_ptr<MeshLib::Mesh> mesh(
194 MeshLib::IO::readMeshFromFile(input_mesh_arg.getValue()));
195
196 if (!mesh)
197 {
198 return EXIT_FAILURE;
199 }
200
201 INFO("Reordering nodes... ");
202 if (!method_arg.isSet() || method_arg.getValue() < 2)
203 {
204 bool const forced = (method_arg.getValue() == 0);
206 const_cast<std::vector<MeshLib::Element*>&>(mesh->getElements()),
207 forced);
208 }
209 else if (method_arg.getValue() == 2)
210 {
212 const_cast<std::vector<MeshLib::Element*>&>(mesh->getElements()));
213 }
214 else if (method_arg.getValue() == 3)
215 {
217 }
218
219 MeshLib::IO::writeMeshToFile(*mesh, output_mesh_arg.getValue());
220
221 INFO("VTU file written.");
222
223 return EXIT_SUCCESS;
224}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
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.
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)

References fixVtkInconsistencies(), INFO(), GitInfoLib::GitInfo::ogs_version, MeshLib::IO::readMeshFromFile(), reorderNonlinearNodes(), reverseNodeOrder(), and MeshLib::IO::writeMeshToFile().

◆ reorderNonlinearNodes()

void reorderNonlinearNodes ( MeshLib::Mesh & mesh)

Orders the base nodes of each elements before its non-linear nodes.

Definition at line 121 of file NodeReordering.cpp.

122{
123 std::vector<MeshLib::Node*> base_nodes;
124 std::vector<MeshLib::Node*> nonlinear_nodes;
125 for (MeshLib::Element const* e : mesh.getElements())
126 {
127 for (unsigned i = 0; i < e->getNumberOfBaseNodes(); i++)
128 {
129 base_nodes.push_back(const_cast<MeshLib::Node*>(e->getNode(i)));
130 }
131 for (unsigned i = e->getNumberOfBaseNodes(); i < e->getNumberOfNodes();
132 i++)
133 {
134 nonlinear_nodes.push_back(
135 const_cast<MeshLib::Node*>(e->getNode(i)));
136 }
137 }
138
139 BaseLib::makeVectorUnique(base_nodes,
141 BaseLib::makeVectorUnique(nonlinear_nodes,
143
144 std::vector<MeshLib::Node*>& allnodes =
145 const_cast<std::vector<MeshLib::Node*>&>(mesh.getNodes());
146 allnodes.clear();
147
148 allnodes.insert(allnodes.end(), base_nodes.begin(), base_nodes.end());
149 allnodes.insert(allnodes.end(), nonlinear_nodes.begin(),
150 nonlinear_nodes.end());
151
152 mesh.resetNodeIDs();
153}
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition Mesh.h:106
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:180
bool idsComparator(T const a, T const b)
Definition Mesh.h:206

References MeshLib::Mesh::getElements(), MeshLib::Mesh::getNodes(), MeshLib::idsComparator(), BaseLib::makeVectorUnique(), and MeshLib::Mesh::resetNodeIDs().

Referenced by main().

◆ reverseNodeOrder()

void reverseNodeOrder ( std::vector< MeshLib::Element * > & elements,
bool const forced )

Reverses order of nodes. In particular, this fixes issues between OGS5 and OGS6 meshes.

Parameters
elementsMesh elements whose nodes should be reordered
forcedIf true, nodes are reordered for all elements, if false it is first checked if the node order is correct according to OGS6 element definitions.

Definition at line 37 of file NodeReordering.cpp.

39{
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)
43 {
44 if (!forced && elements[i]->testElementNodeOrder())
45 {
46 continue;
47 }
48 n_corrected_elements++;
49
50 const unsigned nElemNodes(elements[i]->getNumberOfBaseNodes());
51 std::vector<MeshLib::Node*> nodes(elements[i]->getNodes(),
52 elements[i]->getNodes() + nElemNodes);
53
54 switch (elements[i]->getGeomType())
55 {
57 for (std::size_t j = 0; j < 4; ++j)
58 {
59 elements[i]->setNode(j, nodes[(j + 1) % 4]);
60 }
61 break;
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]);
67 break;
69 for (std::size_t j = 0; j < 3; ++j)
70 {
71 elements[i]->setNode(j, nodes[j + 3]);
72 elements[i]->setNode(j + 3, nodes[j]);
73 }
74 break;
76 for (std::size_t j = 0; j < 4; ++j)
77 {
78 elements[i]->setNode(j, nodes[j + 4]);
79 elements[i]->setNode(j + 4, nodes[j]);
80 }
81 break;
82 default:
83 for (std::size_t j = 0; j < nElemNodes; ++j)
84 {
85 elements[i]->setNode(j, nodes[nElemNodes - j - 1]);
86 }
87 }
88 }
89
90 INFO("Corrected {:d} elements.", n_corrected_elements);
91}

References getNodes(), MeshLib::HEXAHEDRON, INFO(), MeshLib::PRISM, MeshLib::PYRAMID, and MeshLib::TETRAHEDRON.

Referenced by main().