mirror of
https://github.com/NGSolve/netgen.git
synced 2024-12-25 05:20:34 +05:00
Merge branch 'occ_geo_maxh' into 'master'
add maxh property to occ shapes, add TopoDS_Edge.Split method See merge request jschoeberl/netgen!402
This commit is contained in:
commit
97c4fa724b
@ -225,7 +225,6 @@ namespace netgen
|
||||
const MeshingParameters & mparam)
|
||||
{
|
||||
double s0, s1;
|
||||
double maxh = mparam.maxh;
|
||||
int nsubedges = 1;
|
||||
gp_Pnt pnt, oldpnt;
|
||||
double svalue[DIVIDEEDGESECTIONS];
|
||||
@ -871,7 +870,7 @@ namespace netgen
|
||||
|
||||
|
||||
// Philippose - 15/01/2009
|
||||
double maxh = geom.face_maxh[k-1];
|
||||
double maxh = min2(geom.face_maxh[k-1], OCCGeometry::global_shape_properties[TopoDS::Face(geom.fmap(k)).TShape()].maxh);
|
||||
//double maxh = mparam.maxh;
|
||||
// int noldpoints = mesh->GetNP();
|
||||
int noldsurfel = mesh.GetNSE();
|
||||
@ -1151,13 +1150,14 @@ namespace netgen
|
||||
int face_index = geom.fmap.FindIndex(parent_face);
|
||||
|
||||
if(face_index >= 1) localh = min(localh,geom.face_maxh[face_index - 1]);
|
||||
localh = min2(localh, OCCGeometry::global_shape_properties[parent_face.TShape()].maxh);
|
||||
}
|
||||
|
||||
Handle(Geom_Curve) c = BRep_Tool::Curve(e, s0, s1);
|
||||
|
||||
localh = min2(localh, OCCGeometry::global_shape_properties[e.TShape()].maxh);
|
||||
maxedgelen = max (maxedgelen, len);
|
||||
minedgelen = min (minedgelen, len);
|
||||
|
||||
int maxj = max((int) ceil(len/localh), 2);
|
||||
|
||||
for (int j = 0; j <= maxj; j++)
|
||||
|
@ -214,10 +214,12 @@ namespace netgen
|
||||
public:
|
||||
optional<string> name;
|
||||
optional<Vec<3>> col;
|
||||
double maxh = 1e99;
|
||||
void Merge(const ShapeProperties & prop2)
|
||||
{
|
||||
if (prop2.name) name = prop2.name;
|
||||
if (prop2.col) col = prop2.col;
|
||||
maxh = min2(maxh, prop2.maxh);
|
||||
}
|
||||
};
|
||||
|
||||
|
@ -53,6 +53,9 @@
|
||||
#include <ShapeUpgrade_UnifySameDomain.hxx>
|
||||
#include <GeomLProp_SLProps.hxx>
|
||||
|
||||
#include <BOPTools_AlgoTools.hxx>
|
||||
#include <IntTools_Context.hxx>
|
||||
|
||||
#if OCC_VERSION_MAJOR>=7 && OCC_VERSION_MINOR>=4
|
||||
#define OCC_HAVE_DUMP_JSON
|
||||
#endif
|
||||
@ -679,7 +682,15 @@ DLL_HEADER void ExportNgOCCShapes(py::module &m)
|
||||
}, [](const TopoDS_Shape & self, string name) {
|
||||
OCCGeometry::global_shape_properties[self.TShape()].name = name;
|
||||
})
|
||||
|
||||
.def_property("maxh",
|
||||
[](const TopoDS_Shape& self)
|
||||
{
|
||||
return OCCGeometry::global_shape_properties[self.TShape()].maxh;
|
||||
},
|
||||
[](TopoDS_Shape& self, double val)
|
||||
{
|
||||
OCCGeometry::global_shape_properties[self.TShape()].maxh = val;
|
||||
})
|
||||
.def_property("col", [](const TopoDS_Shape & self) {
|
||||
auto it = OCCGeometry::global_shape_properties.find(self.TShape());
|
||||
Vec<3> col(0.2, 0.2, 0.2);
|
||||
@ -1082,6 +1093,42 @@ DLL_HEADER void ExportNgOCCShapes(py::module &m)
|
||||
curve->D1(s1, p, v);
|
||||
return v;
|
||||
})
|
||||
.def("Split", [](const TopoDS_Edge& self, py::args args)
|
||||
{
|
||||
ListOfShapes new_edges;
|
||||
double s0, s1;
|
||||
auto curve = BRep_Tool::Curve(self, s0, s1);
|
||||
double tstart, t, dist;
|
||||
TopoDS_Vertex vstart, vend;
|
||||
vstart = TopExp::FirstVertex(self);
|
||||
IntTools_Context context;
|
||||
tstart = s0;
|
||||
for(auto arg : args)
|
||||
{
|
||||
if(py::isinstance<py::float_>(arg))
|
||||
t = s0 + py::cast<double>(arg) * (s1-s0);
|
||||
else
|
||||
{
|
||||
auto p = py::cast<gp_Pnt>(arg);
|
||||
auto result = context.ComputePE(p, 0., self, t, dist);
|
||||
if(result != 0)
|
||||
throw Exception("Error in finding splitting points on edge!");
|
||||
}
|
||||
auto p = curve->Value(t);
|
||||
vend = BRepBuilderAPI_MakeVertex(p);
|
||||
auto newE = TopoDS::Edge(self.EmptyCopied());
|
||||
BOPTools_AlgoTools::MakeSplitEdge(self, vstart, tstart, vend, t, newE);
|
||||
new_edges.push_back(newE);
|
||||
vstart = vend;
|
||||
tstart = t;
|
||||
}
|
||||
auto newE = TopoDS::Edge(self.EmptyCopied());
|
||||
t = s1;
|
||||
vend = TopExp::LastVertex(self);
|
||||
BOPTools_AlgoTools::MakeSplitEdge(self, vstart, tstart, vend, t, newE);
|
||||
new_edges.push_back(newE);
|
||||
return new_edges;
|
||||
}, "Splits edge at given parameters. Parameters can either be floating values in (0,1), then edge parametrization is used. Or it can be points, then the projection of these points are used for splitting the edge.")
|
||||
;
|
||||
py::class_<TopoDS_Wire, TopoDS_Shape> (m, "Wire")
|
||||
.def(py::init([](const TopoDS_Edge & edge) {
|
||||
|
Loading…
Reference in New Issue
Block a user