mirror of
https://github.com/NGSolve/netgen.git
synced 2024-11-11 16:49:16 +05:00
little cleanup and modernization in geom2d code
This commit is contained in:
parent
eb77dd284d
commit
e2df8a5abc
@ -308,55 +308,50 @@ namespace netgen
|
||||
|
||||
void SplineGeometry2d :: CopyEdgeMesh (int from, int to, Mesh & mesh, Point3dTree & searchtree)
|
||||
{
|
||||
// const int D = 2;
|
||||
|
||||
NgArray<int, PointIndex::BASE> mappoints (mesh.GetNP());
|
||||
NgArray<double, PointIndex::BASE> param (mesh.GetNP());
|
||||
mappoints = -1;
|
||||
Array<PointIndex, PointIndex> mappoints (mesh.GetNP());
|
||||
Array<double, PointIndex> param (mesh.GetNP());
|
||||
mappoints = PointIndex::INVALID;
|
||||
param = 0;
|
||||
|
||||
Point3d pmin, pmax;
|
||||
mesh.GetBox (pmin, pmax);
|
||||
double diam2 = Dist2(pmin, pmax);
|
||||
|
||||
if (printmessage_importance>0)
|
||||
cout << "copy edge, from = " << from << " to " << to << endl;
|
||||
PrintMessage(3, string("Copy edge, from ") + ToString(from) + " to " + ToString(to));
|
||||
|
||||
for (int i = 1; i <= mesh.GetNSeg(); i++)
|
||||
for (const auto& seg : mesh.LineSegments())
|
||||
{
|
||||
const Segment & seg = mesh.LineSegment(i);
|
||||
if (seg.edgenr == from)
|
||||
{
|
||||
mappoints.Elem(seg[0]) = 1;
|
||||
param.Elem(seg[0]) = seg.epgeominfo[0].dist;
|
||||
mappoints[seg[0]] = 1;
|
||||
param[seg[0]] = seg.epgeominfo[0].dist;
|
||||
|
||||
mappoints.Elem(seg[1]) = 1;
|
||||
param.Elem(seg[1]) = seg.epgeominfo[1].dist;
|
||||
mappoints[seg[1]] = 1;
|
||||
param[seg[1]] = seg.epgeominfo[1].dist;
|
||||
}
|
||||
}
|
||||
|
||||
bool mapped = false;
|
||||
for (int i = 1; i <= mappoints.Size(); i++)
|
||||
for (auto i : Range(mappoints))
|
||||
{
|
||||
if (mappoints.Get(i) != -1)
|
||||
if (mappoints[i].IsValid())
|
||||
{
|
||||
Point<2> newp = splines.Get(to)->GetPoint (param.Get(i));
|
||||
Point<2> newp = splines.Get(to)->GetPoint (param[i]);
|
||||
Point<3> newp3 (newp(0), newp(1), 0);
|
||||
|
||||
int npi = -1;
|
||||
PointIndex npi = PointIndex::INVALID;
|
||||
|
||||
for (PointIndex pi = PointIndex::BASE;
|
||||
pi < mesh.GetNP()+PointIndex::BASE; pi++)
|
||||
for (auto pi : Range(mesh.Points()))
|
||||
if (Dist2 (mesh.Point(pi), newp3) < 1e-12 * diam2)
|
||||
npi = pi;
|
||||
|
||||
if (npi == -1)
|
||||
if (!npi.IsValid())
|
||||
{
|
||||
npi = mesh.AddPoint (newp3);
|
||||
searchtree.Insert (newp3, npi);
|
||||
}
|
||||
|
||||
mappoints.Elem(i) = npi;
|
||||
mappoints[i] = npi;
|
||||
|
||||
mesh.GetIdentifications().Add (i, npi, to);
|
||||
mapped = true;
|
||||
@ -375,15 +370,15 @@ namespace netgen
|
||||
Segment nseg;
|
||||
nseg.edgenr = to;
|
||||
nseg.si = GetSpline(to-1).bc; // splines.Get(to)->bc;
|
||||
nseg[0] = mappoints.Get(seg[0]);
|
||||
nseg[1] = mappoints.Get(seg[1]);
|
||||
nseg[0] = mappoints[seg[0]];
|
||||
nseg[1] = mappoints[seg[1]];
|
||||
nseg.domin = GetSpline(to-1).leftdom;
|
||||
nseg.domout = GetSpline(to-1).rightdom;
|
||||
|
||||
nseg.epgeominfo[0].edgenr = to;
|
||||
nseg.epgeominfo[0].dist = param.Get(seg[0]);
|
||||
nseg.epgeominfo[0].dist = param[seg[0]];
|
||||
nseg.epgeominfo[1].edgenr = to;
|
||||
nseg.epgeominfo[1].dist = param.Get(seg[1]);
|
||||
nseg.epgeominfo[1].dist = param[seg[1]];
|
||||
mesh.AddSegment (nseg);
|
||||
}
|
||||
}
|
||||
|
@ -62,33 +62,27 @@ DLL_HEADER void ExportGeom2d(py::module &m)
|
||||
}),
|
||||
py::arg("x"), py::arg("y"), py::arg("maxh") = 1e99, py::arg("hpref")=0, py::arg("name")="")
|
||||
.def("Append", FunctionPointer([](SplineGeometry2d &self, py::list segment, int leftdomain, int rightdomain,
|
||||
py::object bc, py::object copy, double maxh, double hpref)
|
||||
optional<variant<int, string>> bc, optional<int> copy, double maxh,
|
||||
double hpref)
|
||||
{
|
||||
py::extract<std::string> segtype(segment[0]);
|
||||
auto segtype = py::cast<std::string>(segment[0]);
|
||||
|
||||
SplineSegExt * seg;
|
||||
if (segtype().compare("line") == 0)
|
||||
if (segtype == "line")
|
||||
{
|
||||
py::extract<int> point_index1(segment[1]);
|
||||
py::extract<int> point_index2(segment[2]);
|
||||
//point_index1.check()
|
||||
|
||||
LineSeg<2> * l = new LineSeg<2>(self.GetPoint(point_index1()), self.GetPoint(point_index2()));
|
||||
LineSeg<2> * l = new LineSeg<2>(self.GetPoint(py::cast<int>(segment[1])),
|
||||
self.GetPoint(py::cast<int>(segment[2])));
|
||||
seg = new SplineSegExt(*l);
|
||||
}
|
||||
else if (segtype().compare("spline3") == 0)
|
||||
else if (segtype == "spline3")
|
||||
{
|
||||
py::extract<int> point_index1(segment[1]);
|
||||
py::extract<int> point_index2(segment[2]);
|
||||
py::extract<int> point_index3(segment[3]);
|
||||
|
||||
SplineSeg3<2> * seg3 = new SplineSeg3<2>(self.GetPoint(point_index1()), self.GetPoint(point_index2()), self.GetPoint(point_index3()));
|
||||
SplineSeg3<2> * seg3 = new SplineSeg3<2>(self.GetPoint(py::cast<int>(segment[1])),
|
||||
self.GetPoint(py::cast<int>(segment[2])),
|
||||
self.GetPoint(py::cast<int>(segment[3])));
|
||||
seg = new SplineSegExt(*seg3);
|
||||
}
|
||||
else
|
||||
{
|
||||
cout << "Appended segment is not a line or a spline3" << endl;
|
||||
}
|
||||
throw Exception("Appended segment is not a line or a spline3");
|
||||
seg->leftdom = leftdomain;
|
||||
seg->rightdom = rightdomain;
|
||||
seg->hmax = maxh;
|
||||
@ -96,23 +90,27 @@ DLL_HEADER void ExportGeom2d(py::module &m)
|
||||
seg->hpref_right = hpref;
|
||||
seg->reffak = 1;
|
||||
seg->copyfrom = -1;
|
||||
if (py::extract<int>(copy).check())
|
||||
seg->copyfrom = py::extract<int>(copy)()+1;
|
||||
|
||||
if (py::extract<int>(bc).check())
|
||||
seg->bc = py::extract<int>(bc)();
|
||||
else if (py::extract<string>(bc).check())
|
||||
if (copy.has_value())
|
||||
seg->copyfrom = *copy+1;
|
||||
|
||||
if (bc.has_value())
|
||||
{
|
||||
string bcname = py::extract<string>(bc)();
|
||||
seg->bc = self.GetNSplines()+1;
|
||||
self.SetBCName(seg->bc, bcname);
|
||||
if(auto intptr = get_if<int>(&*bc); intptr)
|
||||
seg->bc = *intptr;
|
||||
else
|
||||
{
|
||||
auto bcname = get<string>(*bc);
|
||||
seg->bc = self.GetNSplines() + 1;
|
||||
self.SetBCName(seg->bc, bcname);
|
||||
}
|
||||
}
|
||||
else
|
||||
seg->bc = self.GetNSplines()+1;
|
||||
self.AppendSegment(seg);
|
||||
return self.GetNSplines()-1;
|
||||
}), py::arg("point_indices"), py::arg("leftdomain") = 1, py::arg("rightdomain") = py::int_(0),
|
||||
py::arg("bc")=NGDummyArgument(), py::arg("copy")=NGDummyArgument(), py::arg("maxh")=1e99, py::arg("hpref")=0)
|
||||
py::arg("bc")=nullopt, py::arg("copy")=nullopt, py::arg("maxh")=1e99,
|
||||
py::arg("hpref")=0)
|
||||
|
||||
|
||||
.def("AppendSegment", FunctionPointer([](SplineGeometry2d &self, py::list point_indices, int leftdomain, int rightdomain)
|
||||
@ -132,6 +130,8 @@ DLL_HEADER void ExportGeom2d(py::module &m)
|
||||
seg = new SplineSegExt(*seg3);
|
||||
|
||||
}
|
||||
else
|
||||
throw Exception("Can only append segments with 2 or 3 points!");
|
||||
seg->leftdom = leftdomain;
|
||||
seg->rightdom = rightdomain;
|
||||
seg->hmax = 1e99;
|
||||
|
Loading…
Reference in New Issue
Block a user