DOC: explain the semantics of the precision preference parameters

Describe meaning of attributes in StdMeshers.xml

Fix regression of SALOME_TESTS/Grids/smesh/3D_mesh_NETGEN_02/C4
  SMESH_MesherHelper.cxx

Fix SALOME crash on a case of #http://www.salome-platform.org/forum/forum_10/47659287
  SMESH_MAT2d.cxx
This commit is contained in:
eap 2016-10-06 14:11:45 +03:00
parent e18c7bf133
commit f08f994ee3
5 changed files with 50 additions and 15 deletions

View File

@ -134,7 +134,15 @@ or in later sessions with this module according to the preferences.
when a hypothesis is modified. This allows saving time by omitting
loading data of a large mesh that is planned to be recomputed with other parameters.
- <b>Input fields precision</b>
- <b>Input fields precision</b> - allows to adjust input precision of
different parameters. The semantics of the precision values is
described in detail in <em>Using input widgets</em> chapter of GUI
documentation (Introduction to Salome Platform / Introduction to GUI /
Using input widgets). In brief: \b positive precision value is the
maximum allowed number of digits after the decimal point in the
fixed-point format; \b nagative precision value is the maximum
allowed number of significant digits in mantissa in either the
fixed-point or scientific format.
- <b>Length precision</b> - allows to adjust input precision of coordinates and dimensions.
- <b>Angular precision</b> - allows to adjust input precision of angles.
- <b>Length tolerance precision</b> - allows to adjust input precision of tolerance of coordinates and dimensions.

View File

@ -26,6 +26,31 @@
<!-- GUI customization for MESH component -->
<!-- Attributes of hypotheses/algorithms:
type - string identifier of a hyp.
label-id - hypothesis type name in Create Mesh dialog.
icon-id - not used.
group-id - (optional) integer ID of a group the hyp belongs to in Create Mesh dialog;
by default the hyp is in the last group.
priority - (optional) priority within the group; by default the hyp is last in group.
dim - dimension; defines a tab page in Create Mesh dialog.
context - (optional) allowed context: [LOCAL, GLOBAL, ANY(default)]. LOCAL - the hyp
can be only local (on sub-mesh). GLOBAL - the hyp can be only GLOBAL (on mesh).
auxiliary - (optional) Boolean. Is additional hyp or not. Default is "false".
hypos - list of types of compatible hyps of the algorithm.
opt-hypos = (optional) list of types of compatible ADDITIONAL hyps of the algorithm.
output - geometry of elements generated by the algo. Used to define compatible algos of
different dimensions. Compatible algos have equal geometries in "input" and "output".
input - geometry of elements accepted by algorithm input. Used to define compatible algos of
different dimensions. Compatible algos have equal geometries in "input" and "output".
need-hyp - (optional) Boolean. Does the algo require a hypothesis or not. Default is "false".
need-geom - (optional) Boolean. Can the algo work w/o geometry or not. Default is "true".
support-submeshes - (optional) Boolean. Does an multi-dimensional algo support sub-meshes.
Default is "false".
-->
<meshers>
<meshers-group name ="Standard Meshers"
@ -235,7 +260,7 @@
label-id ="Segments around Vertex"
icon-id ="mesh_algo_regular.png"
hypos ="SegmentLengthAroundVertex"
output ="VERTEX"
output ="NODE"
need-hyp ="true"
dim ="0"/>
@ -246,7 +271,7 @@
priority ="10"
hypos ="Adaptive1D,LocalLength,MaxLength,Arithmetic1D,GeometricProgression,StartEndLength,NumberOfSegments,Deflection1D,AutomaticLength,FixedPoints1D"
opt-hypos="Propagation,PropagOfDistribution,QuadraticMesh"
input ="VERTEX"
input ="NODE"
output ="EDGE"
need-hyp ="true"
dim ="1">
@ -274,7 +299,7 @@
priority ="20"
hypos ="Adaptive1D,LocalLength,MaxLength,Arithmetic1D,GeometricProgression,StartEndLength,NumberOfSegments,Deflection1D,AutomaticLength,FixedPoints1D"
opt-hypos="Propagation,PropagOfDistribution,QuadraticMesh"
input ="VERTEX"
input ="NODE"
output ="EDGE"
need-hyp ="true"
dim ="1">
@ -513,7 +538,7 @@
icon-id ="mesh_algo_regular.png"
group-id="1"
priority="30"
input ="VERTEX"
input ="NODE"
output ="EDGE"
dim ="1">
<python-wrap>

View File

@ -2945,9 +2945,11 @@ bool SMESH_MesherHelper::IsReversedSubMesh (const TopoDS_Face& theFace)
double u0 = GetNodeU( TopoDS::Edge( E ), nn[0], nn[1], &ok );
double u1 = GetNodeU( TopoDS::Edge( E ), nn[1], nn[0], &ok );
// check that the 2 nodes are connected with a segment (IPAL53055)
const SMDS_MeshElement* seg;
if ( SMESHDS_SubMesh* sm = GetMeshDS()->MeshElements( E ))
if ( sm->NbElements() > 0 && !GetMeshDS()->FindEdge( nn[0], nn[1] ))
ok = false;
if (( sm->NbElements() > 0 ) &&
( seg = GetMeshDS()->FindEdge( nn[0], nn[1] )))
ok = sm->Contains( seg );
if ( ok )
{
isReversed = ( u0 > u1 );

View File

@ -1608,7 +1608,6 @@ bool SMESH_Block::LoadMeshBlock(const SMDS_MeshVolume* theVolume,
const int theNode001Index,
vector<const SMDS_MeshNode*>& theOrderedNodes)
{
MESSAGE(" ::LoadMeshBlock()");
init();
SMDS_VolumeTool vTool;
@ -1736,7 +1735,6 @@ bool SMESH_Block::LoadBlockShapes(const TopoDS_Shell& theShell,
const TopoDS_Vertex& theVertex001,
TopTools_IndexedMapOfOrientedShape& theShapeIDMap )
{
MESSAGE(" ::LoadBlockShapes()");
return ( FindBlockShapes( theShell, theVertex000, theVertex001, theShapeIDMap ) &&
LoadBlockShapes( theShapeIDMap ));
}
@ -1752,8 +1750,6 @@ bool SMESH_Block::FindBlockShapes(const TopoDS_Shell& theShell,
const TopoDS_Vertex& theVertex001,
TopTools_IndexedMapOfOrientedShape& theShapeIDMap )
{
MESSAGE(" ::FindBlockShapes()");
// 8 vertices
TopoDS_Shape V000, V100, V010, V110, V001, V101, V011, V111;
// 12 edges

View File

@ -93,8 +93,8 @@ namespace
list< const TVDEdge* > _edges; // MA edges in CCW order within _cell
InSegment( InPoint * p0, InPoint * p1, size_t iE)
: _p0(p0), _p1(p1), _geomEdgeInd(iE) {}
InSegment() : _p0(0), _p1(0), _geomEdgeInd(0) {}
: _p0(p0), _p1(p1), _geomEdgeInd(iE), _cell(0) {}
InSegment() : _p0(0), _p1(0), _geomEdgeInd(0), _cell(0) {}
const InPoint& point0() const { return *_p0; }
const InPoint& point1() const { return *_p1; }
@ -662,7 +662,7 @@ namespace
// get scale to have the same 2d proportions as in 3d
computeProportionScale( face, uvBox, scale );
// make scale to have coordinates precise enough when converted to int
// make 'scale' such that to have coordinates precise enough when converted to int
gp_XY uvMin = uvBox.CornerMin(), uvMax = uvBox.CornerMax();
uvMin.ChangeCoord(1) = uvMin.X() * scale[0];
@ -672,7 +672,7 @@ namespace
double vMax[2] = { Max( Abs( uvMin.X() ), Abs( uvMax.X() )),
Max( Abs( uvMin.Y() ), Abs( uvMax.Y() )) };
int iMax = ( vMax[0] > vMax[1] ) ? 0 : 1;
const double precision = 1e-5;
const double precision = Min( 1e-5, minSegLen * 1e-2 );
double preciScale = Min( vMax[iMax] / precision,
std::numeric_limits<int>::max() / vMax[iMax] );
preciScale /= scale[iMax];
@ -703,6 +703,8 @@ namespace
{
inPoints[ iP++ ] = points[i-1].getInPoint( scale );
inSegments.push_back( InSegment( & inPoints[ iP-2 ], & inPoints[ iP-1 ], iE ));
if ( inPoints[ iP-2 ] == inPoints[ iP-1 ])
return false; // too short segment
}
}
}
@ -716,6 +718,8 @@ namespace
{
inPoints[ iP++ ] = points[i].getInPoint( scale );
inSegments.push_back( InSegment( & inPoints[ iP-2 ], & inPoints[ iP-1 ], iE ));
if ( inPoints[ iP-2 ] == inPoints[ iP-1 ])
return false; // too short segment
}
}
}