From: jfa Date: Fri, 1 Feb 2008 08:07:48 +0000 (+0000) Subject: NPAL17873: SMESH add a triangle instead of make quadrangles only. X-Git-Tag: for_M2008_07022008~5 X-Git-Url: http://git.salome-platform.org/gitweb/?a=commitdiff_plain;h=f27fb702a4d976ee6b422a4b931d1bb5d74a3b96;p=modules%2Fsmesh.git NPAL17873: SMESH add a triangle instead of make quadrangles only. --- diff --git a/doc/salome/gui/SMESH/images/a-averagelength.png b/doc/salome/gui/SMESH/images/a-averagelength.png index dc007eb57..70e2afd26 100755 Binary files a/doc/salome/gui/SMESH/images/a-averagelength.png and b/doc/salome/gui/SMESH/images/a-averagelength.png differ diff --git a/doc/salome/gui/SMESH/input/1d_meshing_hypo.doc b/doc/salome/gui/SMESH/input/1d_meshing_hypo.doc index 6fd23f521..716ddd90d 100644 --- a/doc/salome/gui/SMESH/input/1d_meshing_hypo.doc +++ b/doc/salome/gui/SMESH/input/1d_meshing_hypo.doc @@ -54,9 +54,18 @@ locations and 1D mesh elements are constructed on segments. Average Length hypothesis can be applied for meshing of edges composing your geometrical object. Definition of this hypothesis consists of setting the \b length of segments, which will split these -edges. The points on the edges generated by these segments will -represent nodes of your mesh. Later these nodes will be used for -meshing of the faces abutting to these edges. +edges, and the \b precision of rounding. The points on the edges +generated by these segments will represent nodes of your mesh. +Later these nodes will be used for meshing of the faces abutting to +these edges. + +The \b precision parameter is used to allow rounding a number of +segments, calculated from the edge length and average length of +segment, to the lower integer, if this value outstands from it in +bounds of the precision. Otherwise, the number of segments is rounded +to the higher integer. Use value 0.5 to provide rounding to the +nearest integer, 1.0 for the lower integer, 0.0 for the higher +integer. Default value is 1e-07. \image html image41.gif @@ -155,4 +164,4 @@ minimum and maximum value of this parameter. \image html image148.gif -*/ \ No newline at end of file +*/ diff --git a/idl/SMESH_BasicHypothesis.idl b/idl/SMESH_BasicHypothesis.idl index 9fe9370ee..5978aa052 100644 --- a/idl/SMESH_BasicHypothesis.idl +++ b/idl/SMESH_BasicHypothesis.idl @@ -43,13 +43,33 @@ module StdMeshers /*! * Sets parameter value */ - void SetLength(in double length) + void SetLength(in double length) + raises (SALOME::SALOME_Exception); + + /*! + * Sets parameter value + * + * Precision parameter is used to allow rounding a number of segments, + * calculated from the edge length and average length of segment, + * to the lower integer, if this value outstands from it in bounds of the precision. + * Otherwise, the number of segments is rounded to the higher integer. + * Use value 0.5 to provide rounding to the nearest integer, + * 1.0 for the lower integer, 0.0 for the higher integer. + * Default value is 1e-07. In old studies, restored from file, + * this value will be set to zero, what corresponds to the old behaviour. + */ + void SetPrecision(in double precision) raises (SALOME::SALOME_Exception); /*! * Returns parameter value */ double GetLength(); + + /*! + * Returns parameter value + */ + double GetPrecision(); }; /*! @@ -85,7 +105,7 @@ module StdMeshers /*! * Sets parameter value */ - void SetNumberOfSegments(in long segmentsNumber) + void SetNumberOfSegments(in long segmentsNumber) raises (SALOME::SALOME_Exception); /*! diff --git a/src/StdMeshers/StdMeshers_CompositeSegment_1D.cxx b/src/StdMeshers/StdMeshers_CompositeSegment_1D.cxx index c60cbee61..381caade1 100644 --- a/src/StdMeshers/StdMeshers_CompositeSegment_1D.cxx +++ b/src/StdMeshers/StdMeshers_CompositeSegment_1D.cxx @@ -119,7 +119,7 @@ namespace { return; for ( int iE = 0; iE < side.NbEdges(); ++iE ) { - // set listener and its data + // set listener and its data EventListenerData * listenerData = new EventListenerData(true); const TopoDS_Edge& edge = side.Edge( iE ); SMESH_subMesh * sm = side.GetMesh()->GetSubMesh( edge ); @@ -333,7 +333,7 @@ bool StdMeshers_CompositeSegment_1D::Compute(SMESH_Mesh & aMesh, auto_ptr< BRepAdaptor_CompCurve > C3d ( side->GetCurve3d() ); double f = C3d->FirstParameter(), l = C3d->LastParameter(); list< double > params; - if ( !computeInternalParameters ( *C3d, side->Length(), f, l, params, false )) + if ( !computeInternalParameters ( aMesh, *C3d, side->Length(), f, l, params, false )) return false; // Redistribute parameters near ends @@ -349,7 +349,7 @@ bool StdMeshers_CompositeSegment_1D::Compute(SMESH_Mesh & aMesh, const SMDS_MeshNode * nFirst = SMESH_Algo::VertexNode( VFirst, meshDS ); const SMDS_MeshNode * nLast = SMESH_Algo::VertexNode( VLast, meshDS ); - if (!nFirst) + if (!nFirst) return error(COMPERR_BAD_INPUT_MESH, TComm("No node on vertex ") <ShapeToIndex(VFirst)); if (!nLast) @@ -414,4 +414,3 @@ bool StdMeshers_CompositeSegment_1D::Compute(SMESH_Mesh & aMesh, return true; } - diff --git a/src/StdMeshers/StdMeshers_LocalLength.cxx b/src/StdMeshers/StdMeshers_LocalLength.cxx index 6fd6bbcf7..e3fe6e898 100644 --- a/src/StdMeshers/StdMeshers_LocalLength.cxx +++ b/src/StdMeshers/StdMeshers_LocalLength.cxx @@ -43,6 +43,7 @@ #include #include #include +#include using namespace std; @@ -56,6 +57,7 @@ StdMeshers_LocalLength::StdMeshers_LocalLength(int hypId, int studyId, SMESH_Gen :SMESH_Hypothesis(hypId, studyId, gen) { _length = 1.; + _precision = Precision::Confusion(); _name = "LocalLength"; _param_algo_dim = 1; // is used by SMESH_Regular_1D } @@ -78,12 +80,13 @@ StdMeshers_LocalLength::~StdMeshers_LocalLength() void StdMeshers_LocalLength::SetLength(double length) throw(SALOME_Exception) { - double oldLength = _length; - if (length <= 0) - throw SALOME_Exception(LOCALIZED("length must be positive")); - _length = length; - if (oldLength != _length) - NotifySubMeshesHypothesisModification(); + double oldLength = _length; + if (length <= 0) + throw SALOME_Exception(LOCALIZED("length must be positive")); + _length = length; + const double precision = 1e-7; + if (fabs(oldLength - _length) > precision) + NotifySubMeshesHypothesisModification(); } //============================================================================= @@ -94,7 +97,33 @@ void StdMeshers_LocalLength::SetLength(double length) throw(SALOME_Exception) double StdMeshers_LocalLength::GetLength() const { - return _length; + return _length; +} + +//============================================================================= +/*! + * + */ +//============================================================================= +void StdMeshers_LocalLength::SetPrecision (double thePrecision) throw(SALOME_Exception) +{ + double oldPrecision = _precision; + if (_precision < 0) + throw SALOME_Exception(LOCALIZED("precision cannot be negative")); + _precision = thePrecision; + const double precision = 1e-8; + if (fabs(oldPrecision - _precision) > precision) + NotifySubMeshesHypothesisModification(); +} + +//============================================================================= +/*! + * + */ +//============================================================================= +double StdMeshers_LocalLength::GetPrecision() const +{ + return _precision; } //============================================================================= @@ -105,7 +134,7 @@ double StdMeshers_LocalLength::GetLength() const ostream & StdMeshers_LocalLength::SaveTo(ostream & save) { - save << this->_length; + save << this->_length << " " << this->_precision; return save; } @@ -119,11 +148,23 @@ istream & StdMeshers_LocalLength::LoadFrom(istream & load) { bool isOK = true; double a; + isOK = (load >> a); if (isOK) this->_length = a; else load.clear(ios::badbit | load.rdstate()); + + isOK = (load >> a); + if (isOK) + this->_precision = a; + else + { + load.clear(ios::badbit | load.rdstate()); + // old format, without precision + _precision = 0.; + } + return load; } @@ -190,5 +231,7 @@ bool StdMeshers_LocalLength::SetParametersByMesh(const SMESH_Mesh* theMesh, if ( nbEdges ) _length /= nbEdges; + _precision = Precision::Confusion(); + return nbEdges; } diff --git a/src/StdMeshers/StdMeshers_LocalLength.hxx b/src/StdMeshers/StdMeshers_LocalLength.hxx index 79a85df7f..93fc49a6f 100644 --- a/src/StdMeshers/StdMeshers_LocalLength.hxx +++ b/src/StdMeshers/StdMeshers_LocalLength.hxx @@ -35,15 +35,17 @@ #include "SMESH_Hypothesis.hxx" #include "Utils_SALOME_Exception.hxx" -class STDMESHERS_EXPORT StdMeshers_LocalLength:public SMESH_Hypothesis +class STDMESHERS_EXPORT StdMeshers_LocalLength: public SMESH_Hypothesis { public: StdMeshers_LocalLength(int hypId, int studyId, SMESH_Gen * gen); virtual ~ StdMeshers_LocalLength(); void SetLength(double length) throw(SALOME_Exception); + void SetPrecision(double precision) throw(SALOME_Exception); double GetLength() const; + double GetPrecision() const; virtual std::ostream & SaveTo(std::ostream & save); virtual std::istream & LoadFrom(std::istream & load); @@ -60,6 +62,7 @@ class STDMESHERS_EXPORT StdMeshers_LocalLength:public SMESH_Hypothesis protected: double _length; + double _precision; }; #endif diff --git a/src/StdMeshers/StdMeshers_RadialPrism_3D.cxx b/src/StdMeshers/StdMeshers_RadialPrism_3D.cxx index cc1b4f89e..4bcd5519a 100644 --- a/src/StdMeshers/StdMeshers_RadialPrism_3D.cxx +++ b/src/StdMeshers/StdMeshers_RadialPrism_3D.cxx @@ -160,7 +160,7 @@ bool StdMeshers_RadialPrism_3D::Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& a // to delete helper at exit from Compute() std::auto_ptr helperDeleter( myHelper ); - // get 2 shells + // get 2 shells TopoDS_Solid solid = TopoDS::Solid( aShape ); TopoDS_Shell outerShell = BRepTools::OuterShell( solid ); TopoDS_Shape innerShell; @@ -336,7 +336,7 @@ public: BRepAdaptor_Curve C3D(edge); double f = C3D.FirstParameter(), l = C3D.LastParameter(); list< double > params; - if ( !StdMeshers_Regular_1D::computeInternalParameters( C3D, len, f, l, params, false )) + if ( !StdMeshers_Regular_1D::computeInternalParameters( aMesh, C3D, len, f, l, params, false )) return error("StdMeshers_Regular_1D failed to compute layers distribution"); positions.clear(); diff --git a/src/StdMeshers/StdMeshers_Regular_1D.cxx b/src/StdMeshers/StdMeshers_Regular_1D.cxx index 9d051ae1a..b2a507433 100644 --- a/src/StdMeshers/StdMeshers_Regular_1D.cxx +++ b/src/StdMeshers/StdMeshers_Regular_1D.cxx @@ -146,7 +146,9 @@ bool StdMeshers_Regular_1D::CheckHypothesis const StdMeshers_LocalLength * hyp = dynamic_cast (theHyp); ASSERT(hyp); - _value[ BEG_LENGTH_IND ] = _value[ END_LENGTH_IND ] = hyp->GetLength(); + //_value[ BEG_LENGTH_IND ] = _value[ END_LENGTH_IND ] = hyp->GetLength(); + _value[ BEG_LENGTH_IND ] = hyp->GetLength(); + _value[ END_LENGTH_IND ] = hyp->GetPrecision(); ASSERT( _value[ BEG_LENGTH_IND ] > 0 ); _hypType = LOCAL_LENGTH; aStatus = SMESH_Hypothesis::HYP_OK; @@ -224,7 +226,9 @@ bool StdMeshers_Regular_1D::CheckHypothesis StdMeshers_AutomaticLength * hyp = const_cast (dynamic_cast (theHyp)); ASSERT(hyp); - _value[ BEG_LENGTH_IND ] = _value[ END_LENGTH_IND ] = hyp->GetLength( &aMesh, aShape ); + //_value[ BEG_LENGTH_IND ] = _value[ END_LENGTH_IND ] = hyp->GetLength( &aMesh, aShape ); + _value[ BEG_LENGTH_IND ] = hyp->GetLength( &aMesh, aShape ); + _value[ END_LENGTH_IND ] = Precision::Confusion(); // ?? or set to zero? ASSERT( _value[ BEG_LENGTH_IND ] > 0 ); _hypType = LOCAL_LENGTH; aStatus = SMESH_Hypothesis::HYP_OK; @@ -243,7 +247,7 @@ static bool computeParamByFunc(Adaptor3d_Curve& C3d, double first, double last, // never do this way //OSD::SetSignal( true ); - if( nbSeg<=0 ) + if (nbSeg <= 0) return false; MESSAGE( "computeParamByFunc" ); @@ -251,12 +255,12 @@ static bool computeParamByFunc(Adaptor3d_Curve& C3d, double first, double last, int nbPnt = 1 + nbSeg; vector x(nbPnt, 0.); - if( !buildDistribution( func, 0.0, 1.0, nbSeg, x, 1E-4 ) ) + if (!buildDistribution(func, 0.0, 1.0, nbSeg, x, 1E-4)) return false; MESSAGE( "Points:\n" ); char buf[1024]; - for( int i=0; i<=nbSeg; i++ ) + for ( int i=0; i<=nbSeg; i++ ) { sprintf( buf, "%f\n", float(x[i] ) ); MESSAGE( buf ); @@ -524,7 +528,7 @@ void StdMeshers_Regular_1D::redistributeNearVertices (SMESH_Mesh & theM std::swap( algo._value[ BEG_LENGTH_IND ], algo._value[ END_LENGTH_IND ]); } list params; - if ( algo.computeInternalParameters( theC3d, L, from, to, params, false )) + if ( algo.computeInternalParameters( theMesh, theC3d, L, from, to, params, false )) { if ( isEnd1 ) params.reverse(); while ( 1 + nHalf-- ) @@ -547,12 +551,14 @@ void StdMeshers_Regular_1D::redistributeNearVertices (SMESH_Mesh & theM * */ //============================================================================= -bool StdMeshers_Regular_1D::computeInternalParameters(Adaptor3d_Curve& theC3d, +bool StdMeshers_Regular_1D::computeInternalParameters(SMESH_Mesh & theMesh, + Adaptor3d_Curve& theC3d, double theLength, double theFirstU, double theLastU, list & theParams, - const bool theReverse) + const bool theReverse, + bool theConsiderPropagation) { theParams.clear(); @@ -568,6 +574,38 @@ bool StdMeshers_Regular_1D::computeInternalParameters(Adaptor3d_Curve& theC3d, { // Local Length hypothesis double nbseg = ceil(theLength / _value[ BEG_LENGTH_IND ]); // integer sup + + // NPAL17873: + bool isFound = false; + if (theConsiderPropagation && !_mainEdge.IsNull()) // propagated from some other edge + { + // Advanced processing to assure equal number of segments in case of Propagation + SMESH_subMesh* sm = theMesh.GetSubMeshContaining(_mainEdge); + if (sm) { + bool computed = sm->IsMeshComputed(); + if (!computed) { + if (sm->GetComputeState() == SMESH_subMesh::READY_TO_COMPUTE) { + sm->ComputeStateEngine(SMESH_subMesh::COMPUTE); + computed = sm->IsMeshComputed(); + } + } + if (computed) { + SMESHDS_SubMesh* smds = sm->GetSubMeshDS(); + int nb_segments = smds->NbElements(); + if (nbseg - 1 <= nb_segments && nb_segments <= nbseg + 1) { + isFound = true; + nbseg = nb_segments; + } + } + } + } + if (!isFound) // not found by meshed edge in the propagation chain, use precision + { + double aPrecision = _value[ END_LENGTH_IND ]; + double nbseg_prec = ceil((theLength / _value[ BEG_LENGTH_IND ]) - aPrecision); + if (nbseg_prec == (nbseg - 1)) nbseg--; + } + if (nbseg <= 0) nbseg = 1; // degenerated edge eltSize = theLength / nbseg; @@ -736,14 +774,14 @@ bool StdMeshers_Regular_1D::computeInternalParameters(Adaptor3d_Curve& theC3d, */ //============================================================================= -bool StdMeshers_Regular_1D::Compute(SMESH_Mesh & aMesh, const TopoDS_Shape & aShape) +bool StdMeshers_Regular_1D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape & theShape) { if ( _hypType == NONE ) return false; - SMESHDS_Mesh * meshDS = aMesh.GetMeshDS(); + SMESHDS_Mesh * meshDS = theMesh.GetMeshDS(); - const TopoDS_Edge & EE = TopoDS::Edge(aShape); + const TopoDS_Edge & EE = TopoDS::Edge(theShape); TopoDS_Edge E = TopoDS::Edge(EE.Oriented(TopAbs_FORWARD)); int shapeID = meshDS->ShapeToIndex( E ); @@ -769,10 +807,10 @@ bool StdMeshers_Regular_1D::Compute(SMESH_Mesh & aMesh, const TopoDS_Shape & aSh BRepAdaptor_Curve C3d( E ); double length = EdgeLength( E ); - if ( ! computeInternalParameters( C3d, length, f, l, params, reversed )) { + if ( ! computeInternalParameters( theMesh, C3d, length, f, l, params, reversed, true )) { return false; } - redistributeNearVertices( aMesh, C3d, length, params, VFirst, VLast ); + redistributeNearVertices( theMesh, C3d, length, params, VFirst, VLast ); // edge extrema (indexes : 1 & NbPoints) already in SMDS (TopoDS_Vertex) // only internal nodes receive an edge position with param on curve diff --git a/src/StdMeshers/StdMeshers_Regular_1D.hxx b/src/StdMeshers/StdMeshers_Regular_1D.hxx index 269a034bc..4a22e7253 100644 --- a/src/StdMeshers/StdMeshers_Regular_1D.hxx +++ b/src/StdMeshers/StdMeshers_Regular_1D.hxx @@ -73,12 +73,14 @@ public: protected: - virtual bool computeInternalParameters (Adaptor3d_Curve & theC3d, - double theLength, - double theFirstU, - double theLastU, - std::list< double > & theParameters, - const bool theReverse); + virtual bool computeInternalParameters (SMESH_Mesh & theMesh, + Adaptor3d_Curve & theC3d, + double theLength, + double theFirstU, + double theLastU, + std::list & theParameters, + const bool theReverse, + bool theConsiderPropagation = false); virtual void redistributeNearVertices (SMESH_Mesh & theMesh, Adaptor3d_Curve & theC3d, @@ -101,12 +103,12 @@ protected: BEG_LENGTH_IND = 0, END_LENGTH_IND = 1, DEFLECTION_IND = 0 - }; + }; enum IValueIndex { NB_SEGMENTS_IND = 0, DISTR_TYPE_IND = 1, - CONV_MODE_IND = 2 + CONV_MODE_IND = 2 }; enum VValueIndex { diff --git a/src/StdMeshersGUI/StdMeshersGUI_StdHypothesisCreator.cxx b/src/StdMeshersGUI/StdMeshersGUI_StdHypothesisCreator.cxx index 13f87e560..bfe80bdf1 100644 --- a/src/StdMeshersGUI/StdMeshersGUI_StdHypothesisCreator.cxx +++ b/src/StdMeshersGUI/StdMeshersGUI_StdHypothesisCreator.cxx @@ -388,6 +388,7 @@ QString StdMeshersGUI_StdHypothesisCreator::storeParams() const StdMeshers::StdMeshers_LocalLength::_narrow( hypothesis() ); h->SetLength( params[0].myValue.toDouble() ); + h->SetPrecision( params[1].myValue.toDouble() ); } else if( hypType()=="SegmentLengthAroundVertex" ) { @@ -535,6 +536,9 @@ bool StdMeshersGUI_StdHypothesisCreator::stdParams( ListOfStdParams& p ) const item.myName = tr("SMESH_LOCAL_LENGTH_PARAM"); item.myValue = h->GetLength(); p.append( item ); + item.myName = tr("SMESH_LOCAL_LENGTH_PRECISION"); + item.myValue = h->GetPrecision(); + p.append( item ); } else if( hypType()=="SegmentLengthAroundVertex" ) { @@ -704,12 +708,15 @@ bool StdMeshersGUI_StdHypothesisCreator::stdParams( ListOfStdParams& p ) const */ //================================================================================ -void StdMeshersGUI_StdHypothesisCreator::attuneStdWidget( QWidget* w, const int ) const +void StdMeshersGUI_StdHypothesisCreator::attuneStdWidget (QWidget* w, const int param) const { SMESHGUI_SpinBox* sb = w->inherits( "SMESHGUI_SpinBox" ) ? ( SMESHGUI_SpinBox* )w : 0; if( hypType()=="LocalLength" && sb ) { - sb->RangeStepAndValidator( VALUE_SMALL, VALUE_MAX, 1.0, 6 ); + if (param == 0) // Length + sb->RangeStepAndValidator( VALUE_SMALL, VALUE_MAX, 1.0, 6 ); + else // Precision + sb->RangeStepAndValidator( 0.0, 1.0, 0.05, 6 ); } else if( hypType()=="Arithmetic1D" && sb ) { diff --git a/src/StdMeshersGUI/StdMeshers_msg_en.po b/src/StdMeshersGUI/StdMeshers_msg_en.po index 922119e14..d8cfd1435 100644 --- a/src/StdMeshersGUI/StdMeshers_msg_en.po +++ b/src/StdMeshersGUI/StdMeshers_msg_en.po @@ -37,6 +37,9 @@ msgstr "Average Length" msgid "SMESH_LOCAL_LENGTH_PARAM" msgstr "Length" +msgid "SMESH_LOCAL_LENGTH_PRECISION" +msgstr "Precision" + msgid "SMESH_LOCAL_LENGTH_TITLE" msgstr "Hypothesis Construction" diff --git a/src/StdMeshers_I/StdMeshers_LocalLength_i.cxx b/src/StdMeshers_I/StdMeshers_LocalLength_i.cxx index 4190c1b8c..f360400ca 100644 --- a/src/StdMeshers_I/StdMeshers_LocalLength_i.cxx +++ b/src/StdMeshers_I/StdMeshers_LocalLength_i.cxx @@ -48,15 +48,15 @@ using namespace std; //============================================================================= StdMeshers_LocalLength_i::StdMeshers_LocalLength_i( PortableServer::POA_ptr thePOA, - int theStudyId, - ::SMESH_Gen* theGenImpl ) - : SALOME::GenericObj_i( thePOA ), + int theStudyId, + ::SMESH_Gen* theGenImpl ) + : SALOME::GenericObj_i( thePOA ), SMESH_Hypothesis_i( thePOA ) { MESSAGE( "StdMeshers_LocalLength_i::StdMeshers_LocalLength_i" ); myBaseImpl = new ::StdMeshers_LocalLength( theGenImpl->GetANewId(), - theStudyId, - theGenImpl ); + theStudyId, + theGenImpl ); } //============================================================================= @@ -79,7 +79,6 @@ StdMeshers_LocalLength_i::~StdMeshers_LocalLength_i() * Set length */ //============================================================================= - void StdMeshers_LocalLength_i::SetLength( CORBA::Double theLength ) throw ( SALOME::SALOME_Exception ) { @@ -97,6 +96,30 @@ void StdMeshers_LocalLength_i::SetLength( CORBA::Double theLength ) SMESH::TPythonDump() << _this() << ".SetLength( " << theLength << " )"; } +//============================================================================= +/*! + * StdMeshers_LocalLength_i::SetPrecision + * + * Set length + */ +//============================================================================= +void StdMeshers_LocalLength_i::SetPrecision( CORBA::Double thePrecision ) + throw ( SALOME::SALOME_Exception ) +{ + MESSAGE( "StdMeshers_LocalLength_i::SetPrecision" ); + ASSERT( myBaseImpl ); + try { + this->GetImpl()->SetPrecision( thePrecision ); + } + catch ( SALOME_Exception& S_ex ) { + THROW_SALOME_CORBA_EXCEPTION( S_ex.what(), + SALOME::BAD_PARAM ); + } + + // Update Python script + SMESH::TPythonDump() << _this() << ".SetPrecision( " << thePrecision << " )"; +} + //============================================================================= /*! * StdMeshers_LocalLength_i::GetLength @@ -104,7 +127,6 @@ void StdMeshers_LocalLength_i::SetLength( CORBA::Double theLength ) * Get length */ //============================================================================= - CORBA::Double StdMeshers_LocalLength_i::GetLength() { MESSAGE( "StdMeshers_LocalLength_i::GetLength" ); @@ -112,6 +134,20 @@ CORBA::Double StdMeshers_LocalLength_i::GetLength() return this->GetImpl()->GetLength(); } +//============================================================================= +/*! + * StdMeshers_LocalLength_i::GetPrecision + * + * Get precision + */ +//============================================================================= +CORBA::Double StdMeshers_LocalLength_i::GetPrecision() +{ + MESSAGE( "StdMeshers_LocalLength_i::GetPrecision" ); + ASSERT( myBaseImpl ); + return this->GetImpl()->GetPrecision(); +} + //============================================================================= /*! * StdMeshers_LocalLength_i::GetImpl @@ -119,7 +155,6 @@ CORBA::Double StdMeshers_LocalLength_i::GetLength() * Get implementation */ //============================================================================= - ::StdMeshers_LocalLength* StdMeshers_LocalLength_i::GetImpl() { MESSAGE( "StdMeshers_LocalLength_i::GetImpl" ); @@ -139,4 +174,3 @@ CORBA::Boolean StdMeshers_LocalLength_i::IsDimSupported( SMESH::Dimension type ) { return type == SMESH::DIM_1D; } - diff --git a/src/StdMeshers_I/StdMeshers_LocalLength_i.hxx b/src/StdMeshers_I/StdMeshers_LocalLength_i.hxx index 62c73cbc8..a51997803 100644 --- a/src/StdMeshers_I/StdMeshers_LocalLength_i.hxx +++ b/src/StdMeshers_I/StdMeshers_LocalLength_i.hxx @@ -58,8 +58,14 @@ public: // Set length void SetLength( CORBA::Double theLength ) throw ( SALOME::SALOME_Exception ); + // Set precision + void SetPrecision( CORBA::Double thePrecision ) + throw ( SALOME::SALOME_Exception ); + // Get length CORBA::Double GetLength(); + // Get precision + CORBA::Double GetPrecision(); // Get implementation ::StdMeshers_LocalLength* GetImpl();