mirror of
https://git.dev.opencascade.org/repos/occt.git
synced 2025-04-05 18:16:23 +03:00
0026305: BRepFeat_MakePrism returns inconsistent results && 026315: BRepFeat_MakeRevol fails to create revol from shape
Test-cases for issues #26305 and #26315
This commit is contained in:
parent
c0e32b3c3d
commit
d7988ee19f
@ -54,6 +54,7 @@
|
|||||||
#include <LocOpe_CSIntersector.hxx>
|
#include <LocOpe_CSIntersector.hxx>
|
||||||
#include <LocOpe_PntFace.hxx>
|
#include <LocOpe_PntFace.hxx>
|
||||||
#include <LocOpe_BuildShape.hxx>
|
#include <LocOpe_BuildShape.hxx>
|
||||||
|
#include <ElSLib.hxx>
|
||||||
|
|
||||||
#include <TColGeom_SequenceOfCurve.hxx>
|
#include <TColGeom_SequenceOfCurve.hxx>
|
||||||
|
|
||||||
@ -490,14 +491,9 @@ void BRepFeat::FaceUntil(const TopoDS_Shape& Sbase,
|
|||||||
{
|
{
|
||||||
Bnd_Box B;
|
Bnd_Box B;
|
||||||
BRepBndLib::Add(Sbase,B);
|
BRepBndLib::Add(Sbase,B);
|
||||||
Standard_Real c[6], bnd;
|
Standard_Real x[2], y[2], z[2];
|
||||||
B.Get(c[0],c[2],c[4],c[1],c[3],c[5]);
|
B.Get(x[0],y[0],z[0],x[1],y[1],z[1]);
|
||||||
bnd = c[0];
|
Standard_Real diam = 10.*Sqrt(B.SquareExtent());
|
||||||
for(Standard_Integer i = 1 ; i < 6; i++) {
|
|
||||||
if(c[i] > bnd) bnd = c[i];
|
|
||||||
}
|
|
||||||
bnd = 10*bnd;
|
|
||||||
|
|
||||||
|
|
||||||
Handle(Geom_Surface) s = BRep_Tool::Surface(FUntil);
|
Handle(Geom_Surface) s = BRep_Tool::Surface(FUntil);
|
||||||
Handle(Standard_Type) styp = s->DynamicType();
|
Handle(Standard_Type) styp = s->DynamicType();
|
||||||
@ -507,16 +503,80 @@ void BRepFeat::FaceUntil(const TopoDS_Shape& Sbase,
|
|||||||
}
|
}
|
||||||
Handle(Geom_RectangularTrimmedSurface) str;
|
Handle(Geom_RectangularTrimmedSurface) str;
|
||||||
if (styp == STANDARD_TYPE(Geom_Plane)) {
|
if (styp == STANDARD_TYPE(Geom_Plane)) {
|
||||||
|
gp_Pln aPln = Handle(Geom_Plane)::DownCast(s)->Pln();
|
||||||
|
Standard_Real u, v, umin = RealLast(), umax = -umin,
|
||||||
|
vmin = RealLast(), vmax = -vmin;
|
||||||
|
for(Standard_Integer i = 0 ; i < 2; i++)
|
||||||
|
{
|
||||||
|
for(Standard_Integer j = 0; j < 2; j++)
|
||||||
|
{
|
||||||
|
for(Standard_Integer k = 0; k < 2; k++)
|
||||||
|
{
|
||||||
|
gp_Pnt aP(x[i], y[j], z[k]);
|
||||||
|
ElSLib::Parameters(aPln, aP, u, v);
|
||||||
|
if(u < umin)
|
||||||
|
umin = u;
|
||||||
|
if(u > umax)
|
||||||
|
umax = u;
|
||||||
|
if(v < vmin)
|
||||||
|
vmin = v;
|
||||||
|
if(v > vmax)
|
||||||
|
vmax = v;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
}
|
||||||
|
umin -= diam;
|
||||||
|
umax += diam;
|
||||||
|
vmin -= diam;
|
||||||
|
vmax += diam;
|
||||||
str = new Geom_RectangularTrimmedSurface
|
str = new Geom_RectangularTrimmedSurface
|
||||||
(s, bnd, -bnd, bnd, -bnd, Standard_True, Standard_True);
|
(s, umin, umax, vmin, vmax, Standard_True, Standard_True);
|
||||||
}
|
}
|
||||||
else if (styp == STANDARD_TYPE(Geom_CylindricalSurface)) {
|
else if (styp == STANDARD_TYPE(Geom_CylindricalSurface)) {
|
||||||
|
gp_Cylinder aCyl = Handle(Geom_CylindricalSurface)::DownCast(s)->Cylinder();
|
||||||
|
Standard_Real u, v, vmin = RealLast(), vmax = -vmin;
|
||||||
|
for(Standard_Integer i = 0 ; i < 2; i++)
|
||||||
|
{
|
||||||
|
for(Standard_Integer j = 0; j < 2; j++)
|
||||||
|
{
|
||||||
|
for(Standard_Integer k = 0; k < 2; k++)
|
||||||
|
{
|
||||||
|
gp_Pnt aP(x[i], y[j], z[k]);
|
||||||
|
ElSLib::Parameters(aCyl, aP, u, v);
|
||||||
|
if(v < vmin)
|
||||||
|
vmin = v;
|
||||||
|
if(v > vmax)
|
||||||
|
vmax = v;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
}
|
||||||
|
vmin -= diam;
|
||||||
|
vmax += diam;
|
||||||
str = new Geom_RectangularTrimmedSurface
|
str = new Geom_RectangularTrimmedSurface
|
||||||
(s, bnd, -bnd, Standard_False, Standard_True);
|
(s, vmin, vmax, Standard_False, Standard_True);
|
||||||
}
|
}
|
||||||
else if (styp == STANDARD_TYPE(Geom_ConicalSurface)) {
|
else if (styp == STANDARD_TYPE(Geom_ConicalSurface)) {
|
||||||
|
gp_Cone aCon = Handle(Geom_ConicalSurface)::DownCast(s)->Cone();
|
||||||
|
Standard_Real u, v, vmin = RealLast(), vmax = -vmin;
|
||||||
|
for(Standard_Integer i = 0 ; i < 2; i++)
|
||||||
|
{
|
||||||
|
for(Standard_Integer j = 0; j < 2; j++)
|
||||||
|
{
|
||||||
|
for(Standard_Integer k = 0; k < 2; k++)
|
||||||
|
{
|
||||||
|
gp_Pnt aP(x[i], y[j], z[k]);
|
||||||
|
ElSLib::Parameters(aCon, aP, u, v);
|
||||||
|
if(v < vmin)
|
||||||
|
vmin = v;
|
||||||
|
if(v > vmax)
|
||||||
|
vmax = v;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
}
|
||||||
|
vmin -= diam;
|
||||||
|
vmax += diam;
|
||||||
str = new Geom_RectangularTrimmedSurface
|
str = new Geom_RectangularTrimmedSurface
|
||||||
(s, bnd, -bnd, Standard_False, Standard_True);
|
(s, vmin, vmax, Standard_False, Standard_True);
|
||||||
}
|
}
|
||||||
else {
|
else {
|
||||||
FUntil.Nullify();
|
FUntil.Nullify();
|
||||||
|
@ -260,15 +260,15 @@ void BRepFeat_MakePrism::Perform(const Standard_Real Length)
|
|||||||
if(myLShape.ShapeType() == TopAbs_WIRE) {
|
if(myLShape.ShapeType() == TopAbs_WIRE) {
|
||||||
TopExp_Explorer ex1(VraiPrism, TopAbs_FACE);
|
TopExp_Explorer ex1(VraiPrism, TopAbs_FACE);
|
||||||
for(; ex1.More(); ex1.Next()) {
|
for(; ex1.More(); ex1.Next()) {
|
||||||
TopExp_Explorer ex2(ex1.Current(), TopAbs_WIRE);
|
TopExp_Explorer ex2(ex1.Current(), TopAbs_WIRE);
|
||||||
for(; ex2.More(); ex2.Next()) {
|
for(; ex2.More(); ex2.Next()) {
|
||||||
if(ex2.Current().IsSame(myLShape)) {
|
if(ex2.Current().IsSame(myLShape)) {
|
||||||
FFace = TopoDS::Face(ex1.Current());
|
FFace = TopoDS::Face(ex1.Current());
|
||||||
found = Standard_True;
|
found = Standard_True;
|
||||||
break;
|
break;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
if(found) break;
|
if(found) break;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
@ -276,9 +276,9 @@ void BRepFeat_MakePrism::Perform(const Standard_Real Length)
|
|||||||
for(; exp.More(); exp.Next()) {
|
for(; exp.More(); exp.Next()) {
|
||||||
const TopoDS_Face& ff = TopoDS::Face(exp.Current());
|
const TopoDS_Face& ff = TopoDS::Face(exp.Current());
|
||||||
if(ToFuse(ff, FFace)) {
|
if(ToFuse(ff, FFace)) {
|
||||||
TopTools_DataMapOfShapeListOfShape sl;
|
TopTools_DataMapOfShapeListOfShape sl;
|
||||||
if(!FFace.IsSame(myPbase) && BRepFeat::IsInside(ff, FFace))
|
if(!FFace.IsSame(myPbase) && BRepFeat::IsInside(ff, FFace))
|
||||||
break;
|
break;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
@ -570,15 +570,17 @@ void BRepFeat_MakePrism::Perform(const TopoDS_Shape& From,
|
|||||||
ASI2.Perform(scur);
|
ASI2.Perform(scur);
|
||||||
TopAbs_Orientation OrU, OrF;
|
TopAbs_Orientation OrU, OrF;
|
||||||
TopoDS_Face FFrom, FUntil;
|
TopoDS_Face FFrom, FUntil;
|
||||||
|
Standard_Real ParF, ParU;
|
||||||
if (ASI1.IsDone() && ASI1.NbPoints(1) >=1) {
|
if (ASI1.IsDone() && ASI1.NbPoints(1) >=1) {
|
||||||
if (myFuse == 1) {
|
if (myFuse == 1) {
|
||||||
OrU = ASI1.Point(1,1).Orientation();
|
OrU = ASI1.Point(1,1).Orientation();
|
||||||
}
|
}
|
||||||
else {
|
else {
|
||||||
OrU = ASI1.Point(1,ASI1.NbPoints(1)).Orientation();
|
OrU = ASI1.Point(1,ASI1.NbPoints(1)).Orientation();
|
||||||
}
|
}
|
||||||
if(sens==-1) OrU = TopAbs::Reverse(OrU);
|
if(sens==-1) OrU = TopAbs::Reverse(OrU);
|
||||||
FUntil = ASI1.Point(1,1).Face();
|
FUntil = ASI1.Point(1,1).Face();
|
||||||
|
ParU = ASI1.Point(1,1).Parameter();
|
||||||
}
|
}
|
||||||
else {
|
else {
|
||||||
NotDone();
|
NotDone();
|
||||||
@ -589,12 +591,20 @@ void BRepFeat_MakePrism::Perform(const TopoDS_Shape& From,
|
|||||||
OrF = ASI2.Point(1,1).Orientation();
|
OrF = ASI2.Point(1,1).Orientation();
|
||||||
if(sens==1) OrF = TopAbs::Reverse(OrF);
|
if(sens==1) OrF = TopAbs::Reverse(OrF);
|
||||||
FFrom = ASI2.Point(1,1).Face();
|
FFrom = ASI2.Point(1,1).Face();
|
||||||
|
ParF = ASI2.Point(1,1).Parameter();
|
||||||
}
|
}
|
||||||
else {
|
else {
|
||||||
NotDone();
|
NotDone();
|
||||||
myStatusError = BRepFeat_NoIntersectF;
|
myStatusError = BRepFeat_NoIntersectF;
|
||||||
return;
|
return;
|
||||||
}
|
}
|
||||||
|
if(tran > 0 && (Abs(ParU) < Abs(ParF)))
|
||||||
|
{
|
||||||
|
TopAbs_Orientation Or;
|
||||||
|
Or = OrU;
|
||||||
|
OrU = OrF;
|
||||||
|
OrF = Or;
|
||||||
|
}
|
||||||
TopoDS_Shape Comp;
|
TopoDS_Shape Comp;
|
||||||
BRep_Builder B;
|
BRep_Builder B;
|
||||||
B.MakeCompound(TopoDS::Compound(Comp));
|
B.MakeCompound(TopoDS::Compound(Comp));
|
||||||
|
@ -78,22 +78,22 @@ extern Standard_Boolean BRepFeat_GettraceFEAT();
|
|||||||
#endif
|
#endif
|
||||||
|
|
||||||
static void MajMap(const TopoDS_Shape&, // base
|
static void MajMap(const TopoDS_Shape&, // base
|
||||||
const LocOpe_Revol&,
|
const LocOpe_Revol&,
|
||||||
TopTools_DataMapOfShapeListOfShape&, // myMap
|
TopTools_DataMapOfShapeListOfShape&, // myMap
|
||||||
TopoDS_Shape&, // myFShape
|
TopoDS_Shape&, // myFShape
|
||||||
TopoDS_Shape&); // myLShape
|
TopoDS_Shape&); // myLShape
|
||||||
|
|
||||||
|
|
||||||
static void VerifGluedFaces(const TopoDS_Face& theSkface,
|
static void VerifGluedFaces(const TopoDS_Face& theSkface,
|
||||||
const TopoDS_Shape& thePbase,
|
const TopoDS_Shape& thePbase,
|
||||||
Handle(Geom_Curve)& theBCurve,
|
Handle(Geom_Curve)& theBCurve,
|
||||||
TColGeom_SequenceOfCurve& theCurves,
|
TColGeom_SequenceOfCurve& theCurves,
|
||||||
LocOpe_Revol& theRevol,
|
LocOpe_Revol& theRevol,
|
||||||
TopTools_DataMapOfShapeShape& theMap);
|
TopTools_DataMapOfShapeShape& theMap);
|
||||||
|
|
||||||
|
|
||||||
static Standard_Boolean ToFuse(const TopoDS_Face& ,
|
static Standard_Boolean ToFuse(const TopoDS_Face& ,
|
||||||
const TopoDS_Face&);
|
const TopoDS_Face&);
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
@ -104,11 +104,11 @@ static Standard_Boolean ToFuse(const TopoDS_Face& ,
|
|||||||
//=======================================================================
|
//=======================================================================
|
||||||
|
|
||||||
void BRepFeat_MakeRevol::Init(const TopoDS_Shape& Sbase,
|
void BRepFeat_MakeRevol::Init(const TopoDS_Shape& Sbase,
|
||||||
const TopoDS_Shape& Pbase,
|
const TopoDS_Shape& Pbase,
|
||||||
const TopoDS_Face& Skface,
|
const TopoDS_Face& Skface,
|
||||||
const gp_Ax1& Axis,
|
const gp_Ax1& Axis,
|
||||||
const Standard_Integer Mode,
|
const Standard_Integer Mode,
|
||||||
const Standard_Boolean Modify)
|
const Standard_Boolean Modify)
|
||||||
{
|
{
|
||||||
#ifdef OCCT_DEBUG
|
#ifdef OCCT_DEBUG
|
||||||
Standard_Boolean trc = BRepFeat_GettraceFEAT();
|
Standard_Boolean trc = BRepFeat_GettraceFEAT();
|
||||||
@ -140,7 +140,7 @@ void BRepFeat_MakeRevol::Init(const TopoDS_Shape& Sbase,
|
|||||||
myJustGluer = Standard_False;
|
myJustGluer = Standard_False;
|
||||||
|
|
||||||
//-------------- ifv
|
//-------------- ifv
|
||||||
// mySkface.Nullify();
|
// mySkface.Nullify();
|
||||||
//-------------- ifv
|
//-------------- ifv
|
||||||
|
|
||||||
|
|
||||||
@ -171,7 +171,7 @@ void BRepFeat_MakeRevol::Init(const TopoDS_Shape& Sbase,
|
|||||||
//=======================================================================
|
//=======================================================================
|
||||||
|
|
||||||
void BRepFeat_MakeRevol::Add(const TopoDS_Edge& E,
|
void BRepFeat_MakeRevol::Add(const TopoDS_Edge& E,
|
||||||
const TopoDS_Face& F)
|
const TopoDS_Face& F)
|
||||||
{
|
{
|
||||||
#ifdef OCCT_DEBUG
|
#ifdef OCCT_DEBUG
|
||||||
Standard_Boolean trc = BRepFeat_GettraceFEAT();
|
Standard_Boolean trc = BRepFeat_GettraceFEAT();
|
||||||
@ -237,13 +237,13 @@ void BRepFeat_MakeRevol::Perform(const Standard_Real Angle)
|
|||||||
if(RevolComp) {
|
if(RevolComp) {
|
||||||
/*
|
/*
|
||||||
if (!mySkface.IsNull() || !mySlface.IsEmpty()) {
|
if (!mySkface.IsNull() || !mySlface.IsEmpty()) {
|
||||||
for (exp.Init(mySbase,TopAbs_FACE); exp.More(); exp.Next()) {
|
for (exp.Init(mySbase,TopAbs_FACE); exp.More(); exp.Next()) {
|
||||||
if (exp.Current().IsSame(mySkface)) {
|
if (exp.Current().IsSame(mySkface)) {
|
||||||
angledec = M_PI/5; // pourquoi pas
|
angledec = M_PI/5; // pourquoi pas
|
||||||
if (myFuse) angledec = -angledec;
|
if (myFuse) angledec = -angledec;
|
||||||
break;
|
break;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
*/
|
*/
|
||||||
mySkface.Nullify();
|
mySkface.Nullify();
|
||||||
@ -252,7 +252,7 @@ void BRepFeat_MakeRevol::Perform(const Standard_Real Angle)
|
|||||||
else theRevol.Perform(myPbase, myAxis, Angle, angledec);
|
else theRevol.Perform(myPbase, myAxis, Angle, angledec);
|
||||||
|
|
||||||
TopoDS_Shape VraiRevol = theRevol.Shape();
|
TopoDS_Shape VraiRevol = theRevol.Shape();
|
||||||
|
|
||||||
MajMap(myPbase,theRevol,myMap,myFShape,myLShape);
|
MajMap(myPbase,theRevol,myMap,myFShape,myLShape);
|
||||||
|
|
||||||
myGShape = VraiRevol;
|
myGShape = VraiRevol;
|
||||||
@ -268,32 +268,32 @@ void BRepFeat_MakeRevol::Perform(const Standard_Real Angle)
|
|||||||
}
|
}
|
||||||
|
|
||||||
TopoDS_Face FFace;
|
TopoDS_Face FFace;
|
||||||
|
|
||||||
Standard_Boolean found = Standard_False;
|
Standard_Boolean found = Standard_False;
|
||||||
|
|
||||||
if(!mySkface.IsNull() || !mySlface.IsEmpty()) {
|
if(!mySkface.IsNull() || !mySlface.IsEmpty()) {
|
||||||
if(myLShape.ShapeType() == TopAbs_WIRE) {
|
if(myLShape.ShapeType() == TopAbs_WIRE) {
|
||||||
TopExp_Explorer ex1(VraiRevol, TopAbs_FACE);
|
TopExp_Explorer ex1(VraiRevol, TopAbs_FACE);
|
||||||
for(; ex1.More(); ex1.Next()) {
|
for(; ex1.More(); ex1.Next()) {
|
||||||
TopExp_Explorer ex2(ex1.Current(), TopAbs_WIRE);
|
TopExp_Explorer ex2(ex1.Current(), TopAbs_WIRE);
|
||||||
for(; ex2.More(); ex2.Next()) {
|
for(; ex2.More(); ex2.Next()) {
|
||||||
if(ex2.Current().IsSame(myLShape)) {
|
if(ex2.Current().IsSame(myLShape)) {
|
||||||
FFace = TopoDS::Face(ex1.Current());
|
FFace = TopoDS::Face(ex1.Current());
|
||||||
found = Standard_True;
|
found = Standard_True;
|
||||||
break;
|
break;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
if(found) break;
|
if(found) break;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
TopExp_Explorer exp(mySbase, TopAbs_FACE);
|
TopExp_Explorer exp(mySbase, TopAbs_FACE);
|
||||||
for(; exp.More(); exp.Next()) {
|
for(; exp.More(); exp.Next()) {
|
||||||
const TopoDS_Face& ff = TopoDS::Face(exp.Current());
|
const TopoDS_Face& ff = TopoDS::Face(exp.Current());
|
||||||
if(ToFuse(ff, FFace)) {
|
if(ToFuse(ff, FFace)) {
|
||||||
TopTools_DataMapOfShapeListOfShape sl;
|
TopTools_DataMapOfShapeListOfShape sl;
|
||||||
if(!FFace.IsSame(myPbase) && BRepFeat::IsInside(ff, FFace))
|
if(!FFace.IsSame(myPbase) && BRepFeat::IsInside(ff, FFace))
|
||||||
break;
|
break;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
@ -372,8 +372,8 @@ void BRepFeat_MakeRevol::Perform(const TopoDS_Shape& Until)
|
|||||||
Standard_Boolean Trf = TransformShapeFU(1);
|
Standard_Boolean Trf = TransformShapeFU(1);
|
||||||
ShapeUntilValid();
|
ShapeUntilValid();
|
||||||
|
|
||||||
// Do systematically almost complete revolution
|
// Do systematically almost complete revolution
|
||||||
// BRepSweep_Revol theRevol(myPbase,myAxis,2.*M_PI-10.*Precision::Angular());
|
// BRepSweep_Revol theRevol(myPbase,myAxis,2.*M_PI-10.*Precision::Angular());
|
||||||
LocOpe_Revol theRevol;
|
LocOpe_Revol theRevol;
|
||||||
if(!TourComplet) {
|
if(!TourComplet) {
|
||||||
Angle = 2.*M_PI- 3*M_PI/180.;
|
Angle = 2.*M_PI- 3*M_PI/180.;
|
||||||
@ -387,10 +387,10 @@ void BRepFeat_MakeRevol::Perform(const TopoDS_Shape& Until)
|
|||||||
|
|
||||||
|
|
||||||
if(!Trf) {
|
if(!Trf) {
|
||||||
|
|
||||||
myGShape = VraiRevol;
|
myGShape = VraiRevol;
|
||||||
GeneratedShapeValid();
|
GeneratedShapeValid();
|
||||||
|
|
||||||
TopoDS_Shape Base = theRevol.FirstShape();
|
TopoDS_Shape Base = theRevol.FirstShape();
|
||||||
exp.Init(Base, TopAbs_FACE);
|
exp.Init(Base, TopAbs_FACE);
|
||||||
TopoDS_Face theBase = TopoDS::Face(exp.Current());
|
TopoDS_Face theBase = TopoDS::Face(exp.Current());
|
||||||
@ -430,79 +430,79 @@ void BRepFeat_MakeRevol::Perform(const TopoDS_Shape& Until)
|
|||||||
TopoDS_Shape Cutsh = trP.Shape();
|
TopoDS_Shape Cutsh = trP.Shape();
|
||||||
TopExp_Explorer ex(Cutsh, TopAbs_SOLID);
|
TopExp_Explorer ex(Cutsh, TopAbs_SOLID);
|
||||||
for(; ex.More(); ex.Next()) {
|
for(; ex.More(); ex.Next()) {
|
||||||
TopExp_Explorer ex1(ex.Current(), TopAbs_FACE);
|
TopExp_Explorer ex1(ex.Current(), TopAbs_FACE);
|
||||||
for(; ex1.More(); ex1.Next()) {
|
for(; ex1.More(); ex1.Next()) {
|
||||||
const TopoDS_Face& fac = TopoDS::Face(ex1.Current());
|
const TopoDS_Face& fac = TopoDS::Face(ex1.Current());
|
||||||
if(fac.IsSame(myPbase)) {
|
if(fac.IsSame(myPbase)) {
|
||||||
VraiRevol = ex.Current();
|
VraiRevol = ex.Current();
|
||||||
break;
|
break;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
if(myFuse == 1) {
|
if(myFuse == 1) {
|
||||||
//modified by NIZNHY-PKV Thu Mar 21 18:17:53 2002 f
|
//modified by NIZNHY-PKV Thu Mar 21 18:17:53 2002 f
|
||||||
//BRepAlgo_Fuse f(mySbase, VraiRevol);
|
//BRepAlgo_Fuse f(mySbase, VraiRevol);
|
||||||
//myShape = f.Shape();
|
//myShape = f.Shape();
|
||||||
//UpdateDescendants(f.Builder(), myShape, Standard_False);
|
//UpdateDescendants(f.Builder(), myShape, Standard_False);
|
||||||
BRepAlgoAPI_Fuse f(mySbase, VraiRevol);
|
BRepAlgoAPI_Fuse f(mySbase, VraiRevol);
|
||||||
myShape = f.Shape();
|
myShape = f.Shape();
|
||||||
UpdateDescendants(f, myShape, Standard_False);
|
UpdateDescendants(f, myShape, Standard_False);
|
||||||
//modified by NIZNHY-PKV Thu Mar 21 18:17:57 2002 t
|
//modified by NIZNHY-PKV Thu Mar 21 18:17:57 2002 t
|
||||||
Done();
|
Done();
|
||||||
}
|
}
|
||||||
else if(myFuse == 0) {
|
else if(myFuse == 0) {
|
||||||
//modified by NIZNHY-PKV Thu Mar 21 18:18:23 2002 f
|
//modified by NIZNHY-PKV Thu Mar 21 18:18:23 2002 f
|
||||||
//BRepAlgo_Cut c(mySbase, VraiRevol);
|
//BRepAlgo_Cut c(mySbase, VraiRevol);
|
||||||
//myShape = c.Shape();
|
//myShape = c.Shape();
|
||||||
//UpdateDescendants(c.Builder(), myShape, Standard_False);
|
//UpdateDescendants(c.Builder(), myShape, Standard_False);
|
||||||
BRepAlgoAPI_Cut c(mySbase, VraiRevol);
|
BRepAlgoAPI_Cut c(mySbase, VraiRevol);
|
||||||
myShape = c.Shape();
|
myShape = c.Shape();
|
||||||
UpdateDescendants(c, myShape, Standard_False);
|
UpdateDescendants(c, myShape, Standard_False);
|
||||||
//modified by NIZNHY-PKV Thu Mar 21 18:18:28 2002 t
|
//modified by NIZNHY-PKV Thu Mar 21 18:18:28 2002 t
|
||||||
Done();
|
Done();
|
||||||
}
|
}
|
||||||
else {
|
else {
|
||||||
myShape = VraiRevol;
|
myShape = VraiRevol;
|
||||||
Done();
|
Done();
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
// Loop of control of descendance
|
// Loop of control of descendance
|
||||||
/*
|
/*
|
||||||
TopExp_Explorer expr(mySbase, TopAbs_FACE);
|
TopExp_Explorer expr(mySbase, TopAbs_FACE);
|
||||||
char nom1[20], nom2[20];
|
char nom1[20], nom2[20];
|
||||||
Standard_Integer ii = 0;
|
Standard_Integer ii = 0;
|
||||||
for(; expr.More(); expr.Next()) {
|
for(; expr.More(); expr.Next()) {
|
||||||
ii++;
|
ii++;
|
||||||
sprintf(nom1, "faceinitial_%d", ii);
|
sprintf(nom1, "faceinitial_%d", ii);
|
||||||
DBRep::Set(nom1, expr.Current());
|
DBRep::Set(nom1, expr.Current());
|
||||||
Standard_Integer jj = 0;
|
Standard_Integer jj = 0;
|
||||||
const TopTools_ListOfShape& list = Modified(expr.Current());
|
const TopTools_ListOfShape& list = Modified(expr.Current());
|
||||||
TopTools_ListIteratorOfListOfShape ite(list);
|
TopTools_ListIteratorOfListOfShape ite(list);
|
||||||
for(; ite.More(); ite.Next()) {
|
for(; ite.More(); ite.Next()) {
|
||||||
jj++;
|
jj++;
|
||||||
sprintf(nom2, "facemodifie_%d_%d", ii, jj);
|
sprintf(nom2, "facemodifie_%d_%d", ii, jj);
|
||||||
DBRep::Set(nom2, ite.Value());
|
DBRep::Set(nom2, ite.Value());
|
||||||
}
|
|
||||||
}
|
}
|
||||||
|
}
|
||||||
|
|
||||||
expr.Init(myPbase, TopAbs_EDGE);
|
expr.Init(myPbase, TopAbs_EDGE);
|
||||||
ii=0;
|
ii=0;
|
||||||
for(; expr.More(); expr.Next()) {
|
for(; expr.More(); expr.Next()) {
|
||||||
ii++;
|
ii++;
|
||||||
sprintf(nom1, "edgeinitial_%d", ii);
|
sprintf(nom1, "edgeinitial_%d", ii);
|
||||||
DBRep::Set(nom1, expr.Current());
|
DBRep::Set(nom1, expr.Current());
|
||||||
Standard_Integer jj = 0;
|
Standard_Integer jj = 0;
|
||||||
const TopTools_ListOfShape& list = Generated(expr.Current());
|
const TopTools_ListOfShape& list = Generated(expr.Current());
|
||||||
TopTools_ListIteratorOfListOfShape ite(list);
|
TopTools_ListIteratorOfListOfShape ite(list);
|
||||||
for(; ite.More(); ite.Next()) {
|
for(; ite.More(); ite.Next()) {
|
||||||
jj++;
|
jj++;
|
||||||
sprintf(nom2, "facegeneree_%d_%d", ii, jj);
|
sprintf(nom2, "facegeneree_%d_%d", ii, jj);
|
||||||
DBRep::Set(nom2, ite.Value());
|
DBRep::Set(nom2, ite.Value());
|
||||||
}
|
|
||||||
}
|
}
|
||||||
*/
|
}
|
||||||
}
|
*/
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
//=======================================================================
|
//=======================================================================
|
||||||
@ -511,7 +511,7 @@ void BRepFeat_MakeRevol::Perform(const TopoDS_Shape& Until)
|
|||||||
//=======================================================================
|
//=======================================================================
|
||||||
|
|
||||||
void BRepFeat_MakeRevol::Perform(const TopoDS_Shape& From,
|
void BRepFeat_MakeRevol::Perform(const TopoDS_Shape& From,
|
||||||
const TopoDS_Shape& Until)
|
const TopoDS_Shape& Until)
|
||||||
{
|
{
|
||||||
#ifdef OCCT_DEBUG
|
#ifdef OCCT_DEBUG
|
||||||
Standard_Boolean trc = BRepFeat_GettraceFEAT();
|
Standard_Boolean trc = BRepFeat_GettraceFEAT();
|
||||||
@ -553,7 +553,7 @@ void BRepFeat_MakeRevol::Perform(const TopoDS_Shape& From,
|
|||||||
mySUntil = Until;
|
mySUntil = Until;
|
||||||
Standard_Boolean Trfu = TransformShapeFU(1);
|
Standard_Boolean Trfu = TransformShapeFU(1);
|
||||||
ShapeUntilValid();
|
ShapeUntilValid();
|
||||||
|
|
||||||
if(Trfu != Trff) {
|
if(Trfu != Trff) {
|
||||||
NotDone();
|
NotDone();
|
||||||
myStatusError = BRepFeat_IncTypes;
|
myStatusError = BRepFeat_IncTypes;
|
||||||
@ -570,7 +570,7 @@ void BRepFeat_MakeRevol::Perform(const TopoDS_Shape& From,
|
|||||||
myGShape = VraiRevol;
|
myGShape = VraiRevol;
|
||||||
GeneratedShapeValid();
|
GeneratedShapeValid();
|
||||||
GluedFacesValid();
|
GluedFacesValid();
|
||||||
// VerifGluedFaces(mySkface, theBase, myBCurve, myCurves, theRevol, myGluedF);
|
// VerifGluedFaces(mySkface, theBase, myBCurve, myCurves, theRevol, myGluedF);
|
||||||
|
|
||||||
theRevol.Curves(myCurves);
|
theRevol.Curves(myCurves);
|
||||||
myBCurve = theRevol.BarycCurve();
|
myBCurve = theRevol.BarycCurve();
|
||||||
@ -604,7 +604,8 @@ void BRepFeat_MakeRevol::Perform(const TopoDS_Shape& From,
|
|||||||
pr1 = ElCLib::InPeriod(pr1,PrU-2*M_PI,PrU);
|
pr1 = ElCLib::InPeriod(pr1,PrU-2*M_PI,PrU);
|
||||||
Standard_Real pr2 = ASI2.Point(1,ASI2.NbPoints(1)).Parameter();
|
Standard_Real pr2 = ASI2.Point(1,ASI2.NbPoints(1)).Parameter();
|
||||||
pr2 = ElCLib::InPeriod(pr2,PrU-2*M_PI,PrU);
|
pr2 = ElCLib::InPeriod(pr2,PrU-2*M_PI,PrU);
|
||||||
OrF = OrU;
|
//OrF = OrU;
|
||||||
|
OrF = TopAbs::Reverse(OrU);
|
||||||
FFrom = ASI2.Point(1,1).Face();
|
FFrom = ASI2.Point(1,1).Face();
|
||||||
PrF = Max(pr1, pr2);
|
PrF = Max(pr1, pr2);
|
||||||
}
|
}
|
||||||
@ -631,18 +632,18 @@ void BRepFeat_MakeRevol::Perform(const TopoDS_Shape& From,
|
|||||||
//modified by NIZNHY-PKV Thu Mar 21 18:18:57 2002 t
|
//modified by NIZNHY-PKV Thu Mar 21 18:18:57 2002 t
|
||||||
TopoDS_Shape Cutsh = trP.Shape();
|
TopoDS_Shape Cutsh = trP.Shape();
|
||||||
TopExp_Explorer ex(Cutsh, TopAbs_SOLID);
|
TopExp_Explorer ex(Cutsh, TopAbs_SOLID);
|
||||||
// Standard_Real PrF = BRepFeat::ParametricBarycenter(mySFrom, myBCurve);
|
// Standard_Real PrF = BRepFeat::ParametricBarycenter(mySFrom, myBCurve);
|
||||||
// Standard_Real PrU = BRepFeat::ParametricBarycenter(mySUntil, myBCurve);
|
// Standard_Real PrU = BRepFeat::ParametricBarycenter(mySUntil, myBCurve);
|
||||||
VraiRevol = ex.Current();
|
VraiRevol = ex.Current();
|
||||||
for(; ex.More(); ex.Next()) {
|
for(; ex.More(); ex.Next()) {
|
||||||
Standard_Real PrCur = BRepFeat::
|
Standard_Real PrCur = BRepFeat::
|
||||||
ParametricBarycenter(ex.Current(), myBCurve);
|
ParametricBarycenter(ex.Current(), myBCurve);
|
||||||
if(PrF <= PrCur && PrU >= PrCur) {
|
if(PrF <= PrCur && PrU >= PrCur) {
|
||||||
VraiRevol = ex.Current();
|
VraiRevol = ex.Current();
|
||||||
break;
|
break;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
if(myFuse == 1) {
|
if(myFuse == 1 && !myJustFeat) {
|
||||||
//modified by NIZNHY-PKV Thu Mar 21 18:19:14 2002 f
|
//modified by NIZNHY-PKV Thu Mar 21 18:19:14 2002 f
|
||||||
//BRepAlgo_Fuse f(mySbase, VraiRevol);
|
//BRepAlgo_Fuse f(mySbase, VraiRevol);
|
||||||
//myShape = f.Shape();
|
//myShape = f.Shape();
|
||||||
@ -653,7 +654,7 @@ void BRepFeat_MakeRevol::Perform(const TopoDS_Shape& From,
|
|||||||
//modified by NIZNHY-PKV Thu Mar 21 18:19:18 2002 t
|
//modified by NIZNHY-PKV Thu Mar 21 18:19:18 2002 t
|
||||||
Done();
|
Done();
|
||||||
}
|
}
|
||||||
else if(myFuse == 0) {
|
else if(myFuse == 0 && !myJustFeat) {
|
||||||
//modified by NIZNHY-PKV Thu Mar 21 18:19:46 2002 f
|
//modified by NIZNHY-PKV Thu Mar 21 18:19:46 2002 f
|
||||||
//BRepAlgo_Cut c(mySbase, VraiRevol);
|
//BRepAlgo_Cut c(mySbase, VraiRevol);
|
||||||
//myShape = c.Shape();
|
//myShape = c.Shape();
|
||||||
@ -692,7 +693,7 @@ void BRepFeat_MakeRevol::PerformThruAll()
|
|||||||
//=======================================================================
|
//=======================================================================
|
||||||
|
|
||||||
void BRepFeat_MakeRevol::PerformUntilAngle(const TopoDS_Shape& Until,
|
void BRepFeat_MakeRevol::PerformUntilAngle(const TopoDS_Shape& Until,
|
||||||
const Standard_Real Angle)
|
const Standard_Real Angle)
|
||||||
{
|
{
|
||||||
#ifdef OCCT_DEBUG
|
#ifdef OCCT_DEBUG
|
||||||
Standard_Boolean trc = BRepFeat_GettraceFEAT();
|
Standard_Boolean trc = BRepFeat_GettraceFEAT();
|
||||||
@ -721,8 +722,8 @@ void BRepFeat_MakeRevol::PerformUntilAngle(const TopoDS_Shape& Until,
|
|||||||
Standard_Boolean Trf = TransformShapeFU(1);
|
Standard_Boolean Trf = TransformShapeFU(1);
|
||||||
ShapeUntilValid();
|
ShapeUntilValid();
|
||||||
|
|
||||||
// Produce systematicallt an almost complete revolution
|
// Produce systematicallt an almost complete revolution
|
||||||
// BRepSweep_Revol theRevol(myPbase,myAxis,2.*M_PI-10.*Precision::Angular());
|
// BRepSweep_Revol theRevol(myPbase,myAxis,2.*M_PI-10.*Precision::Angular());
|
||||||
LocOpe_Revol theRevol;
|
LocOpe_Revol theRevol;
|
||||||
theRevol.Perform(myPbase, myAxis, Angle);
|
theRevol.Perform(myPbase, myAxis, Angle);
|
||||||
TopoDS_Shape VraiRevol = theRevol.Shape();
|
TopoDS_Shape VraiRevol = theRevol.Shape();
|
||||||
@ -732,7 +733,7 @@ void BRepFeat_MakeRevol::PerformUntilAngle(const TopoDS_Shape& Until,
|
|||||||
if(Trf) {
|
if(Trf) {
|
||||||
myGShape = VraiRevol;
|
myGShape = VraiRevol;
|
||||||
GeneratedShapeValid();
|
GeneratedShapeValid();
|
||||||
|
|
||||||
TopoDS_Shape Base = theRevol.FirstShape();
|
TopoDS_Shape Base = theRevol.FirstShape();
|
||||||
exp.Init(Base, TopAbs_FACE);
|
exp.Init(Base, TopAbs_FACE);
|
||||||
TopoDS_Face theBase = TopoDS::Face(exp.Current());
|
TopoDS_Face theBase = TopoDS::Face(exp.Current());
|
||||||
@ -773,40 +774,40 @@ void BRepFeat_MakeRevol::PerformUntilAngle(const TopoDS_Shape& Until,
|
|||||||
TopoDS_Shape Cutsh = trP.Shape();
|
TopoDS_Shape Cutsh = trP.Shape();
|
||||||
TopExp_Explorer ex(Cutsh, TopAbs_SOLID);
|
TopExp_Explorer ex(Cutsh, TopAbs_SOLID);
|
||||||
for(; ex.More(); ex.Next()) {
|
for(; ex.More(); ex.Next()) {
|
||||||
TopExp_Explorer ex1(ex.Current(), TopAbs_FACE);
|
TopExp_Explorer ex1(ex.Current(), TopAbs_FACE);
|
||||||
for(; ex1.More(); ex1.Next()) {
|
for(; ex1.More(); ex1.Next()) {
|
||||||
const TopoDS_Face& fac = TopoDS::Face(ex1.Current());
|
const TopoDS_Face& fac = TopoDS::Face(ex1.Current());
|
||||||
if(fac.IsSame(myPbase)) {
|
if(fac.IsSame(myPbase)) {
|
||||||
VraiRevol = ex.Current();
|
VraiRevol = ex.Current();
|
||||||
break;
|
break;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
if(myFuse == 1) {
|
if(myFuse == 1) {
|
||||||
//modified by NIZNHY-PKV Thu Mar 21 18:20:36 2002 f
|
//modified by NIZNHY-PKV Thu Mar 21 18:20:36 2002 f
|
||||||
//BRepAlgo_Fuse f(mySbase, VraiRevol);
|
//BRepAlgo_Fuse f(mySbase, VraiRevol);
|
||||||
//myShape = f.Shape();
|
//myShape = f.Shape();
|
||||||
//UpdateDescendants(f.Builder(), myShape, Standard_False);
|
//UpdateDescendants(f.Builder(), myShape, Standard_False);
|
||||||
BRepAlgoAPI_Fuse f(mySbase, VraiRevol);
|
BRepAlgoAPI_Fuse f(mySbase, VraiRevol);
|
||||||
myShape = f.Shape();
|
myShape = f.Shape();
|
||||||
UpdateDescendants(f, myShape, Standard_False);
|
UpdateDescendants(f, myShape, Standard_False);
|
||||||
//modified by NIZNHY-PKV Thu Mar 21 18:20:40 2002 t
|
//modified by NIZNHY-PKV Thu Mar 21 18:20:40 2002 t
|
||||||
Done();
|
Done();
|
||||||
}
|
}
|
||||||
else if(myFuse == 0) {
|
else if(myFuse == 0) {
|
||||||
//modified by NIZNHY-PKV Thu Mar 21 18:21:07 2002 f
|
//modified by NIZNHY-PKV Thu Mar 21 18:21:07 2002 f
|
||||||
//BRepAlgo_Cut c(mySbase, VraiRevol);
|
//BRepAlgo_Cut c(mySbase, VraiRevol);
|
||||||
//myShape = c.Shape();
|
//myShape = c.Shape();
|
||||||
//UpdateDescendants(c.Builder(), myShape, Standard_False);
|
//UpdateDescendants(c.Builder(), myShape, Standard_False);
|
||||||
BRepAlgoAPI_Cut c(mySbase, VraiRevol);
|
BRepAlgoAPI_Cut c(mySbase, VraiRevol);
|
||||||
myShape = c.Shape();
|
myShape = c.Shape();
|
||||||
UpdateDescendants(c, myShape, Standard_False);
|
UpdateDescendants(c, myShape, Standard_False);
|
||||||
//modified by NIZNHY-PKV Thu Mar 21 18:21:26 2002 t
|
//modified by NIZNHY-PKV Thu Mar 21 18:21:26 2002 t
|
||||||
Done();
|
Done();
|
||||||
}
|
}
|
||||||
else {
|
else {
|
||||||
myShape = VraiRevol;
|
myShape = VraiRevol;
|
||||||
Done();
|
Done();
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
@ -838,15 +839,15 @@ Handle(Geom_Curve) BRepFeat_MakeRevol::BarycCurve()
|
|||||||
//=======================================================================
|
//=======================================================================
|
||||||
|
|
||||||
static void VerifGluedFaces(const TopoDS_Face& theSkface,
|
static void VerifGluedFaces(const TopoDS_Face& theSkface,
|
||||||
const TopoDS_Shape& thePbase,
|
const TopoDS_Shape& thePbase,
|
||||||
Handle(Geom_Curve)& theBCurve,
|
Handle(Geom_Curve)& theBCurve,
|
||||||
TColGeom_SequenceOfCurve& theCurves,
|
TColGeom_SequenceOfCurve& theCurves,
|
||||||
LocOpe_Revol& theRevol,
|
LocOpe_Revol& theRevol,
|
||||||
TopTools_DataMapOfShapeShape& theMap)
|
TopTools_DataMapOfShapeShape& theMap)
|
||||||
{
|
{
|
||||||
Standard_Boolean GluedFaces = Standard_True;
|
Standard_Boolean GluedFaces = Standard_True;
|
||||||
TopoDS_Shape VraiRevol = theRevol.Shape();
|
TopoDS_Shape VraiRevol = theRevol.Shape();
|
||||||
|
|
||||||
TColGeom_SequenceOfCurve scur;
|
TColGeom_SequenceOfCurve scur;
|
||||||
theRevol.Curves(theCurves);
|
theRevol.Curves(theCurves);
|
||||||
theBCurve = theRevol.BarycCurve();
|
theBCurve = theRevol.BarycCurve();
|
||||||
@ -871,13 +872,13 @@ static void VerifGluedFaces(const TopoDS_Face& theSkface,
|
|||||||
for(; ex.More(); ex.Next()) {
|
for(; ex.More(); ex.Next()) {
|
||||||
TopExp_Explorer ex1(ex.Current(), TopAbs_FACE);
|
TopExp_Explorer ex1(ex.Current(), TopAbs_FACE);
|
||||||
for(; ex1.More(); ex1.Next()) {
|
for(; ex1.More(); ex1.Next()) {
|
||||||
const TopoDS_Face& fac1 = TopoDS::Face(ex1.Current());
|
const TopoDS_Face& fac1 = TopoDS::Face(ex1.Current());
|
||||||
TopExp_Explorer ex2(thePbase, TopAbs_FACE);
|
TopExp_Explorer ex2(thePbase, TopAbs_FACE);
|
||||||
for(; ex2.More(); ex2.Next()) {
|
for(; ex2.More(); ex2.Next()) {
|
||||||
const TopoDS_Face& fac2 = TopoDS::Face(ex2.Current());
|
const TopoDS_Face& fac2 = TopoDS::Face(ex2.Current());
|
||||||
if(fac1.IsSame(fac2)) break;
|
if(fac1.IsSame(fac2)) break;
|
||||||
}
|
}
|
||||||
if (ex2.More()) break;
|
if (ex2.More()) break;
|
||||||
}
|
}
|
||||||
if (ex1.More()) continue;
|
if (ex1.More()) continue;
|
||||||
GluedFaces = Standard_False;
|
GluedFaces = Standard_False;
|
||||||
@ -899,10 +900,10 @@ static void VerifGluedFaces(const TopoDS_Face& theSkface,
|
|||||||
//=======================================================================
|
//=======================================================================
|
||||||
|
|
||||||
static void MajMap(const TopoDS_Shape& theB,
|
static void MajMap(const TopoDS_Shape& theB,
|
||||||
const LocOpe_Revol& theP,
|
const LocOpe_Revol& theP,
|
||||||
TopTools_DataMapOfShapeListOfShape& theMap, // myMap
|
TopTools_DataMapOfShapeListOfShape& theMap, // myMap
|
||||||
TopoDS_Shape& theFShape, // myFShape
|
TopoDS_Shape& theFShape, // myFShape
|
||||||
TopoDS_Shape& theLShape) // myLShape
|
TopoDS_Shape& theLShape) // myLShape
|
||||||
{
|
{
|
||||||
TopExp_Explorer exp(theP.FirstShape(),TopAbs_WIRE);
|
TopExp_Explorer exp(theP.FirstShape(),TopAbs_WIRE);
|
||||||
if (exp.More()) {
|
if (exp.More()) {
|
||||||
@ -913,7 +914,7 @@ static void MajMap(const TopoDS_Shape& theB,
|
|||||||
theMap(theFShape).Append(exp.Current());
|
theMap(theFShape).Append(exp.Current());
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
exp.Init(theP.LastShape(),TopAbs_WIRE);
|
exp.Init(theP.LastShape(),TopAbs_WIRE);
|
||||||
if (exp.More()) {
|
if (exp.More()) {
|
||||||
theLShape = exp.Current();
|
theLShape = exp.Current();
|
||||||
@ -941,7 +942,7 @@ static void MajMap(const TopoDS_Shape& theB,
|
|||||||
//=======================================================================
|
//=======================================================================
|
||||||
|
|
||||||
Standard_Boolean ToFuse(const TopoDS_Face& F1,
|
Standard_Boolean ToFuse(const TopoDS_Face& F1,
|
||||||
const TopoDS_Face& F2)
|
const TopoDS_Face& F2)
|
||||||
{
|
{
|
||||||
if (F1.IsNull() || F2.IsNull()) {
|
if (F1.IsNull() || F2.IsNull()) {
|
||||||
return Standard_False;
|
return Standard_False;
|
||||||
|
34
tests/bugs/modalg_6/bug26305_1
Normal file
34
tests/bugs/modalg_6/bug26305_1
Normal file
@ -0,0 +1,34 @@
|
|||||||
|
puts "========"
|
||||||
|
puts "OCC26305"
|
||||||
|
puts "========"
|
||||||
|
puts ""
|
||||||
|
###################################################
|
||||||
|
# BRepFeat_MakePrism returns inconsistent results
|
||||||
|
###################################################
|
||||||
|
|
||||||
|
circle aCircle -9 -9 0 0 0 1 1 0 0 10
|
||||||
|
mkedge anEdge aCircle
|
||||||
|
wire aWire anEdge
|
||||||
|
mkplane aCircle aWire 0
|
||||||
|
|
||||||
|
plane aFromPlane 0 0 -10 0 0 1 1 0 0
|
||||||
|
mkface aFromPlane aFromPlane
|
||||||
|
plane aToPlane 0 0 10 0 0 1 1 0 0
|
||||||
|
mkface aToPlane aToPlane
|
||||||
|
|
||||||
|
featprism aCircle aCircle aCircle 0 0 1 2 1
|
||||||
|
featperform prism aResult aFromPlane aToPlane
|
||||||
|
|
||||||
|
checkshape aResult
|
||||||
|
|
||||||
|
vinit
|
||||||
|
vdisplay aResult
|
||||||
|
vsetdispmode aResult 1
|
||||||
|
vfit
|
||||||
|
|
||||||
|
set bug_info [vreadpixel 300 340 name]
|
||||||
|
if {$bug_info == "BLACK 0"} {
|
||||||
|
puts "ERROR: OCC26305 is reproduced. Prism is incorrect."
|
||||||
|
}
|
||||||
|
|
||||||
|
set only_screen 1
|
34
tests/bugs/modalg_6/bug26305_2
Normal file
34
tests/bugs/modalg_6/bug26305_2
Normal file
@ -0,0 +1,34 @@
|
|||||||
|
puts "========"
|
||||||
|
puts "OCC26305"
|
||||||
|
puts "========"
|
||||||
|
puts ""
|
||||||
|
###################################################
|
||||||
|
# BRepFeat_MakePrism returns inconsistent results
|
||||||
|
###################################################
|
||||||
|
|
||||||
|
circle aCircle -10 -10 0 0 0 1 1 0 0 10
|
||||||
|
mkedge anEdge aCircle
|
||||||
|
wire aWire anEdge
|
||||||
|
mkplane aCircle aWire 0
|
||||||
|
|
||||||
|
plane aFromPlane 0 0 -10 0 0 1 1 0 0
|
||||||
|
mkface aFromPlane aFromPlane
|
||||||
|
plane aToPlane 0 0 10 0 0 1 1 0 0
|
||||||
|
mkface aToPlane aToPlane
|
||||||
|
|
||||||
|
featprism aCircle aCircle aCircle 0 0 1 2 1
|
||||||
|
featperform prism aResult aFromPlane aToPlane
|
||||||
|
|
||||||
|
checkshape aResult
|
||||||
|
|
||||||
|
vinit
|
||||||
|
vdisplay aResult
|
||||||
|
vsetdispmode aResult 1
|
||||||
|
vfit
|
||||||
|
|
||||||
|
set bug_info [vreadpixel 300 340 name]
|
||||||
|
if {$bug_info == "BLACK 0"} {
|
||||||
|
puts "ERROR: OCC26305 is reproduced. Prism is incorrect."
|
||||||
|
}
|
||||||
|
|
||||||
|
set only_screen 1
|
34
tests/bugs/modalg_6/bug26305_3
Normal file
34
tests/bugs/modalg_6/bug26305_3
Normal file
@ -0,0 +1,34 @@
|
|||||||
|
puts "========"
|
||||||
|
puts "OCC26305"
|
||||||
|
puts "========"
|
||||||
|
puts ""
|
||||||
|
###################################################
|
||||||
|
# BRepFeat_MakePrism returns inconsistent results
|
||||||
|
###################################################
|
||||||
|
|
||||||
|
circle aCircle 0 0 0 0 0 1 1 0 0 10
|
||||||
|
mkedge anEdge aCircle
|
||||||
|
wire aWire anEdge
|
||||||
|
mkplane aCircle aWire 0
|
||||||
|
|
||||||
|
plane aFromPlane 0 0 10 0 0 1 1 0 0
|
||||||
|
mkface aFromPlane aFromPlane
|
||||||
|
plane aToPlane 0 0 20 0 0 1 1 0 0
|
||||||
|
mkface aToPlane aToPlane
|
||||||
|
|
||||||
|
featprism aCircle aCircle aCircle 0 0 1 2 1
|
||||||
|
featperform prism aResult aToPlane aFromPlane
|
||||||
|
|
||||||
|
checkshape aResult
|
||||||
|
|
||||||
|
vinit
|
||||||
|
vdisplay aResult
|
||||||
|
vsetdispmode aResult 1
|
||||||
|
vfit
|
||||||
|
|
||||||
|
set bug_info [vreadpixel 300 300 name]
|
||||||
|
if {$bug_info == "BLACK 0"} {
|
||||||
|
puts "ERROR: OCC26305 is reproduced. Prism is incorrect."
|
||||||
|
}
|
||||||
|
|
||||||
|
set only_screen 1
|
29
tests/bugs/modalg_6/bug26315
Normal file
29
tests/bugs/modalg_6/bug26315
Normal file
@ -0,0 +1,29 @@
|
|||||||
|
puts "========"
|
||||||
|
puts "OCC26315"
|
||||||
|
puts "========"
|
||||||
|
puts ""
|
||||||
|
#######################################################
|
||||||
|
# BRepFeat_MakeRevol fails to create revol from shape
|
||||||
|
#######################################################
|
||||||
|
|
||||||
|
circle aCircle 0 0 0 0 0 1 1 0 0 10
|
||||||
|
mkedge anEdge aCircle
|
||||||
|
wire aWire anEdge
|
||||||
|
mkplane aCircle aWire 0
|
||||||
|
|
||||||
|
plane aFromPlane 0 0 -10 0 0 1 1 0 0
|
||||||
|
mkface aFromPlane aFromPlane
|
||||||
|
plane aToPlane 0 0 10 0 0 1 1 0 0
|
||||||
|
mkface aToPlane aToPlane
|
||||||
|
|
||||||
|
featrevol aCircle aCircle aCircle 20 0 0 0 1 0 2 1
|
||||||
|
featperform revol aResult aFromPlane aToPlane
|
||||||
|
|
||||||
|
checkshape aResult
|
||||||
|
|
||||||
|
vinit
|
||||||
|
vdisplay aResult
|
||||||
|
vsetdispmode aResult 1
|
||||||
|
vfit
|
||||||
|
|
||||||
|
set only_screen 1
|
Loading…
x
Reference in New Issue
Block a user