![]() |
OGS
|
2013/13/06 KR Initial implementation
Definition in file NodeReordering.cpp.
#include <tclap/CmdLine.h>#include <Eigen/Dense>#include <algorithm>#include <array>#include <concepts>#include <functional>#include <vector>#include "BaseLib/Algorithm.h"#include "BaseLib/Logging.h"#include "BaseLib/MPI.h"#include "BaseLib/TCLAPArguments.h"#include "InfoLib/GitInfo.h"#include "MeshLib/ElementCoordinatesMappingLocal.h"#include "MeshLib/Elements/Element.h"#include "MeshLib/IO/readMeshFromFile.h"#include "MeshLib/IO/writeMeshToFile.h"#include "MeshLib/Mesh.h"#include "MeshLib/MeshEnums.h"#include "MeshLib/Node.h"#include "MeshLib/Utils/GetSpaceDimension.h"#include "NumLib/Fem/ShapeFunction/ShapeHex20.h"#include "NumLib/Fem/ShapeFunction/ShapeHex8.h"#include "NumLib/Fem/ShapeFunction/ShapeLine2.h"#include "NumLib/Fem/ShapeFunction/ShapeLine3.h"#include "NumLib/Fem/ShapeFunction/ShapePoint1.h"#include "NumLib/Fem/ShapeFunction/ShapePrism15.h"#include "NumLib/Fem/ShapeFunction/ShapePrism6.h"#include "NumLib/Fem/ShapeFunction/ShapePyra13.h"#include "NumLib/Fem/ShapeFunction/ShapePyra5.h"#include "NumLib/Fem/ShapeFunction/ShapeQuad4.h"#include "NumLib/Fem/ShapeFunction/ShapeQuad8.h"#include "NumLib/Fem/ShapeFunction/ShapeQuad9.h"#include "NumLib/Fem/ShapeFunction/ShapeTet10.h"#include "NumLib/Fem/ShapeFunction/ShapeTet4.h"#include "NumLib/Fem/ShapeFunction/ShapeTri3.h"#include "NumLib/Fem/ShapeFunction/ShapeTri6.h"Go to the source code of this file.
Classes | |
| struct | ElementReorderConfigBase |
Functions | |
| template<typename GradShapeFunction, int Dim> | |
| bool | checkJacobianDeterminant (MeshLib::Element const &e, int const mesh_space_dimension, std::array< double, Dim > const &xi, bool const check_reordered=false) |
| template<ShapeFunction ShapeFunc, int Dim> | |
| ElementReorderConfigBase | makeElementConfig (std::array< double, Dim > xi, std::function< void(MeshLib::Element &, const std::vector< MeshLib::Node * > &)> reorder) |
| void | reverseNodeOrder (std::vector< MeshLib::Element * > &elements, int const mesh_space_dimension, bool const forced) |
| Reverses order of nodes. In particular, this fixes issues between (Gmsh or 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[]) |
Variables | |
| static const std::array< ElementReorderConfigBase, static_cast< int >(MeshLib::CellType::enum_length)> | element_configs_array |
| bool checkJacobianDeterminant | ( | MeshLib::Element const & | e, |
| int const | mesh_space_dimension, | ||
| std::array< double, Dim > const & | xi, | ||
| bool const | check_reordered = false ) |
Definition at line 52 of file NodeReordering.cpp.
References MathLib::Point3d::asEigenVector3d(), MeshLib::Element::getID(), MeshLib::ElementCoordinatesMappingLocal::getMappedCoordinates(), MeshLib::Element::getNumberOfNodes(), and OGS_FATAL.
Referenced by makeElementConfig().
| 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 369 of file NodeReordering.cpp.
References MeshLib::PRISM.
Referenced by main().
| int main | ( | int | argc, |
| char * | argv[] ) |
Definition at line 427 of file NodeReordering.cpp.
References MeshLib::enum_length, enum_length, fixVtkInconsistencies(), MeshLib::getSpaceDimension(), MeshLib::HEX20, MeshLib::HEX27, MeshLib::HEX8, INFO(), BaseLib::initOGSLogger(), MeshLib::INVALID, INVALID, MeshLib::LINE2, MeshLib::LINE3, BaseLib::makeLogLevelArg(), OGS_FATAL, GitInfoLib::GitInfo::ogs_version, MeshLib::POINT1, MeshLib::PRISM15, MeshLib::PRISM18, MeshLib::PRISM6, MeshLib::PYRAMID13, MeshLib::PYRAMID5, MeshLib::QUAD4, MeshLib::QUAD8, MeshLib::QUAD9, MeshLib::IO::readMeshFromFile(), reorderNonlinearNodes(), reverseNodeOrder(), MeshLib::TET10, MeshLib::TET4, MeshLib::TRI3, MeshLib::TRI6, and MeshLib::IO::writeMeshToFile().
| ElementReorderConfigBase makeElementConfig | ( | std::array< double, Dim > | xi, |
| std::function< void(MeshLib::Element &, const std::vector< MeshLib::Node * > &)> | reorder ) |
Definition at line 105 of file NodeReordering.cpp.
References checkJacobianDeterminant().
| void reorderNonlinearNodes | ( | MeshLib::Mesh & | mesh | ) |
Orders the base nodes of each elements before its non-linear nodes.
Definition at line 393 of file NodeReordering.cpp.
References MeshLib::Mesh::getElements(), MeshLib::Mesh::getNodes(), MeshLib::idsComparator(), BaseLib::makeVectorUnique(), and MeshLib::Mesh::resetNodeIDs().
Referenced by main().
| void reverseNodeOrder | ( | std::vector< MeshLib::Element * > & | elements, |
| int const | mesh_space_dimension, | ||
| bool const | forced ) |
Reverses order of nodes. In particular, this fixes issues between (Gmsh or OGS5) and OGS6 meshes.
| elements | Mesh elements whose nodes should be reordered |
| forced | If 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 292 of file NodeReordering.cpp.
References MeshLib::CellType2String(), element_configs_array, MeshLib::HEX27, INFO(), MeshLib::INVALID, OGS_FATAL, MeshLib::POINT1, and MeshLib::PRISM18.
Referenced by main().
|
static |
Definition at line 122 of file NodeReordering.cpp.
Referenced by reverseNodeOrder().