#ifdef NG_PYTHON #include <../general/ngpython.hpp> #include #include #include "../meshing/python_mesh.hpp" #ifdef WIN32 #define DLL_HEADER __declspec(dllexport) #endif using namespace netgen; namespace netgen { //extern shared_ptr mesh; extern shared_ptr ng_geometry; } static string stlparameter_description = R"delimiter( STL Specific Meshing Parameters ------------------------------- yangle: float = Angle for edge detection contyangle: float = Edges continue if angle > contyangle edgecornerangle: float = Angle of geometry edge at which the mesher should set a point. )delimiter"; void CreateSTLParametersFromKwargs(STLParameters& stlparam, py::kwargs kwargs) { if(kwargs.contains("yangle")) stlparam.yangle = py::cast(kwargs.attr("pop")("yangle")); if(kwargs.contains("contyangle")) stlparam.contyangle = py::cast(kwargs.attr("pop")("contyangle")); if(kwargs.contains("edgecornerangle")) stlparam.edgecornerangle = py::cast(kwargs.attr("pop")("edgecornerangle")); if(kwargs.contains("chartangle")) stlparam.chartangle = py::cast(kwargs.attr("pop")("chartangle")); if(kwargs.contains("outerchartangle")) stlparam.outerchartangle = py::cast(kwargs.attr("pop")("outerchartangle")); if(kwargs.contains("usesearchtree")) stlparam.usesearchtree = py::cast(kwargs.attr("pop")("usesearchtree")); if(kwargs.contains("resthatlasfac")) { auto val = kwargs.attr("pop")("resthatlasfac"); if(val.is_none()) stlparam.resthatlasenable = false; else { stlparam.resthatlasenable = true; stlparam.resthatlasfac = py::cast(val); } } if(kwargs.contains("atlasminh")) stlparam.atlasminh = py::cast(kwargs.attr("pop")("atlasminh")); if(kwargs.contains("resthsurfcurvfac")) { auto val = kwargs.attr("pop")("resthsurfcurvfac"); if(val.is_none()) stlparam.resthsurfcurvenable = false; else { stlparam.resthsurfcurvenable = true; stlparam.resthsurfcurvfac = py::cast(val); } } if(kwargs.contains("resthchartdistfac")) { auto val = kwargs.attr("pop")("resthchartdistfac"); if(val.is_none()) stlparam.resthchartdistenable = false; else { stlparam.resthchartdistenable = true; stlparam.resthchartdistfac = py::cast(val); } } if(kwargs.contains("resthcloseedgefac")) { auto val = kwargs.attr("pop")("resthcloseedgefac"); if(val.is_none()) stlparam.resthcloseedgeenable = false; else { stlparam.resthcloseedgeenable = true; stlparam.resthcloseedgefac = py::cast(val); } } if(kwargs.contains("resthedgeanglefac")) { auto val = kwargs.attr("pop")("resthedgeanglefac"); if(val.is_none()) stlparam.resthedgeangleenable = false; else { stlparam.resthedgeangleenable = true; stlparam.resthedgeanglefac = py::cast(val); } } if(kwargs.contains("resthsurfmeshcurvfac")) { auto val = kwargs.attr("pop")("resthsurfmeshcurvfac"); if(val.is_none()) stlparam.resthsurfmeshcurvenable = false; else { stlparam.resthsurfmeshcurvenable = true; stlparam.resthsurfmeshcurvfac = py::cast(val); } } if(kwargs.contains("resthlinelengthfac")) { auto val = kwargs.attr("pop")("resthlinelengthfac"); if(val.is_none()) stlparam.resthlinelengthenable = false; else { stlparam.resthlinelengthenable = true; stlparam.resthlinelengthfac = py::cast(val); } } if(kwargs.contains("recalc_h_opt")) stlparam.recalc_h_opt = py::cast(kwargs.attr("pop")("recalc_h_opt")); } DLL_HEADER void ExportSTL(py::module & m) { py::class_, NetgenGeometry> (m,"STLGeometry") .def(py::init<>()) .def(py::init<>([](const string& filename) { ifstream ist(filename); return shared_ptr(STLGeometry::Load(ist)); }), py::arg("filename"), py::call_guard()) .def(NGSPickle()) .def("_visualizationData", [](shared_ptr stl_geo) { std::vector vertices; std::vector trigs; std::vector normals; std::vector min = {std::numeric_limits::max(), std::numeric_limits::max(), std::numeric_limits::max()}; std::vector max = {std::numeric_limits::lowest(), std::numeric_limits::lowest(), std::numeric_limits::lowest()}; std::vector surfnames; surfnames.push_back("stl"); vertices.reserve(stl_geo->GetNT()*3*3); trigs.reserve(stl_geo->GetNT()*4); normals.reserve(stl_geo->GetNT()*3*3); size_t ii = 0; for(int i = 0; i < stl_geo->GetNT(); i++) { auto& trig = stl_geo->GetTriangle(i+1); for(int k = 0; k < 3; k++) { trigs.push_back(ii++); auto& pnt = stl_geo->GetPoint(trig[k]); for (int l = 0; l < 3; l++) { float val = pnt[l]; vertices.push_back(val); min[l] = min2(min[l], val); max[l] = max2(max[l], val); normals.push_back(trig.Normal()[l]); } } trigs.push_back(0); } py::gil_scoped_acquire ac; py::dict res; py::list snames; for(auto name : surfnames) snames.append(py::cast(name)); res["vertices"] = MoveToNumpy(vertices); res["triangles"] = MoveToNumpy(trigs); res["normals"] = MoveToNumpy(normals); res["surfnames"] = snames; res["min"] = MoveToNumpy(min); res["max"] = MoveToNumpy(max); return res; }, py::call_guard()) .def("GenerateMesh", [] (shared_ptr geo, MeshingParameters* pars, py::kwargs kwargs) { MeshingParameters mp; { py::gil_scoped_acquire aq; STLParameters stlparam; if(pars) { if(pars->geometrySpecificParameters.has_value() && (pars->geometrySpecificParameters.type() == typeid(py::kwargs))) { auto mp_kwargs = any_cast(pars->geometrySpecificParameters); py::print("geometry specific kwargs:", mp_kwargs); CreateSTLParametersFromKwargs(stlparam, mp_kwargs); pars->geometrySpecificParameters.reset(); } mp = *pars; } CreateSTLParametersFromKwargs(stlparam, kwargs); CreateMPfromKwargs(mp, kwargs); // this will throw if any kwargs are not passed mp.geometrySpecificParameters = stlparam; } auto mesh = make_shared(); mesh->SetGeometry(geo); ng_geometry = geo; SetGlobalMesh(mesh); geo->GenerateMesh(mesh,mp); return mesh; }, py::arg("mp") = nullptr, py::call_guard(), (meshingparameter_description + stlparameter_description).c_str()) ; m.def("LoadSTLGeometry", [] (const string & filename) { cout << "WARNING: LoadSTLGeometry is deprecated, use the STLGeometry(filename) constructor instead!" << endl; ifstream ist(filename); return shared_ptr(STLGeometry::Load(ist)); },py::call_guard()); } PYBIND11_MODULE(libstl, m) { ExportSTL(m); } #endif