mirror of
https://github.com/NGSolve/netgen.git
synced 2024-11-11 16:49:16 +05:00
csg2d - no bc in vertex, handle maxh
This commit is contained in:
parent
671566ef31
commit
b14178b352
@ -98,7 +98,8 @@ Vertex * Vertex :: Insert(Point<2> p, double lam)
|
||||
current = current->next;
|
||||
|
||||
auto pre = current->prev;
|
||||
vnew->bc = pre->bc;
|
||||
if(lam > -1.0)
|
||||
vnew->info = pre->info;
|
||||
|
||||
pre->next = vnew.get();
|
||||
vnew->prev = pre;
|
||||
@ -1148,7 +1149,7 @@ void CreateResult(Solid2d & sp, Solid2d & sr, bool UNION)
|
||||
auto & vnew = R.AppendVertex(*V);
|
||||
if ((status == EXIT) ^ UNION)
|
||||
{
|
||||
vnew.bc = V->bc;
|
||||
vnew.info = V->info;
|
||||
if(V->spline)
|
||||
vnew.spline = *V->spline;
|
||||
else
|
||||
@ -1166,7 +1167,7 @@ void CreateResult(Solid2d & sp, Solid2d & sr, bool UNION)
|
||||
}
|
||||
else
|
||||
vnew.spline = nullopt;
|
||||
vnew.bc = V->bc;
|
||||
vnew.info = V->info;
|
||||
V->is_intersection = false; // mark visited vertices
|
||||
}
|
||||
if(V == I)
|
||||
@ -1338,11 +1339,8 @@ Solid2d :: Solid2d(const Array<std::variant<Point<2>, EdgeInfo>> & points, strin
|
||||
if(v->info.bc==BC_DEFAULT)
|
||||
v->info.bc = bc;
|
||||
|
||||
v->bc = v->info.bc;
|
||||
if(v->info.control_point)
|
||||
{
|
||||
v->spline = Spline(*v, *v->info.control_point, *v->next);
|
||||
}
|
||||
}
|
||||
|
||||
polys.Append(l);
|
||||
@ -1472,6 +1470,7 @@ shared_ptr<netgen::SplineGeometry2d> CSG2d :: GenerateSplineGeometry()
|
||||
int bc;
|
||||
int p2;
|
||||
double weight;
|
||||
double maxh = 1e99;
|
||||
};
|
||||
|
||||
auto geo = std::make_shared<netgen::SplineGeometry2d>();
|
||||
@ -1575,11 +1574,13 @@ shared_ptr<netgen::SplineGeometry2d> CSG2d :: GenerateSplineGeometry()
|
||||
ls.p2 = pi2;
|
||||
ls.weight = weight;
|
||||
|
||||
if(bcmap.count(p0.bc)==0)
|
||||
bcmap[p0.bc] = bcmap.size()+1;
|
||||
if(bcmap.count(p0.info.bc)==0)
|
||||
bcmap[p0.info.bc] = bcmap.size()+1;
|
||||
|
||||
if(ls.bc==0 || p0.bc != BC_DEFAULT)
|
||||
ls.bc = bcmap[p0.bc];
|
||||
if(ls.bc==0 || p0.info.bc != BC_DEFAULT)
|
||||
ls.bc = bcmap[p0.info.bc];
|
||||
|
||||
ls.maxh = min(ls.maxh, p0.info.maxh);
|
||||
|
||||
if(li!=ri)
|
||||
{
|
||||
@ -1619,7 +1620,7 @@ shared_ptr<netgen::SplineGeometry2d> CSG2d :: GenerateSplineGeometry()
|
||||
seg->bc = ls.bc;
|
||||
seg->reffak = 1;
|
||||
seg->copyfrom = -1;
|
||||
seg->hmax = 1e99;
|
||||
seg->hmax = ls.maxh;
|
||||
geo->AppendSegment(seg);
|
||||
}
|
||||
return geo;
|
||||
|
@ -88,7 +88,7 @@ struct EdgeInfo
|
||||
if(other.bc != BC_DEFAULT)
|
||||
bc = other.bc;
|
||||
if(other.maxh != MAXH_DEFAULT)
|
||||
maxh = other.maxh;
|
||||
maxh = min(maxh, other.maxh);
|
||||
}
|
||||
};
|
||||
|
||||
@ -107,8 +107,6 @@ struct Vertex : Point<2>
|
||||
IntersectionLabel label = NONE; // type of intersection vertex
|
||||
EntryExitLabel enex = NEITHER; // entry/exit "flag"
|
||||
|
||||
string bc = "";
|
||||
|
||||
// In case the edge this - next is curved, store the spline information here
|
||||
optional<Spline> spline = nullopt;
|
||||
EdgeInfo info;
|
||||
@ -413,7 +411,7 @@ struct Loop
|
||||
Vertex & AppendVertex(const Vertex & v)
|
||||
{
|
||||
auto & vnew = Append( static_cast<Point<2>>(v), true );
|
||||
vnew.bc = v.bc;
|
||||
vnew.info = v.info;
|
||||
if(v.spline)
|
||||
vnew.spline = *v.spline;
|
||||
return vnew;
|
||||
@ -545,7 +543,7 @@ struct Loop
|
||||
void SetBC(string bc)
|
||||
{
|
||||
for(auto v : Vertices(ALL))
|
||||
v->bc = bc;
|
||||
v->info.bc = bc;
|
||||
}
|
||||
|
||||
size_t Size() const
|
||||
@ -617,6 +615,28 @@ struct Solid2d
|
||||
{
|
||||
return RotateRad( ang/180.*M_PI );
|
||||
}
|
||||
|
||||
Solid2d & BC(string bc)
|
||||
{
|
||||
for(auto & p : polys)
|
||||
for(auto v : p.Vertices(ALL))
|
||||
v->info.bc = bc;
|
||||
return *this;
|
||||
}
|
||||
|
||||
Solid2d & Maxh(double maxh)
|
||||
{
|
||||
for(auto & p : polys)
|
||||
for(auto v : p.Vertices(ALL))
|
||||
v->info.maxh = maxh;
|
||||
return *this;
|
||||
}
|
||||
|
||||
Solid2d & Mat(string mat)
|
||||
{
|
||||
name = mat;
|
||||
return *this;
|
||||
}
|
||||
};
|
||||
|
||||
|
||||
|
@ -400,10 +400,18 @@ DLL_HEADER void ExportGeom2d(py::module &m)
|
||||
py::class_<Solid2d>(m, "Solid2d")
|
||||
.def(py::init<>())
|
||||
.def(py::init<Array<std::variant<Point<2>, EdgeInfo>>, std::string, std::string>(), py::arg("points"), py::arg("mat")=MAT_DEFAULT, py::arg("bc")=BC_DEFAULT)
|
||||
.def_readwrite("name", &Solid2d::name)
|
||||
.def("__mul__", [](Solid2d & self, Solid2d & other) { return self*other; })
|
||||
.def("__add__", [](Solid2d & self, Solid2d & other) { return self+other; })
|
||||
.def("__sub__", [](Solid2d & self, Solid2d & other) { return self-other; })
|
||||
|
||||
.def(py::self+py::self)
|
||||
.def(py::self-py::self)
|
||||
.def(py::self*py::self)
|
||||
.def(py::self+=py::self)
|
||||
.def(py::self-=py::self)
|
||||
.def(py::self*=py::self)
|
||||
|
||||
.def("Mat", &Solid2d::Mat)
|
||||
.def("BC", &Solid2d::BC)
|
||||
.def("Maxh", &Solid2d::Maxh)
|
||||
|
||||
.def("Copy", [](Solid2d & self) -> Solid2d { return self; })
|
||||
.def("Move", &Solid2d::Move)
|
||||
.def("Scale", &Solid2d::Scale)
|
||||
|
Loading…
Reference in New Issue
Block a user