]> SALOME platform Git repositories - modules/smesh.git/commitdiff
Salome HOME
Add “Grading” parameter to Adaptive 1D hypothesis
authorimn <imn@opencascade.com>
Wed, 16 Jul 2014 07:45:56 +0000 (11:45 +0400)
committerimn <imn@opencascade.com>
Wed, 16 Jul 2014 12:43:20 +0000 (16:43 +0400)
14 files changed:
doc/salome/examples/defining_hypotheses_adaptive1d.py
doc/salome/gui/SMESH/images/adaptive1d.png
doc/salome/gui/SMESH/input/1d_meshing_hypo.doc
idl/SMESH_BasicHypothesis.idl
resources/StdMeshers.xml.in
src/SMESH_SWIG/StdMeshersBuilder.py
src/StdMeshers/StdMeshers_Adaptive1D.cxx
src/StdMeshers/StdMeshers_Adaptive1D.hxx
src/StdMeshersGUI/StdMeshersGUI_StdHypothesisCreator.cxx
src/StdMeshersGUI/StdMeshers_msg_en.ts
src/StdMeshersGUI/StdMeshers_msg_fr.ts
src/StdMeshersGUI/StdMeshers_msg_ja.ts
src/StdMeshers_I/StdMeshers_Adaptive1D_i.cxx
src/StdMeshers_I/StdMeshers_Adaptive1D_i.hxx

index 96e2fabfe1912ce45ff7bcdad8640cd474d92771..0110ebc886acc62a5ce424ee58eddfd7726c5e96 100644 (file)
@@ -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()
 
index 8091c8d0ed7e8da41533a57a36ea75a44e36cd2b..dcabeee0fafdef4545744b1cc09fbe4482ab54ee 100644 (file)
Binary files a/doc/salome/gui/SMESH/images/adaptive1d.png and b/doc/salome/gui/SMESH/images/adaptive1d.png differ
index 1228865f7a643638333442e018cd0ee2646eeba1..d11a8baa04d3df70f38238241399c8713acf9961 100644 (file)
@@ -31,8 +31,8 @@ creation of narrow 2D elements.
 
 - <b>Min size</b> parameter limits the minimal segment size. 
 - <b>Max size</b> parameter defines the length of segments on straight edges. 
-- \b Deflection parameter gives maximal distance of a segment from a curved edge.
-
+- <b>Deflection</b> parameter gives maximal distance of a segment from a curved edge.
+- <b>Grading</b> 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"
 
 <b>See Also</b> a \ref tui_1d_adaptive "sample TUI Script" that uses Adaptive hypothesis.
index 72d1b1a2d2817cd9e7cd4bc128332ffe3a182ff2..a8f9b31b1313fb0b580f6f3ea7f4041b0da732be 100644 (file)
@@ -433,6 +433,13 @@ module StdMeshers
      */
     void SetDeflection(in double deflection) raises (SALOME::SALOME_Exception);
     double GetDeflection();
+    
+    /*!
+     * Sets <grading> parameter value, 
+     * i.e. how much size of adjacent elements can differ
+     */
+    void SetGrading(in double grading) raises (SALOME::SALOME_Exception);
+    double GetGrading();
   };
 
   /*!
index 239ff2cb4919b3e4cdce58a31602899989bb665e..cf8146f54948fb0901589684bb9db9a722ffa234 100644 (file)
         <hypo>GeometricProgression=GeometricProgression(SetStartLength(),SetCommonRatio(),SetReversedEdges())</hypo>
         <hypo>StartEndLength=StartEndLength(SetStartLength(),SetEndLength(),SetReversedEdges())</hypo>
         <hypo>Deflection1D=Deflection1D(SetDeflection())</hypo>
-        <hypo>Adaptive1D=Adaptive(SetMinSize(),SetMaxSize(),SetDeflection())</hypo>
+        <hypo>Adaptive1D=Adaptive(SetMinSize(),SetMaxSize(),SetDeflection(),SetGrading())</hypo>
         <hypo>AutomaticLength=AutomaticLength(SetFineness())</hypo>
         <hypo>FixedPoints1D=FixedPoints1D(SetPoints(),SetNbSegments(),SetReversedEdges())</hypo>
         <hypo>Propagation=Propagation()</hypo>
         <hypo>GeometricProgression=GeometricProgression(SetStartLength(),SetCommonRatio(),SetReversedEdges())</hypo>
         <hypo>StartEndLength=StartEndLength(SetStartLength(),SetEndLength(),SetReversedEdges())</hypo>
         <hypo>Deflection1D=Deflection1D(SetDeflection())</hypo>
-        <hypo>Adaptive1D=Adaptive(SetMinSize(),SetMaxSize(),SetDeflection())</hypo>
+        <hypo>Adaptive1D=Adaptive(SetMinSize(),SetMaxSize(),SetDeflection(),SetGrading())</hypo>
         <hypo>AutomaticLength=AutomaticLength(SetFineness())</hypo>
         <hypo>FixedPoints1D=FixedPoints1D(SetPoints(),SetNbSegments(),SetReversedEdges())</hypo>
         <hypo>Propagation=Propagation()</hypo>
index 1e9c735332a99c79bbfadde90cb2a50c2a3eeb7f..054bd7e7e11b96652aa6515650e35f4c8e61e308 100644 (file)
@@ -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
index ecb30b1aa86f3821e5655b71fee58aa1abf382d9..55d5d00b48efb02fd29c996e0ca3910b5bc96499 100644 (file)
@@ -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<double>::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 );
index a6211808737ba44051c5eef788470c62e75cfe0e..200465026d651a1a0d9d26135507773a4bd558bf 100644 (file)
@@ -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 <grading> 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
 };
 
index 891b2c0c86170a4fdc62f383ba55a3eb9d44303d..904f0bbd73e2a63deeb85b9d5444a1dda4d51e9d 100644 (file)
@@ -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" ))
     {
index 08dcd56718567ff24e3e5fafeca25a544978d25e..abf15afdb1dc8730ad1460d3d9f156c115149e43 100644 (file)
         <source>SMESH_INVALID_FUNCTION</source>
         <translation>Function is invalid</translation>
     </message>
+    <message>
+        <source>SMESH_GRADING1D_PARAM</source>
+        <translation>Grading</translation>
+    </message>
     <message>
         <source>SMESH_LAYERS_DISTRIBUTION</source>
         <translation>1D Hypothesis</translation>
index ef19f469fe10b7de6d2778d0e37de149acceb31e..2a434ca5be90f51303fd2fd798e9eb0774325589 100755 (executable)
         <source>SMESH_INVALID_FUNCTION</source>
         <translation>La fonction n&apos;est pas valide</translation>
     </message>
+    <message>
+        <source>SMESH_GRADING1D_PARAM</source>
+        <translation type="unfinished">Grading</translation>
+    </message>
     <message>
         <source>SMESH_LAYERS_DISTRIBUTION</source>
         <translation>Hypothèse 1D </translation>
index 8a90c5e1f61268b7b0f6cedac65e72de2716db34..5e5d155cbbabc5fe76c078a4ee4328d9af1eeae3 100644 (file)
       <source>SMESH_INVALID_FUNCTION</source>
       <translation>関数が無効です。</translation>
     </message>
+    <message>
+        <source>SMESH_GRADING1D_PARAM</source>
+        <translation type="unfinished">Grading</translation>
+    </message>
     <message>
       <source>SMESH_LAYERS_DISTRIBUTION</source>
       <translation>仮説 1 d</translation>
index c6d82b9c11e41d0fe96fc171f40bc2de96dc1637..bcf6cc5ac5a7c8c6b160c5bad976631a245f3ae5 100644 (file)
@@ -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
index 7dd59e708a64d918715fccc168730fe179b0be3b..fada85ba65f1b36381d43ffa3f8c70137778875b 100644 (file)
@@ -70,6 +70,12 @@ class STDMESHERS_I_EXPORT StdMeshers_Adaptive1D_i:
   void SetDeflection( CORBA::Double theLength ) throw (SALOME::SALOME_Exception);
   CORBA::Double GetDeflection();
 
+  /*!
+   * Sets <grading> 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