mirror of
https://github.com/NGSolve/netgen.git
synced 2025-01-23 19:30:33 +05:00
CSG2d interface (Solid2d ctor, EdgeInfo)
This commit is contained in:
parent
b9487cc07a
commit
ceb57a7c5c
@ -1322,6 +1322,34 @@ Solid2d ClipSolids ( Solid2d s1, Solid2d s2, bool intersect)
|
|||||||
return res;
|
return res;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
Solid2d :: Solid2d(const Array<std::variant<Point<2>, EdgeInfo>> & points, string name_)
|
||||||
|
: name(name_)
|
||||||
|
{
|
||||||
|
Loop l;
|
||||||
|
for (auto & v : points)
|
||||||
|
{
|
||||||
|
if(std::holds_alternative<Point<2>>(v))
|
||||||
|
{
|
||||||
|
l.Append(std::get<0>(v), true);
|
||||||
|
}
|
||||||
|
if(std::holds_alternative<EdgeInfo>(v))
|
||||||
|
{
|
||||||
|
l.first->prev->info.Assign( std::get<1>(v) );
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
for(auto v : l.Vertices(ALL))
|
||||||
|
{
|
||||||
|
v->bc = v->info.bc;
|
||||||
|
if(v->info.control_point)
|
||||||
|
{
|
||||||
|
v->spline = Spline(*v, *v->info.control_point, *v->next);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
polys.Append(l);
|
||||||
|
}
|
||||||
|
|
||||||
Solid2d Solid2d :: operator+(Solid2d & other)
|
Solid2d Solid2d :: operator+(Solid2d & other)
|
||||||
{
|
{
|
||||||
if(polys.Size()==0)
|
if(polys.Size()==0)
|
||||||
|
@ -1,6 +1,8 @@
|
|||||||
#ifndef NETGEN_CSG2D_HPP_INCLUDED
|
#ifndef NETGEN_CSG2D_HPP_INCLUDED
|
||||||
#define NETGEN_CSG2D_HPP_INCLUDED
|
#define NETGEN_CSG2D_HPP_INCLUDED
|
||||||
|
|
||||||
|
#include <variant>
|
||||||
|
|
||||||
#include "geometry2d.hpp"
|
#include "geometry2d.hpp"
|
||||||
|
|
||||||
namespace netgen
|
namespace netgen
|
||||||
@ -58,6 +60,34 @@ enum IteratorType
|
|||||||
ALL
|
ALL
|
||||||
};
|
};
|
||||||
|
|
||||||
|
inline constexpr const double MAXH_DEFAULT{1e99};
|
||||||
|
inline const string BC_DEFAULT{""};
|
||||||
|
inline const string MAT_DEFAULT{""};
|
||||||
|
|
||||||
|
struct EdgeInfo
|
||||||
|
{
|
||||||
|
optional<Point<2>> control_point = nullopt; // for spline segments
|
||||||
|
double maxh = MAXH_DEFAULT;
|
||||||
|
string bc = BC_DEFAULT;
|
||||||
|
|
||||||
|
EdgeInfo() = default;
|
||||||
|
EdgeInfo(Point<2> p) : control_point(p) {}
|
||||||
|
EdgeInfo(double h) : maxh(h) {}
|
||||||
|
EdgeInfo(string s) : bc(s) {}
|
||||||
|
EdgeInfo(optional<Point<2>> p, double h, string s)
|
||||||
|
: control_point(p), maxh(h), bc(s)
|
||||||
|
{}
|
||||||
|
|
||||||
|
void Assign( EdgeInfo other )
|
||||||
|
{
|
||||||
|
if(other.control_point != nullopt)
|
||||||
|
control_point = other.control_point;
|
||||||
|
if(other.bc != BC_DEFAULT)
|
||||||
|
bc = other.bc;
|
||||||
|
if(other.maxh != MAXH_DEFAULT)
|
||||||
|
maxh = other.maxh;
|
||||||
|
}
|
||||||
|
};
|
||||||
|
|
||||||
using Spline = SplineSeg3<2>;
|
using Spline = SplineSeg3<2>;
|
||||||
|
|
||||||
@ -80,6 +110,7 @@ struct Vertex : Point<2>
|
|||||||
|
|
||||||
// In case the edge this - next is curved, store the spline information here
|
// In case the edge this - next is curved, store the spline information here
|
||||||
optional<Spline> spline = nullopt;
|
optional<Spline> spline = nullopt;
|
||||||
|
EdgeInfo info;
|
||||||
|
|
||||||
Vertex * Insert(Point<2> p, double lam = -1.0);
|
Vertex * Insert(Point<2> p, double lam = -1.0);
|
||||||
|
|
||||||
@ -538,6 +569,7 @@ struct Solid2d
|
|||||||
|
|
||||||
Solid2d() = default;
|
Solid2d() = default;
|
||||||
Solid2d(string name_) : name(name_) {}
|
Solid2d(string name_) : name(name_) {}
|
||||||
|
Solid2d(const Array<std::variant<Point<2>, EdgeInfo>> & points, string name_);
|
||||||
|
|
||||||
Solid2d operator+(Solid2d & other);
|
Solid2d operator+(Solid2d & other);
|
||||||
Solid2d operator*(Solid2d & other);
|
Solid2d operator*(Solid2d & other);
|
||||||
|
@ -399,7 +399,7 @@ DLL_HEADER void ExportGeom2d(py::module &m)
|
|||||||
|
|
||||||
py::class_<Solid2d>(m, "Solid2d")
|
py::class_<Solid2d>(m, "Solid2d")
|
||||||
.def(py::init<>())
|
.def(py::init<>())
|
||||||
.def(py::init<std::string>())
|
.def(py::init<Array<std::variant<Point<2>, EdgeInfo>>, std::string>(), py::arg("points"), py::arg("mat")=MAT_DEFAULT)
|
||||||
.def_readwrite("name", &Solid2d::name)
|
.def_readwrite("name", &Solid2d::name)
|
||||||
.def("__mul__", [](Solid2d & self, Solid2d & other) { return self*other; })
|
.def("__mul__", [](Solid2d & self, Solid2d & other) { return self*other; })
|
||||||
.def("__add__", [](Solid2d & self, Solid2d & other) { return self+other; })
|
.def("__add__", [](Solid2d & self, Solid2d & other) { return self+other; })
|
||||||
@ -421,6 +421,13 @@ DLL_HEADER void ExportGeom2d(py::module &m)
|
|||||||
.def("Add", &CSG2d::Add)
|
.def("Add", &CSG2d::Add)
|
||||||
;
|
;
|
||||||
|
|
||||||
|
py::class_<EdgeInfo>(m, "EdgeInfo")
|
||||||
|
.def(py::init<>())
|
||||||
|
.def(py::init<const Point<2>&>(), py::arg("control_point"))
|
||||||
|
.def(py::init<double>(), py::arg("maxh"))
|
||||||
|
.def(py::init<string>(), py::arg("bc"))
|
||||||
|
.def(py::init<optional<Point<2>>, double, string>(), py::arg("control_point")=nullopt, py::arg("maxh")=MAXH_DEFAULT, py::arg("bc")=BC_DEFAULT)
|
||||||
|
;
|
||||||
}
|
}
|
||||||
|
|
||||||
PYBIND11_MODULE(libgeom2d, m) {
|
PYBIND11_MODULE(libgeom2d, m) {
|
||||||
|
@ -1,4 +1,4 @@
|
|||||||
from .libngpy._geom2d import SplineGeometry, Solid2d, CSG2d, Rectangle, Circle
|
from .libngpy._geom2d import SplineGeometry, Solid2d, CSG2d, Rectangle, Circle, EdgeInfo
|
||||||
from .meshing import meshsize
|
from .meshing import meshsize
|
||||||
|
|
||||||
unit_square = SplineGeometry()
|
unit_square = SplineGeometry()
|
||||||
@ -137,4 +137,10 @@ SplineGeometry.AddSegment = lambda *args, **kwargs : SplineGeometry.Append(*args
|
|||||||
SplineGeometry.AddPoint = lambda *args, **kwargs : SplineGeometry.AppendPoint(*args, **kwargs)
|
SplineGeometry.AddPoint = lambda *args, **kwargs : SplineGeometry.AppendPoint(*args, **kwargs)
|
||||||
SplineGeometry.CreatePML = CreatePML
|
SplineGeometry.CreatePML = CreatePML
|
||||||
|
|
||||||
|
bc = lambda s : EdgeInfo(bc=s)
|
||||||
|
maxh = lambda h : EdgeInfo(maxh=h)
|
||||||
|
def cp(p_or_px, py_or_none = None):
|
||||||
|
if py_or_none is None:
|
||||||
|
return EdgeInfo(control_point=p)
|
||||||
|
else:
|
||||||
|
return EdgeInfo(control_point=(p_or_px,py_or_none))
|
||||||
|
Loading…
Reference in New Issue
Block a user