1
0
mirror of https://git.dev.opencascade.org/repos/occt.git synced 2025-08-09 13:22:24 +03:00

0030008: BRepMesh does not respect angular deflection in internal area of bspline surface

1. Check whether the mesh satisfies the required angular deflection has been amended. Namely normals (to the surface) in the ends of any not "frontier" link are made collinear (with the given angular tolerance).

2. New parameters AngleInterior and DeflectionInterior have been added in IMeshTools_Parameters structure.

3. In case of thin long faces with internal edges, add points of internal edges to control parameters using grabParamsOfInternalEdges() in order to avoid aberrations on its ends. Disable addition of parameters from boundary edges in case of BSpline surface. Deviation can be controlled through the deflection parameter.

4. Grab parameters from edges in case if there is just a single interval on BSpline surface along U and V direction.
This commit is contained in:
nbv
2018-10-19 16:38:02 +03:00
committed by bugmaster
parent 7bd071edb1
commit 46478ffe32
26 changed files with 354 additions and 107 deletions

View File

@@ -44,6 +44,11 @@ public:
protected:
//! Initializes U and V parameters lists using CN continuity intervals.
virtual Standard_Boolean initParameters() const Standard_OVERRIDE
{
return Standard_True;
}
};
#endif

View File

@@ -14,41 +14,46 @@
// commercial license or contractual agreement.
#include <BRepMesh_Deflection.hxx>
#include <BRepMesh_ShapeTool.hxx>
#include <IMeshTools_Parameters.hxx>
#include <IMeshData_Edge.hxx>
#include <IMeshData_Wire.hxx>
#include <IMeshData_Face.hxx>
#include <BRep_Tool.hxx>
#include <Bnd_Box.hxx>
#include <BRepBndLib.hxx>
#include <TopExp_Explorer.hxx>
#include <TopoDS_Vertex.hxx>
#include <BRepMesh_ShapeTool.hxx>
#include <IMeshData_Edge.hxx>
#include <IMeshData_Wire.hxx>
#include <IMeshTools_Parameters.hxx>
#include <TopExp.hxx>
#include <TopoDS_Vertex.hxx>
//=======================================================================
//function : RelativeEdgeDeflection
//purpose :
//=======================================================================
Standard_Real BRepMesh_Deflection::RelativeEdgeDeflection(
const TopoDS_Edge& theEdge,
const Standard_Real theDeflection,
Standard_Real BRepMesh_Deflection::ComputeAbsoluteDeflection(
const TopoDS_Shape& theShape,
const Standard_Real theRelativeDeflection,
const Standard_Real theMaxShapeSize,
Standard_Real& theAdjustmentCoefficient)
{
theAdjustmentCoefficient = 1.;
Standard_Real aEdgeDeflection = theDeflection;
if (theEdge.IsNull())
if (theShape.IsNull())
{
return aEdgeDeflection;
return theRelativeDeflection;
}
Bnd_Box aBox;
BRepBndLib::Add (theEdge, aBox, Standard_False);
BRepMesh_ShapeTool::BoxMaxDimension (aBox, aEdgeDeflection);
BRepBndLib::Add (theShape, aBox, Standard_False);
Standard_Real aShapeSize = theRelativeDeflection;
BRepMesh_ShapeTool::BoxMaxDimension (aBox, aShapeSize);
// Adjust resulting value in relation to the total size
theAdjustmentCoefficient = theMaxShapeSize / (2 * aEdgeDeflection);
Standard_Real aX1, aY1, aZ1, aX2, aY2, aZ2;
aBox.Get(aX1, aY1, aZ1, aX2, aY2, aZ2);
const Standard_Real aMaxShapeSize = (theMaxShapeSize > 0.0) ? theMaxShapeSize :
Max(aX2 - aX1, Max(aY2 - aY1, aZ2 - aZ1));
theAdjustmentCoefficient = aMaxShapeSize / (2 * aShapeSize);
if (theAdjustmentCoefficient < 0.5)
{
theAdjustmentCoefficient = 0.5;
@@ -58,7 +63,7 @@ Standard_Real BRepMesh_Deflection::RelativeEdgeDeflection(
theAdjustmentCoefficient = 2.;
}
return (theAdjustmentCoefficient * aEdgeDeflection * theDeflection);
return (theAdjustmentCoefficient * aShapeSize * theRelativeDeflection);
}
//=======================================================================
@@ -75,8 +80,9 @@ void BRepMesh_Deflection::ComputeDeflection (
if (theParameters.Relative)
{
Standard_Real aScale;
aLinDeflection = RelativeEdgeDeflection (theDEdge->GetEdge (),
theParameters.Deflection, theMaxShapeSize, aScale);
aLinDeflection = ComputeAbsoluteDeflection(theDEdge->GetEdge(),
theParameters.Deflection,
theMaxShapeSize, aScale);
// Is it OK?
aAngDeflection = theParameters.Angle * aScale;
@@ -144,7 +150,15 @@ void BRepMesh_Deflection::ComputeDeflection (
const IMeshData::IFaceHandle& theDFace,
const IMeshTools_Parameters& theParameters)
{
Standard_Real aFaceDeflection = 0.;
Standard_Real aDeflection = theParameters.DeflectionInterior;
if (theParameters.Relative)
{
Standard_Real aScale;
aDeflection = ComputeAbsoluteDeflection(theDFace->GetFace(),
aDeflection, -1.0, aScale);
}
Standard_Real aFaceDeflection = 0.0;
if (theDFace->WiresNb () > 0)
{
for (Standard_Integer aWireIt = 0; aWireIt < theDFace->WiresNb(); ++aWireIt)
@@ -154,10 +168,8 @@ void BRepMesh_Deflection::ComputeDeflection (
aFaceDeflection /= theDFace->WiresNb ();
}
else
{
aFaceDeflection = theParameters.Deflection;
}
aFaceDeflection = Max(aDeflection, aFaceDeflection);
theDFace->SetDeflection (Max(2.* BRepMesh_ShapeTool::MaxFaceTolerance(
theDFace->GetFace()), aFaceDeflection));

View File

@@ -31,16 +31,17 @@ class BRepMesh_Deflection : public Standard_Transient
{
public:
//! Returns relative deflection for edge with respect to shape size.
//! @param theEdge edge for which relative deflection should be computed.
//! @param theDeflection absolute deflection.
//! @param theMaxShapeSize maximum size of a shape.
//! Returns absolute deflection for theShape with respect to the
//! relative deflection and theMaxShapeSize.
//! @param theShape shape for that the deflection should be computed.
//! @param theRelativeDeflection relative deflection.
//! @param theMaxShapeSize maximum size of the whole shape.
//! @param theAdjustmentCoefficient coefficient of adjustment between maximum
//! size of shape and calculated relative deflection.
//! @return relative deflection for the edge.
Standard_EXPORT static Standard_Real RelativeEdgeDeflection (
const TopoDS_Edge& theEdge,
const Standard_Real theDeflection,
//! @return absolute deflection for the shape.
Standard_EXPORT static Standard_Real ComputeAbsoluteDeflection (
const TopoDS_Shape& theShape,
const Standard_Real theRelativeDeflection,
const Standard_Real theMaxShapeSize,
Standard_Real& theAdjustmentCoefficient);

View File

@@ -18,6 +18,7 @@
#include <BRepMesh_DelaunayNodeInsertionMeshAlgo.hxx>
#include <BRepMesh_GeomTool.hxx>
#include <GeomLib.hxx>
//! Extends node insertion Delaunay meshing algo in order to control
//! deflection of generated trianges. Splits triangles failing the check.
@@ -314,16 +315,56 @@ private:
const gp_XY aMidPnt2d = (theNodesInfo[i].Point2d +
theNodesInfo[j].Point2d) / 2.;
usePoint(aMidPnt2d, LineDeviation(theNodesInfo[i].Point, theNodesInfo[j].Point));
if (!usePoint (aMidPnt2d, LineDeviation (theNodesInfo[i].Point,
theNodesInfo[j].Point)))
{
if (!checkLinkEndsForAngularDeviation(theNodesInfo[i],
theNodesInfo[j],
aMidPnt2d))
{
myControlNodes->Append(aMidPnt2d);
}
}
}
}
}
//! Checks the given point (located between the given nodes)
//! for specified angular deviation.
Standard_Boolean checkLinkEndsForAngularDeviation(const TriangleNodeInfo& theNodeInfo1,
const TriangleNodeInfo& theNodeInfo2,
const gp_XY& /*theMidPoint*/)
{
gp_Dir aNorm1, aNorm2;
const Handle(Geom_Surface)& aSurf =
this->getDFace()->GetSurface()->ChangeSurface().Surface().Surface();
if ((GeomLib::NormEstim(aSurf, theNodeInfo1.Point2d, Precision::Confusion(), aNorm1) == 0) &&
(GeomLib::NormEstim(aSurf, theNodeInfo2.Point2d, Precision::Confusion(), aNorm2) == 0))
{
Standard_Real anAngle = aNorm1.Angle(aNorm2);
if (anAngle > this->getParameters().AngleInterior)
return Standard_False;
}
#if 0
else if (GeomLib::NormEstim(aSurf, theMidPoint, Precision::Confusion(), aNorm1) != 0)
{
// It is better to consider the singular point as a node of triangulation.
// However, it leads to hangs up meshing some faces (including faces with
// degenerated edges). E.g. tests "mesh standard_incmesh Q6".
// So, this code fragment is better to implement in the future.
return Standard_False;
}
#endif
return Standard_True;
}
//! Computes deflection of the given point and caches it for
//! insertion in case if it overflows deflection.
//! @return True if point has been cached for insertion.
template<class DeflectionFunctor>
inline void usePoint(
inline Standard_Boolean usePoint(
const gp_XY& thePnt2d,
const DeflectionFunctor& theDeflectionFunctor)
{
@@ -332,7 +373,10 @@ private:
if (!checkDeflectionOfPointAndUpdateCache(thePnt2d, aPnt, theDeflectionFunctor.SquareDeviation(aPnt)))
{
myControlNodes->Append(thePnt2d);
return Standard_True;
}
return Standard_False;
}
//! Checks the given point for specified linear deflection.

View File

@@ -46,7 +46,8 @@ protected:
virtual std::pair<Standard_Integer, Standard_Integer> getCellsCount (const Standard_Integer theVerticesNb) Standard_OVERRIDE
{
return BRepMesh_GeomTool::CellsCount (this->getDFace()->GetSurface(), theVerticesNb,
this->getParameters().Deflection, &this->getRangeSplitter());
this->getDFace()->GetDeflection(),
&this->getRangeSplitter());
}
//! Perfroms processing of generated mesh. Generates surface nodes and inserts them into structure.

View File

@@ -51,7 +51,7 @@ public: //! @name mesher API
Standard_EXPORT BRepMesh_IncrementalMesh(const TopoDS_Shape& theShape,
const IMeshTools_Parameters& theParameters);
//! Performs meshing ot the shape.
//! Performs meshing ot the shape.
Standard_EXPORT virtual void Perform() Standard_OVERRIDE;
public: //! @name accessing to parameters.
@@ -85,12 +85,23 @@ private:
//! Initializes specific parameters
inline void initParameters()
{
if (myParameters.DeflectionInterior < Precision::Confusion())
{
myParameters.DeflectionInterior = myParameters.Deflection;
}
if (myParameters.MinSize < Precision::Confusion())
{
myParameters.MinSize =
Max(IMeshTools_Parameters::RelMinSize() * myParameters.Deflection,
Max(IMeshTools_Parameters::RelMinSize() * Min(myParameters.Deflection,
myParameters.DeflectionInterior),
Precision::Confusion());
}
if (myParameters.AngleInterior < Precision::Angular())
{
myParameters.AngleInterior = 2.0 * myParameters.Angle;
}
}
public: //! @name plugin API

View File

@@ -70,7 +70,8 @@ Handle (IMeshData_Model) BRepMesh_ModelBuilder::Perform (
}
else
{
aModel->SetMaxSize(theParameters.Deflection);
aModel->SetMaxSize(Max(theParameters.Deflection,
theParameters.DeflectionInterior));
}
Handle (IMeshTools_ShapeVisitor) aVisitor =

View File

@@ -14,15 +14,14 @@
// commercial license or contractual agreement.
#include <BRepMesh_NURBSRangeSplitter.hxx>
#include <GeomAdaptor_Curve.hxx>
#include <IMeshTools_Parameters.hxx>
#include <IMeshData_Wire.hxx>
#include <IMeshData_Edge.hxx>
#include <IMeshData_PCurve.hxx>
#include <GeomAbs_IsoType.hxx>
#include <BRepMesh_GeomTool.hxx>
#include <NCollection_Handle.hxx>
#include <algorithm>
#include <BRepMesh_GeomTool.hxx>
#include <GeomAdaptor_Curve.hxx>
#include <GeomLib.hxx>
#include <IMeshData_Edge.hxx>
#include <IMeshData_Wire.hxx>
#include <NCollection_Handle.hxx>
namespace
{
@@ -55,13 +54,6 @@ namespace
{
myParameters = theParameters;
if (myParameters.MinSize <= Precision::Confusion())
{
myParameters.MinSize =
Max(IMeshTools_Parameters::RelMinSize() * myParameters.Deflection,
Precision::Confusion());
}
Standard_Integer aStartIndex, aEndIndex;
if (myIsoU)
{
@@ -106,8 +98,18 @@ namespace
const Standard_Real aSqDist = BRepMesh_GeomTool::SquareDeflectionOfSegment(
myPrevControlPnt, myCurrControlPnt, aMidPnt);
const Standard_Real aSqMaxDeflection = myDFace->GetDeflection() * myDFace->GetDeflection();
if ((aSqDist > aSqMaxDeflection) &&
Standard_Real anAngle = 0.0;
if ((myPrevControlVec.SquareMagnitude() > Precision::SquareConfusion()) &&
(myCurrControlVec.SquareMagnitude() > Precision::SquareConfusion()))
{
anAngle = myPrevControlVec.Angle(myCurrControlVec);
}
const Standard_Real aSqMaxDeflection = myDFace->GetDeflection() *
myDFace->GetDeflection();
if (((aSqDist > aSqMaxDeflection) || (anAngle > myParameters.AngleInterior)) &&
aSqDist > myParameters.MinSize * myParameters.MinSize)
{
// insertion
@@ -120,9 +122,8 @@ namespace
// internals in order to prevent movement of triangle body
// outside the surface in case of highly curved ones, e.g.
// BSpline springs.
if (aSqDist < aSqMaxDeflection &&
myControlParams->Length() > 3 &&
theIndex < myControlParams->Length())
if (((aSqDist < aSqMaxDeflection) || (anAngle < myParameters.AngleInterior)) &&
myControlParams->Length() > 3 && theIndex < myControlParams->Length())
{
// Remove too dense points
const Standard_Real aTmpParam = myControlParams->Value(theIndex + 1);
@@ -162,7 +163,7 @@ namespace
// Lets check parameters for angular deflection.
if (myPrevControlVec.SquareMagnitude() < gp::Resolution() ||
aTmpVec.SquareMagnitude() < gp::Resolution() ||
myPrevControlVec.Angle(aTmpVec) < myParameters.Angle)
myPrevControlVec.Angle(aTmpVec) < myParameters.AngleInterior)
{
// For current Iso line we can remove this parameter.
myControlParamsToRemove->Add(myCurrControlParam);
@@ -257,6 +258,32 @@ namespace
return isAdded;
}
//! Checks whether intervals should be split.
//! Returns true in case if it is impossible to compute normal
//! directly on intervals, false is returned elsewhere.
Standard_Boolean toSplitIntervals (const Handle (Geom_Surface)& theSurf,
const TColStd_Array1OfReal (&theIntervals)[2])
{
Standard_Integer aIntervalU = theIntervals[0].Lower ();
for (; aIntervalU <= theIntervals[0].Upper (); ++aIntervalU)
{
const Standard_Real aParamU = theIntervals[0].Value(aIntervalU);
Standard_Integer aIntervalV = theIntervals[1].Lower ();
for (; aIntervalV <= theIntervals[1].Upper (); ++aIntervalV)
{
gp_Dir aNorm;
const Standard_Real aParamV = theIntervals[1].Value(aIntervalV);
if (GeomLib::NormEstim (theSurf, gp_Pnt2d (aParamU, aParamV), Precision::Confusion (), aNorm) != 0)
{
return Standard_True;
}
// TODO: do not split intervals if there is no normal in the middle of interval.
}
}
return Standard_False;
}
}
//=======================================================================
@@ -287,7 +314,10 @@ void BRepMesh_NURBSRangeSplitter::AdjustRange()
Handle(IMeshData::ListOfPnt2d) BRepMesh_NURBSRangeSplitter::GenerateSurfaceNodes(
const IMeshTools_Parameters& theParameters) const
{
initParameters();
if (!initParameters())
{
return Handle(IMeshData::ListOfPnt2d)();
}
const std::pair<Standard_Real, Standard_Real>& aRangeU = GetRangeU();
const std::pair<Standard_Real, Standard_Real>& aRangeV = GetRangeV();
@@ -357,7 +387,7 @@ Handle(IMeshData::ListOfPnt2d) BRepMesh_NURBSRangeSplitter::GenerateSurfaceNodes
// Function: initParameters
// Purpose :
//=======================================================================
void BRepMesh_NURBSRangeSplitter::initParameters() const
Standard_Boolean BRepMesh_NURBSRangeSplitter::initParameters() const
{
const Handle(BRepAdaptor_HSurface)& aSurface = GetSurface();
@@ -375,21 +405,79 @@ void BRepMesh_NURBSRangeSplitter::initParameters() const
aSurface->UIntervals(aIntervals[0], aContinuity);
aSurface->VIntervals(aIntervals[1], aContinuity);
Standard_Boolean isSplitIntervals =
(aIntervalsNb.first > 1 || aIntervalsNb.second > 1);
const Standard_Boolean isSplitIntervals = toSplitIntervals (
aSurface->ChangeSurface().Surface().Surface(), aIntervals);
if (!isSplitIntervals &&
(aSurface->GetType() == GeomAbs_BezierSurface ||
aSurface->GetType() == GeomAbs_BSplineSurface))
if (!initParamsFromIntervals(aIntervals[0], GetRangeU(), isSplitIntervals,
const_cast<IMeshData::IMapOfReal&>(GetParametersU())))
{
isSplitIntervals = (aSurface->NbUPoles() > 2 && aSurface->NbVPoles() > 2);
//if (!grabParamsOfEdges (Edge_Frontier, Param_U))
{
return Standard_False;
}
}
initParamsFromIntervals(aIntervals[0], GetRangeU(), isSplitIntervals,
const_cast<IMeshData::IMapOfReal&>(GetParametersU()));
if (!initParamsFromIntervals(aIntervals[1], GetRangeV(), isSplitIntervals,
const_cast<IMeshData::IMapOfReal&>(GetParametersV())))
{
//if (!grabParamsOfEdges (Edge_Frontier, Param_V))
{
return Standard_False;
}
}
initParamsFromIntervals(aIntervals[1], GetRangeV(), isSplitIntervals,
const_cast<IMeshData::IMapOfReal&>(GetParametersV()));
return grabParamsOfEdges(Edge_Internal, Param_U | Param_V);
}
//=======================================================================
//function : grabParamsOfInternalEdges
//purpose :
//=======================================================================
Standard_Boolean BRepMesh_NURBSRangeSplitter::grabParamsOfEdges (
const EdgeType theEdgeType,
const Standard_Integer theParamDimensionFlag) const
{
if ((theParamDimensionFlag & (Param_U | Param_V)) == 0)
{
return Standard_False;
}
const IMeshData::IFaceHandle& aDFace = GetDFace ();
for (Standard_Integer aWireIt = 0; aWireIt < aDFace->WiresNb (); ++aWireIt)
{
const IMeshData::IWireHandle& aDWire = aDFace->GetWire (aWireIt);
for (Standard_Integer aEdgeIt = 0; aEdgeIt < aDWire->EdgesNb (); ++aEdgeIt)
{
const IMeshData::IEdgePtr& aDEdge = aDWire->GetEdge (aEdgeIt);
for (Standard_Integer aPCurveIt = 0; aPCurveIt < aDEdge->PCurvesNb (); ++aPCurveIt)
{
const IMeshData::IPCurveHandle& aDPCurve = aDEdge->GetPCurve (aPCurveIt);
if (aDPCurve->GetFace () == aDFace)
{
if (theEdgeType == Edge_Internal && !aDPCurve->IsInternal ())
{
continue;
}
for (Standard_Integer aPointIt = 0; aPointIt < aDPCurve->ParametersNb (); ++aPointIt)
{
const gp_Pnt2d& aPnt2d = aDPCurve->GetPoint (aPointIt);
if (theParamDimensionFlag & Param_U)
{
const_cast<IMeshData::IMapOfReal&>(GetParametersU ()).Add (aPnt2d.X ());
}
if (theParamDimensionFlag & Param_V)
{
const_cast<IMeshData::IMapOfReal&>(GetParametersV ()).Add (aPnt2d.Y ());
}
}
}
}
}
}
return Standard_True;
}
//=======================================================================
@@ -411,16 +499,13 @@ Handle(IMeshData::SequenceOfReal) BRepMesh_NURBSRangeSplitter::computeGrainAndFi
aMinDiff /= theDelta;
}
const Standard_Real aMinSize =
theParameters.MinSize > Precision::Confusion() ? theParameters.MinSize :
Max(IMeshTools_Parameters::RelMinSize() * theParameters.Deflection,
Precision::Confusion());
aMinDiff = Max(aMinSize, aMinDiff);
aMinDiff = Max(theParameters.MinSize, aMinDiff);
const Standard_Real aDiffMaxLim = 0.1 * theRangeDiff;
const Standard_Real aDiffMinLim = Max(0.005 * theRangeDiff, 2. * theTol2d);
const Standard_Real aDiff = Max(aMinSize, Min(aDiffMaxLim, aDiffMinLim));
const Standard_Real aDiffMinLim = Max(0.005 * theRangeDiff,
2. * theTol2d);
const Standard_Real aDiff = Max(theParameters.MinSize,
Min(aDiffMaxLim, aDiffMinLim));
return filterParameters(theSourceParams, aMinDiff, aDiff, theAllocator);
}

View File

@@ -46,7 +46,7 @@ public:
protected:
//! Initializes U and V parameters lists using CN continuity intervals.
Standard_EXPORT virtual void initParameters() const;
Standard_EXPORT virtual Standard_Boolean initParameters() const;
private:
@@ -66,6 +66,23 @@ private:
const Standard_Real theFilterDist,
const Handle(NCollection_IncAllocator)& theAllocator) const;
enum EdgeType
{
Edge_Internal,
Edge_Frontier
};
enum ParamDimension
{
Param_U = 0x1,
Param_V = 0x2
};
//! Finds edges of discrete face and uses its points
//! as auxiliary control parameters for generation of nodes.
Standard_Boolean grabParamsOfEdges (const EdgeType theEdgeType,
const Standard_Integer theParamDimensionFlag) const;
private:
GeomAbs_SurfaceType mySurfaceType;

View File

@@ -26,6 +26,8 @@ struct IMeshTools_Parameters {
:
Angle(0.5),
Deflection(0.001),
AngleInterior(-1.0),
DeflectionInterior(-1.0),
MinSize (-1.0),
InParallel (Standard_False),
Relative (Standard_False),
@@ -42,12 +44,18 @@ struct IMeshTools_Parameters {
return 0.1;
}
//! Angular deflection
//! Angular deflection used to tessellate the boundary edges
Standard_Real Angle;
//! Deflection
//!Linear deflection used to tessellate the boundary edges
Standard_Real Deflection;
//! Angular deflection used to tessellate the face interior
Standard_Real AngleInterior;
//! Linear deflection used to tessellate the face interior
Standard_Real DeflectionInterior;
//! Minimal allowed size of mesh element
Standard_Real MinSize;

View File

@@ -118,8 +118,11 @@ static Standard_Integer incrementalmesh(Draw_Interpretor& di, Standard_Integer n
Builds triangular mesh for the shape\n\
usage: incmesh Shape LinearDeflection [options]\n\
options:\n\
-a val angular deflection in deg\n\
-a val angular deflection for edges in deg\n\
(default ~28.64 deg = 0.5 rad)\n\n\
-ai val angular deflection inside of faces in deg\n\
(default ~57.29 deg = 1 rad)\n\n\
-di val Linear deflection used to tessellate the face interior.\n\
-min minimum size parameter limiting size of triangle's\n\
edges to prevent sinking into amplification in case\n\
of distorted curves and surfaces\n\n\
@@ -142,6 +145,8 @@ options:\n\
}
IMeshTools_Parameters aMeshParams;
aMeshParams.Deflection = aMeshParams.DeflectionInterior =
Max(Draw::Atof(argv[2]), Precision::Confusion());
if (nbarg > 3)
{
@@ -168,8 +173,16 @@ options:\n\
{
aMeshParams.Angle = aVal * M_PI / 180.;
}
else if (aOpt == "-ai")
{
aMeshParams.AngleInterior = aVal * M_PI / 180.;
}
else if (aOpt == "-min")
aMeshParams.MinSize = aVal;
else if (aOpt == "-di")
{
aMeshParams.DeflectionInterior = aVal;
}
else
--i;
}