From 8224f3cd2d4487dd51b1aebb7eb1a3bd5ef030eb Mon Sep 17 00:00:00 2001 From: Matthias Hochsteger Date: Tue, 4 Oct 2022 12:26:02 +0200 Subject: [PATCH] New debug parameter to write mesh on error, python export --- libsrc/meshing/meshfunc.cpp | 14 +++++++++++++- libsrc/meshing/meshtype.cpp | 1 + libsrc/meshing/meshtype.hpp | 2 ++ libsrc/meshing/python_mesh.cpp | 18 ++++++++++++++++++ 4 files changed, 34 insertions(+), 1 deletion(-) diff --git a/libsrc/meshing/meshfunc.cpp b/libsrc/meshing/meshfunc.cpp index 3e4aab21..a59f1205 100644 --- a/libsrc/meshing/meshfunc.cpp +++ b/libsrc/meshing/meshfunc.cpp @@ -359,6 +359,8 @@ namespace netgen if (mesh.HasOpenQuads()) { + if(debugparam.write_mesh_on_error) + md.mesh->Save("open_quads_"+ToString(md.domain)+".vol.gz"); PrintSysError ("mesh has still open quads"); throw NgException ("Stop meshing since too many attempts"); // return MESHING3_GIVEUP; @@ -423,7 +425,11 @@ namespace netgen if (cntsteps > mp.maxoutersteps) - throw NgException ("Stop meshing since too many attempts"); + { + if(debugparam.write_mesh_on_error) + md.mesh->Save("meshing_error_domain_"+ToString(md.domain)+".vol.gz"); + throw NgException ("Stop meshing since too many attempts in domain " + ToString(md.domain)); + } PrintMessage (1, "start tetmeshing"); @@ -504,6 +510,8 @@ namespace netgen bool res = (mesh.CheckConsistentBoundary() != 0); if (res) { + if(debugparam.write_mesh_on_error) + md.mesh->Save("inconsistent_surface_domain_"+ToString(md.domain)+".vol.gz"); PrintError ("Surface mesh not consistent"); throw NgException ("Stop meshing since surface mesh not consistent"); } @@ -592,7 +600,11 @@ namespace netgen { if (mp.checkoverlappingboundary) if (md[i].mesh->CheckOverlappingBoundary()) + { + if(debugparam.write_mesh_on_error) + md[i].mesh->Save("overlapping_mesh_domain_"+ToString(md[i].domain)+".vol.gz"); throw NgException ("Stop meshing since boundary mesh is overlapping"); + } if(md[i].mesh->GetGeometry()->GetGeomType() == Mesh::GEOM_OCC) FillCloseSurface( md[i] ); diff --git a/libsrc/meshing/meshtype.cpp b/libsrc/meshing/meshtype.cpp index d49d99cd..068c1ed8 100644 --- a/libsrc/meshing/meshtype.cpp +++ b/libsrc/meshing/meshtype.cpp @@ -2898,6 +2898,7 @@ namespace netgen haltsegment = 0; haltsegmentp1 = 0; haltsegmentp2 = 0; + write_mesh_on_error = false; }; diff --git a/libsrc/meshing/meshtype.hpp b/libsrc/meshing/meshtype.hpp index 6722bde0..28d1a71e 100644 --- a/libsrc/meshing/meshtype.hpp +++ b/libsrc/meshing/meshtype.hpp @@ -1415,6 +1415,8 @@ namespace netgen /// int haltfacenr; /// + bool write_mesh_on_error; + /// DebugParameters (); }; diff --git a/libsrc/meshing/python_mesh.cpp b/libsrc/meshing/python_mesh.cpp index c79013a3..dd7b7910 100644 --- a/libsrc/meshing/python_mesh.cpp +++ b/libsrc/meshing/python_mesh.cpp @@ -1567,6 +1567,24 @@ project_boundaries : Optional[str] = None return old; })); + py::class_ (m, "_DebugParameters") + .def_readwrite("debugoutput", &DebugParameters::debugoutput) + .def_readwrite("slowchecks", &DebugParameters::slowchecks) + .def_readwrite("haltsuccess", &DebugParameters::haltsuccess) + .def_readwrite("haltnosuccess", &DebugParameters::haltnosuccess) + .def_readwrite("haltlargequalclass", &DebugParameters::haltlargequalclass) + .def_readwrite("haltsegment", &DebugParameters::haltsegment) + .def_readwrite("haltnode", &DebugParameters::haltnode) + .def_readwrite("haltsegmentp1", &DebugParameters::haltsegmentp1) + .def_readwrite("haltsegmentp2", &DebugParameters::haltsegmentp2) + .def_readwrite("haltexistingline", &DebugParameters::haltexistingline) + .def_readwrite("haltoverlap", &DebugParameters::haltoverlap) + .def_readwrite("haltface", &DebugParameters::haltface) + .def_readwrite("haltfacenr", &DebugParameters::haltfacenr) + .def_readwrite("write_mesh_on_error", &DebugParameters::write_mesh_on_error) + ; + + m.attr("debugparam") = py::cast(&debugparam); m.def("ReadCGNSFile", &ReadCGNSFile, py::arg("filename"), py::arg("base")=1, "Read mesh and solution vectors from CGNS file"); m.def("WriteCGNSFile", &WriteCGNSFile, py::arg("mesh"), py::arg("filename"), py::arg("names"), py::arg("values"), py::arg("locations"),