mirror of
https://github.com/NGSolve/netgen.git
synced 2024-12-24 04:50:34 +05:00
New debug parameter to write mesh on error, python export
This commit is contained in:
parent
7712429cc2
commit
8224f3cd2d
@ -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] );
|
||||
|
@ -2898,6 +2898,7 @@ namespace netgen
|
||||
haltsegment = 0;
|
||||
haltsegmentp1 = 0;
|
||||
haltsegmentp2 = 0;
|
||||
write_mesh_on_error = false;
|
||||
};
|
||||
|
||||
|
||||
|
@ -1415,6 +1415,8 @@ namespace netgen
|
||||
///
|
||||
int haltfacenr;
|
||||
///
|
||||
bool write_mesh_on_error;
|
||||
///
|
||||
DebugParameters ();
|
||||
};
|
||||
|
||||
|
@ -1567,6 +1567,24 @@ project_boundaries : Optional[str] = None
|
||||
return old;
|
||||
}));
|
||||
|
||||
py::class_<DebugParameters> (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"),
|
||||
|
Loading…
Reference in New Issue
Block a user