mirror of
https://github.com/NGSolve/netgen.git
synced 2024-12-25 05:20:34 +05:00
features for 2d meshing
This commit is contained in:
parent
5bcf28d196
commit
7165c90fb1
@ -47,7 +47,7 @@ DLL_HEADER void ExportGeom2d()
|
|||||||
return self.geompoints.Size()-1;
|
return self.geompoints.Size()-1;
|
||||||
}),
|
}),
|
||||||
(bp::arg("self"), bp::arg("x"), bp::arg("y"), bp::arg("maxh") = 1e99))
|
(bp::arg("self"), bp::arg("x"), bp::arg("y"), bp::arg("maxh") = 1e99))
|
||||||
.def("Append", FunctionPointer([](SplineGeometry2d &self, bp::list segment, int leftdomain, int rightdomain, bp::object bc)
|
.def("Append", FunctionPointer([](SplineGeometry2d &self, bp::list segment, int leftdomain, int rightdomain, bp::object bc, double maxh)
|
||||||
{
|
{
|
||||||
bp::extract<std::string> segtype(segment[0]);
|
bp::extract<std::string> segtype(segment[0]);
|
||||||
|
|
||||||
@ -76,7 +76,7 @@ DLL_HEADER void ExportGeom2d()
|
|||||||
}
|
}
|
||||||
seg->leftdom = leftdomain;
|
seg->leftdom = leftdomain;
|
||||||
seg->rightdom = rightdomain;
|
seg->rightdom = rightdomain;
|
||||||
seg->hmax = 1e99;
|
seg->hmax = maxh;
|
||||||
seg->reffak = 1;
|
seg->reffak = 1;
|
||||||
seg->copyfrom = -1;
|
seg->copyfrom = -1;
|
||||||
if (bp::extract<int>(bc).check())
|
if (bp::extract<int>(bc).check())
|
||||||
@ -93,7 +93,8 @@ DLL_HEADER void ExportGeom2d()
|
|||||||
seg->bc = self.GetNSplines()+1;
|
seg->bc = self.GetNSplines()+1;
|
||||||
self.AppendSegment(seg);
|
self.AppendSegment(seg);
|
||||||
}), (bp::arg("self"), bp::arg("point_indices"), bp::arg("leftdomain") = 1, bp::arg("rightdomain") = 0,
|
}), (bp::arg("self"), bp::arg("point_indices"), bp::arg("leftdomain") = 1, bp::arg("rightdomain") = 0,
|
||||||
bp::arg("bc")=bp::object()))
|
bp::arg("bc")=bp::object(), bp::arg("maxh")=1e99
|
||||||
|
))
|
||||||
|
|
||||||
|
|
||||||
.def("AppendSegment", FunctionPointer([](SplineGeometry2d &self, bp::list point_indices, int leftdomain, int rightdomain)
|
.def("AppendSegment", FunctionPointer([](SplineGeometry2d &self, bp::list point_indices, int leftdomain, int rightdomain)
|
||||||
|
@ -20,6 +20,30 @@ pnums = [unit_square.AppendPoint(*p) for p in pnts]
|
|||||||
for l1,l2,bc in lines:
|
for l1,l2,bc in lines:
|
||||||
unit_square.Append( ["line", pnums[l1], pnums[l2]], bc=bc)
|
unit_square.Append( ["line", pnums[l1], pnums[l2]], bc=bc)
|
||||||
|
|
||||||
|
|
||||||
|
def MakeRectangle (geo, p1, p2, **args):
|
||||||
|
p1x, p1y = p1
|
||||||
|
p2x, p2y = p2
|
||||||
|
p1x,p2x = min(p1x,p2x), max(p1x, p2x)
|
||||||
|
p1y,p2y = min(p1y,p2y), max(p1y, p2y)
|
||||||
|
pts = [geo.AppendPoint(*p) for p in [(p1x,p1y), (p2x, p1y), (p2x, p2y), (p1x, p2y)]]
|
||||||
|
for p1,p2 in [(0,1), (1, 2), (2, 3), (3, 0)]:
|
||||||
|
geo.Append( ["line", pts[p1], pts[p2]], **args)
|
||||||
|
|
||||||
|
def MakeCircle (geo, c, r, **args):
|
||||||
|
cx,cy = c
|
||||||
|
pts = [geo.AppendPoint(*p) for p in [(cx,cy-r), (cx+r,cy-r), (cx+r,cy), (cx+r,cy+r), \
|
||||||
|
(cx,cy+r), (cx-r,cy+r), (cx-r,cy), (cx-r,cy-r)]]
|
||||||
|
for p1,p2,p3 in [(0,1,2), (2,3,4), (4, 5, 6), (6, 7, 0)]:
|
||||||
|
geo.Append( ["spline3", pts[p1], pts[p2], pts[p3]], **args)
|
||||||
|
|
||||||
|
|
||||||
|
SplineGeometry.AddCircle = lambda geo, c, r, **args : MakeCircle(geo, c, r, **args)
|
||||||
|
SplineGeometry.AddRectangle = lambda geo, p1, p2, **args : MakeRectangle(geo, p1, p2, **args)
|
||||||
|
SplineGeometry.AddSegment = SplineGeometry.Append
|
||||||
|
SplineGeometry.AddPoint = SplineGeometry.AppendPoint
|
||||||
|
|
||||||
|
|
||||||
__all__ = ['SplineGeometry', 'unit_square']
|
__all__ = ['SplineGeometry', 'unit_square']
|
||||||
|
|
||||||
|
|
||||||
|
Loading…
Reference in New Issue
Block a user