Merge branch 'curved_splinesurface' into 'master'

Curved splinesurface

See merge request jschoeberl/netgen!110
This commit is contained in:
Joachim Schöberl 2018-12-30 14:33:55 +00:00
commit 459f168126
4 changed files with 47 additions and 6 deletions

View File

@ -1236,7 +1236,15 @@ namespace netgen
return max2(bb/(aa*aa),aa/(bb*bb)); return max2(bb/(aa*aa),aa/(bb*bb));
} }
int EllipticCylinder :: IsIdentic(const Surface& s2, int& inv, double eps) const
{
const EllipticCylinder* ps2 = dynamic_cast<const EllipticCylinder*>(&s2);
if (!ps2) return 0;
if((vl - ps2->vl).Length() > eps || (vs - ps2->vs).Length() > eps || (a-ps2->a).Length() > eps)
return 0;
return 1;
}
Point<3> EllipticCylinder :: GetSurfacePoint () const Point<3> EllipticCylinder :: GetSurfacePoint () const
{ {

View File

@ -303,6 +303,7 @@ namespace netgen
const Box<3> & bbox, const Box<3> & bbox,
double facets) const; double facets) const;
virtual int IsIdentic (const Surface & s2, int & inv, double eps) const;
virtual double MaxCurvature () const; virtual double MaxCurvature () const;

View File

@ -235,13 +235,17 @@ DLL_HEADER void ExportCSG(py::module &m)
return self.GetNP()-1; return self.GetNP()-1;
}), }),
py::arg("x"),py::arg("y"),py::arg("z"),py::arg("hpref")=false) py::arg("x"),py::arg("y"),py::arg("z"),py::arg("hpref")=false)
.def("AddSegment", FunctionPointer .def("AddSegment", [] (SplineSurface & self, int i1, int i2, string bcname, double maxh)
([] (SplineSurface & self, int i1, int i2, string bcname, double maxh)
{ {
auto seg = make_shared<LineSeg<3>>(self.GetPoint(i1),self.GetPoint(i2)); auto seg = make_shared<LineSeg<3>>(self.GetPoint(i1),self.GetPoint(i2));
self.AppendSegment(seg,bcname,maxh); self.AppendSegment(seg,bcname,maxh);
}), },
py::arg("pnt1"),py::arg("pnt2"),py::arg("bcname")="default", py::arg("maxh")=-1.) 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.)
; ;
py::class_<SPSolid, shared_ptr<SPSolid>> (m, "Solid") py::class_<SPSolid, shared_ptr<SPSolid>> (m, "Solid")
@ -471,7 +475,12 @@ However, when r = 0, the top part becomes a point(tip) and meshing fails!
self.GetTopLevelObject(tlonr) -> SetBCProp(surf->GetBase()->GetBCProperty()); self.GetTopLevelObject(tlonr) -> SetBCProp(surf->GetBase()->GetBCProperty());
self.GetTopLevelObject(tlonr) -> SetBCName(surf->GetBase()->GetBCName()); self.GetTopLevelObject(tlonr) -> SetBCName(surf->GetBase()->GetBCName());
self.GetTopLevelObject(tlonr) -> SetMaxH(surf->GetBase()->GetMaxH()); self.GetTopLevelObject(tlonr) -> SetMaxH(surf->GetBase()->GetMaxH());
for(auto p : surf->GetPoints()) Array<Point<3>> non_midpoints;
for(auto spline : surf->GetSplines())
{
non_midpoints.Append(spline->GetPoint(0));
}
for(auto p : non_midpoints)
self.AddUserPoint(p); self.AddUserPoint(p);
self.AddSplineSurface(surf); self.AddSplineSurface(surf);
}), }),

View File

@ -61,9 +61,32 @@ void SplineSurface :: AppendPoint(const Point<3> & p, const double reffac, const
} }
cuttings->Append(plane); cuttings->Append(plane);
} }
else
{
auto spline3 = dynamic_cast<SplineSeg3<3>*>(spline.get());
if(spline3)
{
auto p1 = Point<3>(spline3->StartPI());
Project(p1);
auto p2 = Point<3>(spline3->TangentPoint());
Project(p2);
auto p3 = Point<3>(spline3->EndPI());
Project(p3);
Vec<3> v1 = p2-p1;
Vec<3> v2 = p2-p3;
Point<3> mid = p1 - v2;
cout << "mid point = " << mid << endl;
cout << "v1 = " << v1 << endl;
cout << "v2 = " << v2 << endl;
auto cyl = make_shared<EllipticCylinder>(mid, v1, v2);
if(maxh[i] > 0)
cyl->SetMaxH(maxh[i]);
cuttings->Append(cyl);
}
else else
throw NgException("Spline type not implemented for SplineSurface!"); throw NgException("Spline type not implemented for SplineSurface!");
} }
}
all_cuts = cuttings; all_cuts = cuttings;
return cuttings; return cuttings;
} }