From: imn Date: Wed, 16 Jul 2014 07:45:56 +0000 (+0400) Subject: Add “Grading” parameter to Adaptive 1D hypothesis X-Git-Tag: V7_5_0a1~56 X-Git-Url: http://git.salome-platform.org/gitweb/?a=commitdiff_plain;h=80fe1ddefc561a7a571ac08807f7d173f45d8080;p=modules%2Fsmesh.git Add “Grading” parameter to Adaptive 1D hypothesis --- diff --git a/doc/salome/examples/defining_hypotheses_adaptive1d.py b/doc/salome/examples/defining_hypotheses_adaptive1d.py index 96e2fabfe..0110ebc88 100644 --- a/doc/salome/examples/defining_hypotheses_adaptive1d.py +++ b/doc/salome/examples/defining_hypotheses_adaptive1d.py @@ -19,12 +19,14 @@ shape = geompy.MakeCut( shape, tool, theName="shape" ) # Parameters of Adaptive hypothesis. minSize and maxSize are such that they do not limit # size of segments because size of geometrical features lies within [2.-100.] range, hence # size of segments is defined by deflection parameter and size of geometrical features only. +# grading is defined how much size of adjacent elements can differ. minSize = 0.1 maxSize = 200 deflection = 0.05 +grading = 0.7 mesh = smesh.Mesh( shape ) -mesh.Segment().Adaptive( minSize, maxSize, deflection ) +mesh.Segment().Adaptive( minSize, maxSize, deflection, grading ) mesh.Triangle().MaxElementArea( 300 ) mesh.Compute() diff --git a/doc/salome/gui/SMESH/images/adaptive1d.png b/doc/salome/gui/SMESH/images/adaptive1d.png index 8091c8d0e..dcabeee0f 100644 Binary files a/doc/salome/gui/SMESH/images/adaptive1d.png and b/doc/salome/gui/SMESH/images/adaptive1d.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 1228865f7..d11a8baa0 100644 --- a/doc/salome/gui/SMESH/input/1d_meshing_hypo.doc +++ b/doc/salome/gui/SMESH/input/1d_meshing_hypo.doc @@ -31,8 +31,8 @@ creation of narrow 2D elements. - Min size parameter limits the minimal segment size. - Max size parameter defines the length of segments on straight edges. -- \b Deflection parameter gives maximal distance of a segment from a curved edge. - +- Deflection parameter gives maximal distance of a segment from a curved edge. +- Grading parameter defines how much size of adjacent elements can differ. \image html adaptive1d_sample_mesh.png "Adaptive hypothesis and Netgen 2D algorithm - the size of mesh segments reflects the size of geometrical features" See Also a \ref tui_1d_adaptive "sample TUI Script" that uses Adaptive hypothesis. diff --git a/idl/SMESH_BasicHypothesis.idl b/idl/SMESH_BasicHypothesis.idl index 72d1b1a2d..a8f9b31b1 100644 --- a/idl/SMESH_BasicHypothesis.idl +++ b/idl/SMESH_BasicHypothesis.idl @@ -433,6 +433,13 @@ module StdMeshers */ void SetDeflection(in double deflection) raises (SALOME::SALOME_Exception); double GetDeflection(); + + /*! + * Sets parameter value, + * i.e. how much size of adjacent elements can differ + */ + void SetGrading(in double grading) raises (SALOME::SALOME_Exception); + double GetGrading(); }; /*! diff --git a/resources/StdMeshers.xml.in b/resources/StdMeshers.xml.in index 239ff2cb4..cf8146f54 100644 --- a/resources/StdMeshers.xml.in +++ b/resources/StdMeshers.xml.in @@ -236,7 +236,7 @@ GeometricProgression=GeometricProgression(SetStartLength(),SetCommonRatio(),SetReversedEdges()) StartEndLength=StartEndLength(SetStartLength(),SetEndLength(),SetReversedEdges()) Deflection1D=Deflection1D(SetDeflection()) - Adaptive1D=Adaptive(SetMinSize(),SetMaxSize(),SetDeflection()) + Adaptive1D=Adaptive(SetMinSize(),SetMaxSize(),SetDeflection(),SetGrading()) AutomaticLength=AutomaticLength(SetFineness()) FixedPoints1D=FixedPoints1D(SetPoints(),SetNbSegments(),SetReversedEdges()) Propagation=Propagation() @@ -262,7 +262,7 @@ GeometricProgression=GeometricProgression(SetStartLength(),SetCommonRatio(),SetReversedEdges()) StartEndLength=StartEndLength(SetStartLength(),SetEndLength(),SetReversedEdges()) Deflection1D=Deflection1D(SetDeflection()) - Adaptive1D=Adaptive(SetMinSize(),SetMaxSize(),SetDeflection()) + Adaptive1D=Adaptive(SetMinSize(),SetMaxSize(),SetDeflection(),SetGrading()) AutomaticLength=AutomaticLength(SetFineness()) FixedPoints1D=FixedPoints1D(SetPoints(),SetNbSegments(),SetReversedEdges()) Propagation=Propagation() diff --git a/src/SMESH_SWIG/StdMeshersBuilder.py b/src/SMESH_SWIG/StdMeshersBuilder.py index 1e9c73533..054bd7e7e 100644 --- a/src/SMESH_SWIG/StdMeshersBuilder.py +++ b/src/SMESH_SWIG/StdMeshersBuilder.py @@ -183,20 +183,23 @@ class StdMeshersBuilder_Segment(Mesh_Algorithm): # @param minSize defines the minimal allowed segment length # @param maxSize defines the maximal allowed segment length # @param deflection defines the maximal allowed distance from a segment to an edge + # @param grading defines how much size of adjacent elements can differ # @param UseExisting if ==true - searches for an existing hypothesis created with # the same parameters, else (default) - creates a new one # @return an instance of StdMeshers_Adaptive1D hypothesis # @ingroup l3_hypos_1dhyps - def Adaptive(self, minSize, maxSize, deflection, UseExisting=False): + def Adaptive(self, minSize, maxSize, deflection, grading, UseExisting=False): from salome.smesh.smeshBuilder import IsEqual compFun = lambda hyp, args: ( IsEqual(hyp.GetMinSize(), args[0]) and \ IsEqual(hyp.GetMaxSize(), args[1]) and \ - IsEqual(hyp.GetDeflection(), args[2])) - hyp = self.Hypothesis("Adaptive1D", [minSize, maxSize, deflection], + IsEqual(hyp.GetDeflection(), args[2]) and \ + IsEqual(hyp.GetGrading(), args[3])) + hyp = self.Hypothesis("Adaptive1D", [minSize, maxSize, deflection, grading], UseExisting=UseExisting, CompareMethod=compFun) hyp.SetMinSize(minSize) hyp.SetMaxSize(maxSize) hyp.SetDeflection(deflection) + hyp.SetGrading(grading) return hyp ## Defines "Arithmetic1D" hypothesis to cut an edge in several segments with a length diff --git a/src/StdMeshers/StdMeshers_Adaptive1D.cxx b/src/StdMeshers/StdMeshers_Adaptive1D.cxx index ecb30b1aa..55d5d00b4 100644 --- a/src/StdMeshers/StdMeshers_Adaptive1D.cxx +++ b/src/StdMeshers/StdMeshers_Adaptive1D.cxx @@ -942,6 +942,7 @@ StdMeshers_Adaptive1D::StdMeshers_Adaptive1D(int hypId, myMinSize = 1e-10; myMaxSize = 1e+10; myDeflection = 1e-2; + myGrading = 1e-2; myAlgo = NULL; _name = "Adaptive1D"; _param_algo_dim = 1; // is used by SMESH_Regular_1D @@ -968,6 +969,20 @@ void StdMeshers_Adaptive1D::SetDeflection(double value) } } //======================================================================= +//function : SetGrading +//purpose : +void StdMeshers_Adaptive1D::SetGrading(double value) + throw(SALOME_Exception) +{ + if (value <= std::numeric_limits::min() ) + throw SALOME_Exception("Grading must be greater that zero"); + if (myGrading != value) + { + myGrading = value; + NotifySubMeshesHypothesisModification(); + } +} +//======================================================================= //function : SetMinSize //purpose : Sets minimal allowed segment length void StdMeshers_Adaptive1D::SetMinSize(double minSize) @@ -1002,7 +1017,7 @@ void StdMeshers_Adaptive1D::SetMaxSize(double maxSize) //purpose : Persistence ostream & StdMeshers_Adaptive1D::SaveTo(ostream & save) { - save << myMinSize << " " << myMaxSize << " " << myDeflection; + save << myMinSize << " " << myMaxSize << " " << myDeflection << " " << myGrading; save << " " << -1 << " " << -1; // preview addition of parameters return save; } @@ -1012,7 +1027,7 @@ ostream & StdMeshers_Adaptive1D::SaveTo(ostream & save) istream & StdMeshers_Adaptive1D::LoadFrom(istream & load) { int dummyParam; - bool isOK = (load >> myMinSize >> myMaxSize >> myDeflection >> dummyParam >> dummyParam); + bool isOK = (load >> myMinSize >> myMaxSize >> myDeflection >> myGrading >> dummyParam >> dummyParam); if (!isOK) load.clear(ios::badbit | load.rdstate()); return load; @@ -1082,6 +1097,7 @@ bool StdMeshers_Adaptive1D::SetParametersByDefaults(const TDefaults& dflts, myMinSize = dflts._elemLength / 10; myMaxSize = dflts._elemLength * 2; myDeflection = myMinSize / 7; + myGrading = 0.7; return true; } @@ -1145,7 +1161,7 @@ bool AdaptiveAlgo::Compute(SMESH_Mesh & theMesh, myMesh = &theMesh; SMESH_MesherHelper helper( theMesh ); - const double grading = 0.7; + const double grading = myHyp->GetGrading(); TopTools_IndexedMapOfShape edgeMap, faceMap; TopExp::MapShapes( theShape, TopAbs_EDGE, edgeMap ); diff --git a/src/StdMeshers/StdMeshers_Adaptive1D.hxx b/src/StdMeshers/StdMeshers_Adaptive1D.hxx index a62118087..200465026 100644 --- a/src/StdMeshers/StdMeshers_Adaptive1D.hxx +++ b/src/StdMeshers/StdMeshers_Adaptive1D.hxx @@ -59,6 +59,13 @@ class STDMESHERS_EXPORT StdMeshers_Adaptive1D : public SMESH_Hypothesis void SetDeflection(double value) throw(SALOME_Exception); double GetDeflection() const { return myDeflection; } + /*! + * Sets parameter value, + * i.e. how much size of adjacent elements can differ + */ + void SetGrading(double value) throw(SALOME_Exception); + double GetGrading() const { return myGrading; } + virtual std::ostream & SaveTo(std::ostream & save); virtual std::istream & LoadFrom(std::istream & load); @@ -83,7 +90,7 @@ class STDMESHERS_EXPORT StdMeshers_Adaptive1D : public SMESH_Hypothesis protected: - double myMinSize, myMaxSize, myDeflection; + double myMinSize, myMaxSize, myDeflection, myGrading; SMESH_Algo* myAlgo; // StdMeshers_AdaptiveAlgo_1D implemented in cxx file }; diff --git a/src/StdMeshersGUI/StdMeshersGUI_StdHypothesisCreator.cxx b/src/StdMeshersGUI/StdMeshersGUI_StdHypothesisCreator.cxx index 891b2c0c8..904f0bbd7 100644 --- a/src/StdMeshersGUI/StdMeshersGUI_StdHypothesisCreator.cxx +++ b/src/StdMeshersGUI/StdMeshersGUI_StdHypothesisCreator.cxx @@ -609,6 +609,8 @@ QString StdMeshersGUI_StdHypothesisCreator::storeParams() const h->SetMaxSize( params[1].myValue.toDouble() ); h->SetVarParameter( params[0].text(), "SetDeflection" ); h->SetDeflection( params[2].myValue.toDouble() ); + h->SetVarParameter( params[0].text(), "SetGrading" ); + h->SetGrading( params[3].myValue.toDouble() ); } else if( hypType()=="AutomaticLength" ) { @@ -1048,6 +1050,11 @@ bool StdMeshersGUI_StdHypothesisCreator::stdParams( ListOfStdParams& p ) const if(!initVariableName( hyp, item, "SetDeflection" )) item.myValue = h->GetDeflection(); p.append( item ); + + item.myName = tr( "SMESH_GRADING1D_PARAM" ); + if(!initVariableName( hyp, item, "SetGrading" )) + item.myValue = h->GetGrading(); + p.append( item ); } else if( hypType()=="AutomaticLength" ) { @@ -1421,7 +1428,10 @@ void StdMeshersGUI_StdHypothesisCreator::attuneStdWidget (QWidget* w, const int) } else if( hypType()=="Adaptive1D" ) { - sb->RangeStepAndValidator( VALUE_SMALL, VALUE_MAX, 1.0, "length_precision" ); + if (sb->objectName() == tr("SMESH_GRADING1D_PARAM")) + sb->RangeStepAndValidator( 0.0, 2.0, 0.1, "length_precision" ); + else + sb->RangeStepAndValidator( VALUE_SMALL, VALUE_MAX, 1.0, "length_precision" ); } else if( hypType().startsWith( "ViscousLayers" )) { diff --git a/src/StdMeshersGUI/StdMeshers_msg_en.ts b/src/StdMeshersGUI/StdMeshers_msg_en.ts index 08dcd5671..abf15afdb 100644 --- a/src/StdMeshersGUI/StdMeshers_msg_en.ts +++ b/src/StdMeshersGUI/StdMeshers_msg_en.ts @@ -150,6 +150,10 @@ SMESH_INVALID_FUNCTION Function is invalid + + SMESH_GRADING1D_PARAM + Grading + SMESH_LAYERS_DISTRIBUTION 1D Hypothesis diff --git a/src/StdMeshersGUI/StdMeshers_msg_fr.ts b/src/StdMeshersGUI/StdMeshers_msg_fr.ts index ef19f469f..2a434ca5b 100755 --- a/src/StdMeshersGUI/StdMeshers_msg_fr.ts +++ b/src/StdMeshersGUI/StdMeshers_msg_fr.ts @@ -135,6 +135,10 @@ SMESH_INVALID_FUNCTION La fonction n'est pas valide + + SMESH_GRADING1D_PARAM + Grading + SMESH_LAYERS_DISTRIBUTION Hypothèse 1D diff --git a/src/StdMeshersGUI/StdMeshers_msg_ja.ts b/src/StdMeshersGUI/StdMeshers_msg_ja.ts index 8a90c5e1f..5e5d155cb 100644 --- a/src/StdMeshersGUI/StdMeshers_msg_ja.ts +++ b/src/StdMeshersGUI/StdMeshers_msg_ja.ts @@ -150,6 +150,10 @@ SMESH_INVALID_FUNCTION 関数が無効です。 + + SMESH_GRADING1D_PARAM + Grading + SMESH_LAYERS_DISTRIBUTION 仮説 1 d diff --git a/src/StdMeshers_I/StdMeshers_Adaptive1D_i.cxx b/src/StdMeshers_I/StdMeshers_Adaptive1D_i.cxx index c6d82b9c1..bcf6cc5ac 100644 --- a/src/StdMeshers_I/StdMeshers_Adaptive1D_i.cxx +++ b/src/StdMeshers_I/StdMeshers_Adaptive1D_i.cxx @@ -152,6 +152,37 @@ CORBA::Double StdMeshers_Adaptive1D_i::GetDeflection() return this->GetImpl()->GetDeflection(); } +//======================================================================= +//function : SetGrading +//purpose : Sets how much size of adjacent elements can differ. +//======================================================================= + +void StdMeshers_Adaptive1D_i::SetGrading( CORBA::Double theValue ) + throw ( SALOME::SALOME_Exception ) +{ + ASSERT( myBaseImpl ); + try { + this->GetImpl()->SetGrading( theValue ); + } + catch ( SALOME_Exception& S_ex ) { + THROW_SALOME_CORBA_EXCEPTION( S_ex.what(), SALOME::BAD_PARAM ); + } + + // Update Python script + SMESH::TPythonDump() << _this() << ".SetGrading( " << SMESH::TVar(theValue) << " )"; +} + +//======================================================================= +//function : GetGrading +//purpose : Returns grading +//======================================================================= + +CORBA::Double StdMeshers_Adaptive1D_i::GetGrading() +{ + ASSERT( myBaseImpl ); + return this->GetImpl()->GetGrading(); +} + //======================================================================= //function : GetImpl //purpose : Get implementation diff --git a/src/StdMeshers_I/StdMeshers_Adaptive1D_i.hxx b/src/StdMeshers_I/StdMeshers_Adaptive1D_i.hxx index 7dd59e708..fada85ba6 100644 --- a/src/StdMeshers_I/StdMeshers_Adaptive1D_i.hxx +++ b/src/StdMeshers_I/StdMeshers_Adaptive1D_i.hxx @@ -70,6 +70,12 @@ class STDMESHERS_I_EXPORT StdMeshers_Adaptive1D_i: void SetDeflection( CORBA::Double theLength ) throw (SALOME::SALOME_Exception); CORBA::Double GetDeflection(); + /*! + * Sets parameter value, + // * i.e. a maximal allowed distance between a segment and an edge. + */ + void SetGrading( CORBA::Double theLength ) throw (SALOME::SALOME_Exception); + CORBA::Double GetGrading(); /*! * Returns implementation