netgen/libsrc/csg/python_csg.cpp

790 lines
32 KiB
C++
Raw Normal View History

2014-08-31 19:05:24 +06:00
#ifdef NG_PYTHON
2014-08-30 06:15:59 +06:00
2014-10-08 20:09:03 +06:00
#include <../general/ngpython.hpp>
#include <core/python_ngcore.hpp>
2014-08-30 06:15:59 +06:00
#include <csg.hpp>
#include "../meshing/python_mesh.hpp"
2014-08-30 06:15:59 +06:00
using namespace netgen;
2014-12-19 19:03:36 +05:00
namespace netgen
{
extern shared_ptr<NetgenGeometry> ng_geometry;
}
2014-08-30 06:15:59 +06:00
2016-08-26 17:33:57 +05:00
2014-08-30 06:15:59 +06:00
2014-09-25 16:50:48 +06:00
// a shadow solid tree using shared pointers.
class SPSolid
{
shared_ptr<SPSolid> s1, s2;
Solid * solid;
2014-10-01 02:36:02 +06:00
int bc = -1;
2015-08-08 22:10:48 +05:00
string bcname = "";
2014-10-18 20:41:20 +06:00
double maxh = -1;
2014-10-05 22:53:19 +06:00
string material;
2014-09-25 16:50:48 +06:00
bool owner;
2017-09-03 01:13:36 +05:00
double red = 0, green = 0, blue = 1;
2014-10-29 00:56:11 +05:00
bool transp = false;
2014-09-25 16:50:48 +06:00
public:
2019-07-15 12:07:35 +05:00
enum optyp { TERM, SECTION, UNION, SUB, EXISTING };
2014-09-25 16:50:48 +06:00
SPSolid (Solid * as) : solid(as), owner(true), op(TERM) { ; }
2019-07-15 12:07:35 +05:00
SPSolid (Solid * as, int /*dummy*/)
: solid(as), owner(false), op(EXISTING) { ; }
2014-09-26 02:23:31 +06:00
~SPSolid ()
{
2014-10-01 02:36:02 +06:00
; // if (owner) delete solid;
2014-09-26 02:23:31 +06:00
}
2014-09-25 16:50:48 +06:00
SPSolid (optyp aop, shared_ptr<SPSolid> as1, shared_ptr<SPSolid> as2)
: s1(as1), s2(as2), owner(true), op(aop)
{
if (aop == UNION)
solid = new Solid (Solid::UNION, s1->GetSolid(), s2->GetSolid());
else if (aop == SECTION)
solid = new Solid (Solid::SECTION, s1->GetSolid(), s2->GetSolid());
2014-09-25 20:42:36 +06:00
else if (aop == SUB)
2014-09-26 02:23:31 +06:00
solid = new Solid (Solid::SUB, s1->GetSolid()); // , s2->GetSolid());
2014-09-25 16:50:48 +06:00
}
Solid * GetSolid() { return solid; }
2014-10-01 18:07:09 +06:00
const Solid * GetSolid() const { return solid; }
2014-09-25 16:50:48 +06:00
void GiveUpOwner()
{
owner = false;
if (s1) s1 -> GiveUpOwner();
if (s2) s2 -> GiveUpOwner();
}
void AddSurfaces(CSGeometry & geom)
{
if (op == TERM)
geom.AddSurfaces (solid->GetPrimitive());
if (s1) s1 -> AddSurfaces (geom);
if (s2) s2 -> AddSurfaces (geom);
}
2014-10-05 22:53:19 +06:00
void SetMaterial (string mat) { material = mat; }
string GetMaterial ()
{
if (!material.empty()) return material;
if (s1)
{
string s1mat = s1->GetMaterial();
if (!s1mat.empty()) return s1mat;
}
if (s2)
{
string s2mat = s2->GetMaterial();
if (!s2mat.empty()) return s2mat;
}
return material;
}
2014-10-01 02:36:02 +06:00
void SetBC(int abc)
{
if (bc == -1)
{
bc = abc;
if (s1) s1 -> SetBC(bc);
if (s2) s2 -> SetBC(bc);
if (op == TERM)
{
Primitive * prim = solid -> GetPrimitive();
for (int i = 0; i < prim->GetNSurfaces(); i++)
prim->GetSurface(i).SetBCProperty (abc);
2014-12-09 19:42:53 +05:00
// cout << "set " << prim->GetNSurfaces() << " surfaces to bc " << bc << endl;
2014-10-01 02:36:02 +06:00
}
}
}
2014-10-18 20:41:20 +06:00
2015-08-08 22:10:48 +05:00
void SetBCName(string name)
{
if (bcname == "")
{
bcname = name;
if (s1) s1 -> SetBCName(name);
if (s2) s2 -> SetBCName(name);
if (op == TERM)
{
Primitive * prim = solid -> GetPrimitive();
for (int i = 0; i < prim->GetNSurfaces(); i++)
prim->GetSurface(i).SetBCName (name);
// cout << "set " << prim->GetNSurfaces() << " surfaces to bc " << bc << endl;
}
}
}
2014-12-05 20:58:49 +05:00
void SetMaxH(double amaxh)
2014-10-18 20:41:20 +06:00
{
if (maxh == -1)
{
maxh = amaxh;
if (s1) s1 -> SetMaxH(maxh);
if (s2) s2 -> SetMaxH(maxh);
if (op == TERM)
{
Primitive * prim = solid -> GetPrimitive();
for (int i = 0; i < prim->GetNSurfaces(); i++)
prim->GetSurface(i).SetMaxH (maxh);
}
}
}
2014-10-29 00:56:11 +05:00
void SetColor(double ared, double agreen, double ablue)
{
red = ared;
green = agreen;
blue = ablue;
2014-10-29 00:56:11 +05:00
}
double GetRed() const { return red; }
double GetGreen() const { return green; }
double GetBlue() const { return blue; }
void SetTransparent() { transp = true; }
bool IsTransparent() { return transp; }
2014-10-18 20:41:20 +06:00
2014-09-25 16:50:48 +06:00
private:
optyp op;
};
2014-08-30 06:15:59 +06:00
2014-10-01 18:07:09 +06:00
inline ostream & operator<< (ostream & ost, const SPSolid & sol)
{
ost << *sol.GetSolid();
return ost;
}
2014-08-30 06:15:59 +06:00
namespace netgen
{
2016-11-04 16:14:52 +05:00
extern CSGeometry * ParseCSG (istream & istr, CSGeometry *instance=nullptr);
2014-08-30 06:15:59 +06:00
}
2016-10-30 19:01:52 +05:00
2016-11-04 16:14:52 +05:00
DLL_HEADER void ExportCSG(py::module &m)
2014-08-30 06:15:59 +06:00
{
2016-11-04 16:14:52 +05:00
py::class_<SplineGeometry<2>> (m, "SplineCurve2d")
.def(py::init<>())
2016-09-06 19:21:05 +05:00
.def ("AddPoint", FunctionPointer
([] (SplineGeometry<2> & self, double x, double y)
{
self.geompoints.Append (GeomPoint<2> (Point<2> (x,y)));
return self.geompoints.Size()-1;
}))
.def ("AddSegment", FunctionPointer
([] (SplineGeometry<2> & self, int i1, int i2)
{
self.splines.Append (new LineSeg<2> (self.geompoints[i1], self.geompoints[i2]));
}))
.def ("AddSegment", FunctionPointer
([] (SplineGeometry<2> & self, int i1, int i2, int i3)
{
self.splines.Append (new SplineSeg3<2> (self.geompoints[i1], self.geompoints[i2], self.geompoints[i3]));
}))
;
2016-11-06 22:25:38 +05:00
py::class_<SplineGeometry<3>,shared_ptr<SplineGeometry<3>>> (m,"SplineCurve3d")
.def(py::init<>())
2016-10-04 22:30:57 +05:00
.def ("AddPoint", FunctionPointer
([] (SplineGeometry<3> & self, double x, double y, double z)
{
self.geompoints.Append (GeomPoint<3> (Point<3> (x,y,z)));
return self.geompoints.Size()-1;
}))
.def ("AddSegment", FunctionPointer
([] (SplineGeometry<3> & self, int i1, int i2)
{
self.splines.Append (new LineSeg<3> (self.geompoints[i1], self.geompoints[i2]));
}))
.def ("AddSegment", FunctionPointer
([] (SplineGeometry<3> & self, int i1, int i2, int i3)
{
self.splines.Append (new SplineSeg3<3> (self.geompoints[i1], self.geompoints[i2], self.geompoints[i3]));
}))
;
2016-11-06 22:25:38 +05:00
py::class_<SplineSurface, shared_ptr<SplineSurface>> (m, "SplineSurface",
"A surface for co dim 2 integrals on the splines")
.def(py::init([](shared_ptr<SPSolid> base, py::list cuts)
2016-10-27 18:41:08 +05:00
{
auto primitive = dynamic_cast<OneSurfacePrimitive*> (base->GetSolid()->GetPrimitive());
2019-07-09 13:39:16 +05:00
auto acuts = make_shared<NgArray<shared_ptr<OneSurfacePrimitive>>>();
2016-11-06 22:25:38 +05:00
for(int i = 0; i<py::len(cuts);i++)
2016-10-27 18:41:08 +05:00
{
2016-11-06 22:25:38 +05:00
py::extract<shared_ptr<SPSolid>> sps(cuts[i]);
2016-10-27 18:41:08 +05:00
if(!sps.check())
throw NgException("Cut must be SurfacePrimitive in constructor of SplineSurface!");
auto sp = dynamic_cast<OneSurfacePrimitive*>(sps()->GetSolid()->GetPrimitive());
if(sp)
2017-02-27 15:32:42 +05:00
acuts->Append(shared_ptr<OneSurfacePrimitive>(sp));
2016-10-27 18:41:08 +05:00
else
throw Exception("Cut must be SurfacePrimitive in constructor of SplineSurface!");
2016-10-27 18:41:08 +05:00
}
if(!primitive)
throw Exception("Base is not a SurfacePrimitive in constructor of SplineSurface!");
return make_shared<SplineSurface>(shared_ptr<OneSurfacePrimitive>(primitive),acuts);
2016-11-09 20:10:04 +05:00
}),py::arg("base"), py::arg("cuts")=py::list())
2016-10-04 22:30:57 +05:00
.def("AddPoint", FunctionPointer
([] (SplineSurface & self, double x, double y, double z, bool hpref)
{
self.AppendPoint(Point<3>(x,y,z),hpref);
return self.GetNP()-1;
}),
2016-11-06 22:25:38 +05:00
py::arg("x"),py::arg("y"),py::arg("z"),py::arg("hpref")=false)
.def("AddSegment", [] (SplineSurface & self, int i1, int i2, string bcname, double maxh)
2016-10-04 22:30:57 +05:00
{
auto seg = make_shared<LineSeg<3>>(self.GetPoint(i1),self.GetPoint(i2));
self.AppendSegment(seg,bcname,maxh);
},
2016-11-06 22:25:38 +05:00
py::arg("pnt1"),py::arg("pnt2"),py::arg("bcname")="default", py::arg("maxh")=-1.)
.def("AddSegment", [] (SplineSurface& self, int i1, int i2, int i3, string bcname, double maxh)
{
auto seg = make_shared<SplineSeg3<3>>(self.GetPoint(i1), self.GetPoint(i2), self.GetPoint(i3));
self.AppendSegment(seg, bcname, maxh);
}, py::arg("pnt1"),py::arg("pnt2"), py::arg("pnt3"),py::arg("bcname")="default", py::arg("maxh")=-1.)
2016-10-04 22:30:57 +05:00
;
2016-09-06 19:21:05 +05:00
2016-11-04 16:14:52 +05:00
py::class_<SPSolid, shared_ptr<SPSolid>> (m, "Solid")
2014-10-01 18:07:09 +06:00
.def ("__str__", &ToString<SPSolid>)
2014-12-09 17:01:39 +05:00
.def ("__add__", FunctionPointer( [] ( shared_ptr<SPSolid> self, shared_ptr<SPSolid> other) { return make_shared<SPSolid> (SPSolid::UNION, self, other); }))
.def ("__mul__", FunctionPointer( [] ( shared_ptr<SPSolid> self, shared_ptr<SPSolid> other) { return make_shared<SPSolid> (SPSolid::SECTION, self, other); }))
.def ("__sub__", FunctionPointer
([] ( shared_ptr<SPSolid> self, shared_ptr<SPSolid> other)
{ return make_shared<SPSolid> (SPSolid::SECTION, self,
make_shared<SPSolid> (SPSolid::SUB, other, nullptr)); }))
2014-10-01 18:07:09 +06:00
2014-10-01 02:36:02 +06:00
.def ("bc", FunctionPointer([](shared_ptr<SPSolid> & self, int nr) -> shared_ptr<SPSolid>
{ self->SetBC(nr); return self; }))
2015-08-08 22:10:48 +05:00
.def ("bc", FunctionPointer([](shared_ptr<SPSolid> & self, string name) -> shared_ptr<SPSolid>
{ self->SetBCName(name); return self; }))
2014-10-18 20:41:20 +06:00
.def ("maxh", FunctionPointer([](shared_ptr<SPSolid> & self, double maxh) -> shared_ptr<SPSolid>
{ self->SetMaxH(maxh); return self; }))
2014-10-05 22:53:19 +06:00
.def ("mat", FunctionPointer([](shared_ptr<SPSolid> & self, string mat) -> shared_ptr<SPSolid>
{ self->SetMaterial(mat); return self; }))
.def ("mat", &SPSolid::GetMaterial)
2016-11-04 16:14:52 +05:00
.def("col", FunctionPointer([](shared_ptr<SPSolid> & self, py::list rgb) -> shared_ptr<SPSolid>
2014-12-09 17:01:39 +05:00
{
2016-11-04 16:14:52 +05:00
py::extract<double> red(rgb[0]);
py::extract<double> green(rgb[1]);
py::extract<double> blue(rgb[2]);
2014-12-09 17:01:39 +05:00
self->SetColor(red(),green(),blue());
return self;
}))
.def("transp", FunctionPointer([](shared_ptr<SPSolid> & self)->shared_ptr < SPSolid > { self->SetTransparent(); return self; }))
2014-10-05 22:53:19 +06:00
;
2014-09-25 16:50:48 +06:00
2016-11-04 16:14:52 +05:00
m.def ("Sphere", FunctionPointer([](Point<3> c, double r)
2014-09-25 16:50:48 +06:00
{
Sphere * sp = new Sphere (c, r);
Solid * sol = new Solid (sp);
return make_shared<SPSolid> (sol);
}));
2017-12-07 13:43:32 +05:00
m.def ("Ellipsoid", FunctionPointer([](Point<3> m, Vec<3> a, Vec<3> b, Vec<3> c)
{
Ellipsoid * ell = new Ellipsoid (m, a, b, c);
Solid * sol = new Solid (ell);
return make_shared<SPSolid> (sol);
}));
2016-11-04 16:14:52 +05:00
m.def ("Plane", FunctionPointer([](Point<3> p, Vec<3> n)
2014-09-25 16:50:48 +06:00
{
Plane * sp = new Plane (p,n);
Solid * sol = new Solid (sp);
return make_shared<SPSolid> (sol);
}));
2016-11-04 16:14:52 +05:00
m.def ("Cone", FunctionPointer([](Point<3> a, Point<3> b, double ra, double rb)
2015-09-01 02:47:27 +05:00
{
Cone * cyl = new Cone (a, b, ra, rb);
Solid * sol = new Solid (cyl);
return make_shared<SPSolid> (sol);
}));
2016-11-04 16:14:52 +05:00
m.def ("Cylinder", FunctionPointer([](Point<3> a, Point<3> b, double r)
2014-09-26 02:23:31 +06:00
{
Cylinder * cyl = new Cylinder (a, b, r);
Solid * sol = new Solid (cyl);
return make_shared<SPSolid> (sol);
}));
2016-11-04 16:14:52 +05:00
m.def ("OrthoBrick", FunctionPointer([](Point<3> p1, Point<3> p2)
2014-09-25 16:50:48 +06:00
{
OrthoBrick * brick = new OrthoBrick (p1,p2);
Solid * sol = new Solid (brick);
return make_shared<SPSolid> (sol);
}));
m.def ("Torus", FunctionPointer([](Point<3> c, Vec<3> n, double R, double r)
{
Torus * torus = new Torus (c,n,R,r);
Solid * sol = new Solid (torus);
return make_shared<SPSolid> (sol);
}));
2016-11-04 16:14:52 +05:00
m.def ("Revolution", FunctionPointer([](Point<3> p1, Point<3> p2,
2016-09-06 19:21:05 +05:00
const SplineGeometry<2> & spline)
{
Revolution * rev = new Revolution (p1, p2, spline);
Solid * sol = new Solid(rev);
return make_shared<SPSolid> (sol);
}));
m.def ("Extrusion", FunctionPointer([](const SplineGeometry<3> & path,
const SplineGeometry<2> & profile,
Vec<3> n)
{
Extrusion * extr = new Extrusion (path,profile,n);
Solid * sol = new Solid(extr);
return make_shared<SPSolid> (sol);
}));
2018-02-26 19:22:35 +05:00
m.def("EllipticCone", [](const Point<3>& a, const Vec<3>& v, const Vec<3>& w,
double h, double r)
{
auto ellcone = new EllipticCone(a,v,w,h,r);
auto sol = new Solid(ellcone);
return make_shared<SPSolid>(sol);
}, py::arg("a"), py::arg("vl"), py::arg("vs"), py::arg("h"), py::arg("r"),
R"raw_string(
An elliptic cone, given by the point 'a' at the base of the cone along the main axis,
the vectors v and w of the long and short axis of the ellipse, respectively,
the height of the cone, h, and ratio of base long axis length to top long axis length, r
Note: The elliptic cone has to be truncated by planes similar to a cone or an elliptic cylinder.
When r =1, the truncated elliptic cone becomes an elliptic cylinder.
When r tends to zero, the truncated elliptic cone tends to a full elliptic cone.
However, when r = 0, the top part becomes a point(tip) and meshing fails!
)raw_string");
m.def("Polyhedron", [](py::list points, py::list faces)
{
auto poly = new Polyhedra();
for(auto p : points)
poly->AddPoint(py::cast<Point<3>>(p));
int fnr = 0;
for(auto face : faces)
{
auto lface = py::cast<py::list>(face);
if(py::len(lface) == 3)
poly->AddFace(py::cast<int>(lface[0]),
py::cast<int>(lface[1]),
py::cast<int>(lface[2]),
fnr++);
else if(py::len(lface) == 4)
{
poly->AddFace(py::cast<int>(lface[0]),
py::cast<int>(lface[1]),
py::cast<int>(lface[2]),
fnr);
poly->AddFace(py::cast<int>(lface[0]),
py::cast<int>(lface[2]),
py::cast<int>(lface[3]),
fnr++);
}
}
return make_shared<SPSolid>(new Solid(poly));
});
2014-09-25 16:50:48 +06:00
2016-11-04 16:14:52 +05:00
m.def ("Or", FunctionPointer([](shared_ptr<SPSolid> s1, shared_ptr<SPSolid> s2)
2014-09-25 16:50:48 +06:00
{
return make_shared<SPSolid> (SPSolid::UNION, s1, s2);
}));
2016-11-04 16:14:52 +05:00
m.def ("And", FunctionPointer([](shared_ptr<SPSolid> s1, shared_ptr<SPSolid> s2)
2014-09-25 16:50:48 +06:00
{
return make_shared<SPSolid> (SPSolid::SECTION, s1, s2);
}));
2014-08-30 06:15:59 +06:00
py::class_<CSGeometry, NetgenGeometry, shared_ptr<CSGeometry>> (m, "CSGeometry")
2016-11-04 16:14:52 +05:00
.def(py::init<>())
2018-11-29 22:35:30 +05:00
.def(py::init([](const string& filename)
{
ifstream ist (filename);
auto geo = make_shared<CSGeometry>();
ParseCSG(ist, geo.get());
geo->FindIdenticSurfaces(1e-8 * geo->MaxSize());
return geo;
}), py::arg("filename"))
.def(NGSPickle<CSGeometry>())
2014-09-26 02:23:31 +06:00
.def("Save", FunctionPointer([] (CSGeometry & self, string filename)
{
cout << "save geometry to file " << filename << endl;
self.Save (filename);
}))
2016-11-04 16:14:52 +05:00
.def("Add",
2017-09-01 12:12:50 +05:00
[] (CSGeometry & self, shared_ptr<SPSolid> solid, py::list bcmod, double maxh,
py::tuple col, bool transparent, int layer)
2014-10-05 22:53:19 +06:00
{
solid->AddSurfaces (self);
solid->GiveUpOwner();
int tlonr = self.SetTopLevelObject (solid->GetSolid());
self.GetTopLevelObject(tlonr) -> SetMaterial(solid->GetMaterial());
2014-12-09 17:01:39 +05:00
self.GetTopLevelObject(tlonr) -> SetRGB(solid->GetRed(),solid->GetGreen(),solid->GetBlue());
2017-09-03 01:13:36 +05:00
// self.GetTopLevelObject(tlonr)->SetTransparent(solid->IsTransparent());
self.GetTopLevelObject(tlonr)->SetTransparent(transparent);
self.GetTopLevelObject(tlonr)->SetMaxH(maxh);
self.GetTopLevelObject(tlonr)->SetLayer(layer);
2014-12-11 14:04:49 +05:00
2017-09-01 12:12:50 +05:00
// cout << "rgb = " << py::len(rgb) << endl;
if (py::len(col)==3)
self.GetTopLevelObject(tlonr) -> SetRGB(py::cast<double>(col[0]),
py::cast<double>(col[1]),
py::cast<double>(col[2]));
2014-12-11 14:04:49 +05:00
// bcmod is list of tuples ( solid, bcnr )
2016-11-04 16:14:52 +05:00
for (int i = 0; i < py::len(bcmod); i++)
2014-12-11 14:04:49 +05:00
{
2016-11-04 16:14:52 +05:00
py::tuple tup = py::extract<py::tuple> (bcmod[i]) ();
auto mod_solid = py::extract<shared_ptr<SPSolid>> (tup[0]) ();
2015-08-26 18:44:37 +05:00
int mod_nr = -1;
string * bcname = nullptr;
2016-11-04 16:14:52 +05:00
py::object val = tup[1];
if (py::extract<int>(val).check()) mod_nr = py::extract<int> (val)();
if (py::extract<string>(val).check()) bcname = new string ( py::extract<string> (val)());
2015-08-26 18:44:37 +05:00
2019-07-09 13:39:16 +05:00
NgArray<int> si;
2014-12-11 14:04:49 +05:00
mod_solid -> GetSolid() -> GetSurfaceIndices (si);
2014-12-19 19:03:36 +05:00
// cout << "change bc on surfaces: " << si << " to " << mod_nr << endl;
2014-12-11 14:04:49 +05:00
for (int j = 0; j < si.Size(); j++)
{
CSGeometry::BCModification bcm;
2015-08-26 18:44:37 +05:00
bcm.bcname = bcname ? new string (*bcname) : nullptr;
2014-12-11 14:04:49 +05:00
bcm.tlonr = tlonr;
bcm.si = si[j];
2014-12-19 19:03:36 +05:00
bcm.bcnr = mod_nr;
2014-12-11 14:04:49 +05:00
self.bcmodifications.Append (bcm);
}
2015-08-26 18:44:37 +05:00
delete bcname;
2014-12-11 14:04:49 +05:00
}
return tlonr;
2016-11-04 16:14:52 +05:00
},
2017-09-01 12:12:50 +05:00
py::arg("solid"), py::arg("bcmod")=py::list(), py::arg("maxh")=1e99,
py::arg("col")=py::tuple(), py::arg("transparent")=false, py::arg("layer")=1
2014-12-11 14:04:49 +05:00
)
2014-09-25 16:50:48 +06:00
2015-11-19 16:30:54 +05:00
.def("AddSurface", FunctionPointer
([] (CSGeometry & self, shared_ptr<SPSolid> surface, shared_ptr<SPSolid> solid)
{
solid->AddSurfaces (self);
solid->GiveUpOwner();
Surface & surf = surface->GetSolid()->GetPrimitive()->GetSurface();
int tlonr = self.SetTopLevelObject (solid->GetSolid(), &surf);
// self.GetTopLevelObject(tlonr) -> SetMaterial(solid->GetMaterial());
self.GetTopLevelObject(tlonr) -> SetBCProp(surf.GetBCProperty());
self.GetTopLevelObject(tlonr) -> SetBCName(surf.GetBCName());
2015-11-19 16:30:54 +05:00
self.GetTopLevelObject(tlonr) -> SetRGB(solid->GetRed(),solid->GetGreen(),solid->GetBlue());
self.GetTopLevelObject(tlonr)->SetTransparent(solid->IsTransparent());
}),
2016-11-04 16:14:52 +05:00
py::arg("surface"), py::arg("solid")
2015-11-19 16:30:54 +05:00
)
2016-10-04 22:30:57 +05:00
.def("AddSplineSurface", FunctionPointer
([] (CSGeometry & self, shared_ptr<SplineSurface> surf)
{
2016-10-27 18:41:08 +05:00
auto cuttings = surf->CreateCuttingSurfaces();
auto spsol = make_shared<SPSolid>(new Solid(surf.get()));
2016-10-27 18:41:08 +05:00
for(auto cut : (*cuttings)){
spsol = make_shared<SPSolid>(SPSolid::SECTION,spsol,make_shared<SPSolid>(new Solid(cut.get())));
2016-10-04 22:30:57 +05:00
}
spsol->AddSurfaces(self);
int tlonr = self.SetTopLevelObject(spsol->GetSolid(), surf.get());
self.GetTopLevelObject(tlonr) -> SetBCProp(surf->GetBase()->GetBCProperty());
self.GetTopLevelObject(tlonr) -> SetBCName(surf->GetBase()->GetBCName());
self.GetTopLevelObject(tlonr) -> SetMaxH(surf->GetBase()->GetMaxH());
2019-07-09 13:39:16 +05:00
NgArray<Point<3>> non_midpoints;
for(auto spline : surf->GetSplines())
{
non_midpoints.Append(spline->GetPoint(0));
}
for(auto p : non_midpoints)
2016-10-27 18:41:08 +05:00
self.AddUserPoint(p);
self.AddSplineSurface(surf);
2016-10-04 22:30:57 +05:00
}),
2016-11-06 22:25:38 +05:00
py::arg("SplineSurface"))
2019-01-10 13:41:42 +05:00
.def("SingularFace", [] (CSGeometry & self, shared_ptr<SPSolid> sol, shared_ptr<SPSolid> surfaces, double factor)
{
int tlonum = -1;
for (int i = 0; i < self.GetNTopLevelObjects(); i++)
if (self.GetTopLevelObject(i)->GetSolid() == sol->GetSolid())
tlonum = i;
if (tlonum == -1) throw NgException("not a top-level-object");
if (!surfaces) surfaces = sol;
2019-01-25 14:41:49 +05:00
auto singface = new SingularFace(tlonum+1, surfaces->GetSolid(), factor);
2019-01-10 13:41:42 +05:00
self.singfaces.Append(singface);
}, py::arg("solid"), py::arg("surfaces")=nullptr, py::arg("factor")=0.25)
.def("SingularEdge", [] (CSGeometry & self, shared_ptr<SPSolid> s1,shared_ptr<SPSolid> s2, double factor)
{
auto singedge = new SingularEdge(1, -1, self, s1->GetSolid(), s2->GetSolid(), factor);
self.singedges.Append (singedge);
})
.def("SingularPoint", [] (CSGeometry & self, shared_ptr<SPSolid> s1,shared_ptr<SPSolid> s2,
shared_ptr<SPSolid> s3, double factor)
{
auto singpoint = new SingularPoint(1, s1->GetSolid(), s2->GetSolid(), s3->GetSolid(), factor);
self.singpoints.Append (singpoint);
})
.def("CloseSurfaces", FunctionPointer
2016-11-04 16:14:52 +05:00
([] (CSGeometry & self, shared_ptr<SPSolid> s1, shared_ptr<SPSolid> s2, py::list aslices )
{
2019-07-09 13:39:16 +05:00
NgArray<int> si1, si2;
s1->GetSolid()->GetSurfaceIndices (si1);
s2->GetSolid()->GetSurfaceIndices (si2);
Flags flags;
try
{
2016-11-04 16:14:52 +05:00
int n = py::len(aslices);
2019-08-06 17:16:13 +05:00
Array<double> slices(n);
for(int i=0; i<n; i++)
{
2016-11-04 16:14:52 +05:00
slices[i]= py::extract<double>(aslices[i])();
}
flags.SetFlag("slices", slices);
}
2016-11-04 16:14:52 +05:00
catch( py::error_already_set const & ) {
cout << "caught python error:" << endl;
PyErr_Print();
}
const TopLevelObject * domain = nullptr;
self.AddIdentification
(new CloseSurfaceIdentification
(self.GetNIdentifications()+1, self,
self.GetSurface (si1[0]), self.GetSurface (si2[0]),
domain,
flags));
}),
2016-11-04 16:14:52 +05:00
py::arg("solid1"), py::arg("solid2"), py::arg("slices")
)
.def("CloseSurfaces", FunctionPointer
2017-08-27 17:52:11 +05:00
([] (CSGeometry & self, shared_ptr<SPSolid> s1, shared_ptr<SPSolid> s2,
int reflevels, shared_ptr<SPSolid> domain_solid)
{
2019-07-09 13:39:16 +05:00
NgArray<int> si1, si2;
s1->GetSolid()->GetSurfaceIndices (si1);
s2->GetSolid()->GetSurfaceIndices (si2);
cout << "surface ids1 = " << si1 << endl;
cout << "surface ids2 = " << si2 << endl;
Flags flags;
2018-01-14 12:01:57 +05:00
const TopLevelObject * domain = nullptr;
if (domain_solid)
domain = self.GetTopLevelObject(domain_solid->GetSolid());
2017-08-27 17:52:11 +05:00
self.AddIdentification
(new CloseSurfaceIdentification
(self.GetNIdentifications()+1, self,
self.GetSurface (si1[0]), self.GetSurface (si2[0]),
domain,
flags));
}),
2017-08-27 17:52:11 +05:00
py::arg("solid1"), py::arg("solid2"), py::arg("reflevels")=2, py::arg("domain")=nullptr
)
.def("PeriodicSurfaces", FunctionPointer
([] (CSGeometry & self, shared_ptr<SPSolid> s1, shared_ptr<SPSolid> s2,
Transformation<3> trafo)
{
2019-07-09 13:39:16 +05:00
NgArray<int> si1, si2;
s1->GetSolid()->GetSurfaceIndices (si1);
s2->GetSolid()->GetSurfaceIndices (si2);
cout << "identify surfaces " << si1[0] << " and " << si2[0] << endl;
self.AddIdentification
(new PeriodicIdentification
(self.GetNIdentifications()+1, self,
self.GetSurface (si1[0]), self.GetSurface (si2[0]),
trafo));
}),
py::arg("solid1"), py::arg("solid2"),
py::arg("trafo")=Transformation<3>(Vec<3>(0,0,0))
)
2020-01-13 20:41:06 +05:00
.def("NameEdge", [] (CSGeometry & self, shared_ptr<SPSolid> s1, shared_ptr<SPSolid> s2, string name)
{
Array<Surface*> surfs1, surfs2;
s1->GetSolid()->ForEachSurface( [&surfs1] (Surface * s, bool inv) { surfs1.Append(s); });
s2->GetSolid()->ForEachSurface( [&surfs2] (Surface * s, bool inv) { surfs2.Append(s); });
for (auto s1 : surfs1)
for (auto s2 : surfs2)
self.named_edges[tuple(s1,s2)] = name;
})
2020-01-15 15:56:23 +05:00
.def("AddPoint", [] (CSGeometry & self, Point<3> p, variant<int,string> index) -> CSGeometry&
2017-09-22 19:55:10 +05:00
{
2020-01-15 15:56:23 +05:00
if (auto pint = std::get_if<int> (&index))
self.AddUserPoint(CSGeometry::UserPoint(p, *pint));
if (auto pstr = std::get_if<string> (&index))
self.AddUserPoint(CSGeometry::UserPoint(p, *pstr));
2017-09-22 19:55:10 +05:00
return self;
})
.def("GetTransparent", FunctionPointer
([] (CSGeometry & self, int tlonr)
{
return self.GetTopLevelObject(tlonr)->GetTransparent();
}),
2016-11-04 16:14:52 +05:00
py::arg("tlonr")
)
.def("SetTransparent", FunctionPointer
([] (CSGeometry & self, int tlonr, bool transparent)
{
self.GetTopLevelObject(tlonr)->SetTransparent(transparent);
}),
2016-11-04 16:14:52 +05:00
py::arg("tlonr"), py::arg("transparent")
)
.def("GetVisible", FunctionPointer
([] (CSGeometry & self, int tlonr)
{
return self.GetTopLevelObject(tlonr)->GetVisible();
}),
2016-11-04 16:14:52 +05:00
py::arg("tlonr")
)
.def("SetVisible", FunctionPointer
([] (CSGeometry & self, int tlonr, bool visible)
{
self.GetTopLevelObject(tlonr)->SetVisible(visible);
}),
2016-11-04 16:14:52 +05:00
py::arg("tlonr"), py::arg("visible")
)
.def("SetBoundingBox", FunctionPointer
([] (CSGeometry & self, Point<3> pmin, Point<3> pmax)
{
self.SetBoundingBox(Box<3> (pmin, pmax));
}),
2016-11-04 16:14:52 +05:00
py::arg("pmin"), py::arg("pmax")
)
.def("Draw", FunctionPointer
([] (shared_ptr<CSGeometry> self)
{
self->FindIdenticSurfaces(1e-6);
self->CalcTriangleApproximation(0.01, 20);
ng_geometry = self;
2016-11-04 16:14:52 +05:00
})
)
2019-07-15 12:07:35 +05:00
.def("GetSolids", [](CSGeometry& self)
{
py::list lst;
for(auto i : Range(self.GetSolids().Size()))
lst.append(make_shared<SPSolid>(self.GetSolids()[i], 1234));
return lst;
})
2016-11-04 16:14:52 +05:00
.def_property_readonly ("ntlo", &CSGeometry::GetNTopLevelObjects)
.def("_visualizationData", [](shared_ptr<CSGeometry> csg_geo)
{
std::vector<float> vertices;
std::vector<int> trigs;
std::vector<float> normals;
std::vector<float> min = {std::numeric_limits<float>::max(),
std::numeric_limits<float>::max(),
std::numeric_limits<float>::max()};
std::vector<float> max = {std::numeric_limits<float>::lowest(),
std::numeric_limits<float>::lowest(),
std::numeric_limits<float>::lowest()};
std::vector<string> surfnames;
for (int i = 0; i < csg_geo->GetNSurf(); i++)
{
auto surf = csg_geo->GetSurface(i);
surfnames.push_back(surf->GetBCName());
}
csg_geo->FindIdenticSurfaces(1e-6);
csg_geo->CalcTriangleApproximation(0.01,100);
auto nto = csg_geo->GetNTopLevelObjects();
size_t np = 0;
size_t ntrig = 0;
for (int i = 0; i < nto; i++){
np += csg_geo->GetTriApprox(i)->GetNP();
ntrig += csg_geo->GetTriApprox(i)->GetNT();
}
vertices.reserve(np*3);
trigs.reserve(ntrig*4);
normals.reserve(np*3);
int offset_points = 0;
for (int i = 0; i < nto; i++)
{
auto triapprox = csg_geo->GetTriApprox(i);
for (int j = 0; j < triapprox->GetNP(); j++)
for(int k = 0; k < 3; k++) {
float val = triapprox->GetPoint(j)[k];
vertices.push_back(val);
min[k] = min2(min[k], val);
max[k] = max2(max[k],val);
normals.push_back(triapprox->GetNormal(j)[k]);
}
for (int j = 0; j < triapprox->GetNT(); j++)
{
for(int k = 0; k < 3; k++)
trigs.push_back(triapprox->GetTriangle(j)[k]+offset_points);
trigs.push_back(triapprox->GetTriangle(j).SurfaceIndex());
}
offset_points += triapprox->GetNP();
}
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<py::gil_scoped_release>())
.def("GenerateMesh", [](shared_ptr<CSGeometry> geo,
MeshingParameters* pars, py::kwargs kwargs)
2014-08-30 06:15:59 +06:00
{
MeshingParameters mp;
if(pars) mp = *pars;
{
py::gil_scoped_acquire aq;
CreateMPfromKwargs(mp, kwargs);
}
auto mesh = make_shared<Mesh>();
SetGlobalMesh (mesh);
mesh->SetGeometry(geo);
2015-08-08 22:10:48 +05:00
ng_geometry = geo;
geo->FindIdenticSurfaces(1e-8 * geo->MaxSize());
2019-11-26 21:08:21 +05:00
auto result = geo->GenerateMesh (mesh, mp);
if(result != 0)
throw Exception("Meshing failed!");
return mesh;
}, py::arg("mp") = nullptr,
meshingparameter_description.c_str(),
py::call_guard<py::gil_scoped_release>())
2014-10-01 02:36:02 +06:00
;
2014-12-09 17:01:39 +05:00
2016-11-04 16:14:52 +05:00
m.def("Save", FunctionPointer
2014-12-09 17:01:39 +05:00
([](const Mesh & self, const string & filename, const CSGeometry & geom)
{
ostream * outfile;
if (filename.substr (filename.length()-3, 3) == ".gz")
outfile = new ogzstream (filename.c_str());
else
outfile = new ofstream (filename.c_str());
self.Save (*outfile);
*outfile << endl << endl << "endmesh" << endl << endl;
geom.SaveToMeshFile (*outfile);
delete outfile;
2018-03-13 02:38:21 +05:00
}),py::call_guard<py::gil_scoped_release>())
2014-12-09 17:01:39 +05:00
;
2016-11-04 16:14:52 +05:00
m.def("ZRefinement", FunctionPointer
([](Mesh & mesh, CSGeometry & geom)
{
ZRefinementOptions opt;
opt.minref = 5;
ZRefinement (mesh, &geom, opt);
2018-03-13 02:38:21 +05:00
}),py::call_guard<py::gil_scoped_release>())
;
2014-08-30 06:15:59 +06:00
}
PYBIND11_MODULE(libcsg, m) {
2016-11-04 16:14:52 +05:00
ExportCSG(m);
2014-08-30 06:15:59 +06:00
}
2014-08-31 19:05:24 +06:00
#endif
2014-08-30 06:15:59 +06:00