Salome HOME
Implementation of Quadrangle (Mapping) for faces built on 3 edges (0018911 from Mantis).
[modules/smesh.git] / src / StdMeshers / StdMeshers_NumberOfSegments.cxx
index 0b9d333556f0052ec19e485c38d81457adbcbfe6..e03a921e85db59d00dd55e409503ec1082e36980 100644 (file)
@@ -1,46 +1,55 @@
-//  SMESH SMESH : implementaion of SMESH idl descriptions
+//  Copyright (C) 2007-2008  CEA/DEN, EDF R&D, OPEN CASCADE
+//
+//  Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
+//  CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
+//
+//  This library is free software; you can redistribute it and/or
+//  modify it under the terms of the GNU Lesser General Public
+//  License as published by the Free Software Foundation; either
+//  version 2.1 of the License.
 //
-//  Copyright (C) 2003  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
-//  CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS 
-// 
-//  This library is free software; you can redistribute it and/or 
-//  modify it under the terms of the GNU Lesser General Public 
-//  License as published by the Free Software Foundation; either 
-//  version 2.1 of the License. 
-// 
-//  This library is distributed in the hope that it will be useful, 
-//  but WITHOUT ANY WARRANTY; without even the implied warranty of 
-//  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU 
-//  Lesser General Public License for more details. 
-// 
-//  You should have received a copy of the GNU Lesser General Public 
-//  License along with this library; if not, write to the Free Software 
-//  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA 
-// 
-//  See http://www.opencascade.org/SALOME/ or email : webmaster.salome@opencascade.org 
+//  This library is distributed in the hope that it will be useful,
+//  but WITHOUT ANY WARRANTY; without even the implied warranty of
+//  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+//  Lesser General Public License for more details.
 //
+//  You should have received a copy of the GNU Lesser General Public
+//  License along with this library; if not, write to the Free Software
+//  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
 //
+//  See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
 //
+//  SMESH SMESH : implementaion of SMESH idl descriptions
 //  File   : StdMeshers_NumberOfSegments.cxx
 //           Moved here from SMESH_NumberOfSegments.cxx
 //  Author : Paul RASCLE, EDF
 //  Module : SMESH
-//  $Header$
-
-using namespace std;
+//
 #include "StdMeshers_NumberOfSegments.hxx"
+
 #include "StdMeshers_Distribution.hxx"
-#include <Standard_ErrorHandler.hxx>
-#include <TCollection_AsciiString.hxx>
+#include "SMESHDS_SubMesh.hxx"
+#include "SMESH_Mesh.hxx"
+
 #include <ExprIntrp_GenExp.hxx>
-#include <Expr_NamedUnknown.hxx>
-#include <CASCatch_CatchSignals.hxx>
-#include <CASCatch_Failure.hxx> 
-#include <CASCatch_ErrorHandler.hxx>
-#include <OSD.hxx>
 #include <Expr_Array1OfNamedUnknown.hxx>
+#include <Expr_NamedUnknown.hxx>
 #include <TColStd_Array1OfReal.hxx>
+#include <TCollection_AsciiString.hxx>
+#include <TopExp.hxx>
+#include <TopTools_IndexedMapOfShape.hxx>
+
+#if (OCC_VERSION_MAJOR << 16 | OCC_VERSION_MINOR << 8 | OCC_VERSION_MAINTENANCE) > 0x060100
+#define NO_CAS_CATCH
+#endif
+
+#include <Standard_Failure.hxx>
+
+#ifdef NO_CAS_CATCH
+#include <Standard_ErrorHandler.hxx>
+#endif
 
+using namespace std;
 
 const double PRECISION = 1e-7;
 
@@ -50,10 +59,11 @@ const double PRECISION = 1e-7;
  */
 //=============================================================================
 
-StdMeshers_NumberOfSegments::StdMeshers_NumberOfSegments(int hypId, int studyId,
-       SMESH_Gen * gen)
+StdMeshers_NumberOfSegments::StdMeshers_NumberOfSegments(int         hypId,
+                                                         int         studyId,
+                                                         SMESH_Gen * gen)
   : SMESH_Hypothesis(hypId, studyId, gen),
-    _numberOfSegments(1),
+    _numberOfSegments(15),//issue 19923
     _distrType(DT_Regular),
     _scaleFactor(1.),
     _convMode(1)  //cut negative by default
@@ -77,17 +87,20 @@ StdMeshers_NumberOfSegments::~StdMeshers_NumberOfSegments()
  *  
  */
 //=============================================================================
-const std::vector<double>& StdMeshers_NumberOfSegments::BuildDistributionExpr( const char* expr, int nbSeg, int conv )
-throw ( SALOME_Exception )
+const vector<double>&
+StdMeshers_NumberOfSegments::BuildDistributionExpr( const char* expr,int nbSeg,int conv )
+  throw ( SALOME_Exception )
 {
   if( !buildDistribution( TCollection_AsciiString( ( Standard_CString )expr ), conv, 0.0, 1.0, nbSeg, _distr, 1E-4 ) )
     _distr.resize( 0 );
   return _distr;
 }
 
-const std::vector<double>& StdMeshers_NumberOfSegments::BuildDistributionTab( const std::vector<double>& tab,
-                                                                             int nbSeg, int conv )
-throw ( SALOME_Exception )
+const vector<double>&
+StdMeshers_NumberOfSegments::BuildDistributionTab( const vector<double>& tab,
+                                                   int nbSeg,
+                                                   int conv )
+  throw ( SALOME_Exception )
 {
   if( !buildDistribution( tab, conv, 0.0, 1.0, nbSeg, _distr, 1E-4 ) )
     _distr.resize( 0 );
@@ -103,14 +116,13 @@ throw ( SALOME_Exception )
 void StdMeshers_NumberOfSegments::SetNumberOfSegments(int segmentsNumber)
 throw(SALOME_Exception)
 {
-       int oldNumberOfSegments = _numberOfSegments;
-       if (segmentsNumber <= 0)
-               throw
-                       SALOME_Exception(LOCALIZED("number of segments must be positive"));
-       _numberOfSegments = segmentsNumber;
-
-       if (oldNumberOfSegments != _numberOfSegments)
-               NotifySubMeshesHypothesisModification();
+  int oldNumberOfSegments = _numberOfSegments;
+  if (segmentsNumber <= 0)
+    throw SALOME_Exception(LOCALIZED("number of segments must be positive"));
+  _numberOfSegments = segmentsNumber;
+
+  if (oldNumberOfSegments != _numberOfSegments)
+    NotifySubMeshesHypothesisModification();
 }
 
 //=============================================================================
@@ -121,7 +133,7 @@ throw(SALOME_Exception)
 
 int StdMeshers_NumberOfSegments::GetNumberOfSegments() const
 {
-       return _numberOfSegments;
+  return _numberOfSegments;
 }
 
 //================================================================================
@@ -164,11 +176,12 @@ void StdMeshers_NumberOfSegments::SetScaleFactor(double scaleFactor)
   throw(SALOME_Exception)
 {
   if (_distrType != DT_Scale)
-    throw SALOME_Exception(LOCALIZED("not a scale distribution"));
+    _distrType = DT_Scale;
+    //throw SALOME_Exception(LOCALIZED("not a scale distribution"));
   if (scaleFactor < PRECISION)
     throw SALOME_Exception(LOCALIZED("scale factor must be positive"));
-  if (fabs(scaleFactor - 1.0) < PRECISION)
-    throw SALOME_Exception(LOCALIZED("scale factor must not be equal to 1"));
+  //if (fabs(scaleFactor - 1.0) < PRECISION)
+  //  throw SALOME_Exception(LOCALIZED("scale factor must not be equal to 1"));
 
   if (fabs(_scaleFactor - scaleFactor) > PRECISION)
   {
@@ -197,11 +210,12 @@ double StdMeshers_NumberOfSegments::GetScaleFactor() const
  */
 //================================================================================
 
-void StdMeshers_NumberOfSegments::SetTableFunction(const std::vector<double>& table)
+void StdMeshers_NumberOfSegments::SetTableFunction(const vector<double>& table)
   throw(SALOME_Exception)
 {
   if (_distrType != DT_TabFunc)
-    throw SALOME_Exception(LOCALIZED("not a table function distribution"));
+    _distrType = DT_TabFunc;
+  //throw SALOME_Exception(LOCALIZED("not a table function distribution"));
   if ( (table.size() % 2) != 0 )
     throw SALOME_Exception(LOCALIZED("odd size of vector of table function"));
 
@@ -209,24 +223,19 @@ void StdMeshers_NumberOfSegments::SetTableFunction(const std::vector<double>& ta
   double prev = -PRECISION;
   bool isSame = table.size() == _table.size();
 
-  OSD::SetSignal( true );
-  CASCatch_CatchSignals aCatchSignals;
-  aCatchSignals.Activate();
-
   bool pos = false;
   for (i=0; i < table.size()/2; i++) {
     double par = table[i*2];
     double val = table[i*2+1];
     if( _convMode==0 )
     {
-      CASCatch_TRY
-      {
+      try {
+#ifdef NO_CAS_CATCH
+        OCC_CATCH_SIGNALS;
+#endif
        val = pow( 10.0, val );
-      }
-      CASCatch_CATCH(CASCatch_Failure)
-      {
-       aCatchSignals.Deactivate();
-       Handle(CASCatch_Failure) aFail = CASCatch_Failure::Caught();
+      } catch(Standard_Failure) {
+       Handle(Standard_Failure) aFail = Standard_Failure::Caught();
        throw SALOME_Exception( LOCALIZED( "invalid value"));
        return;
       }
@@ -251,7 +260,6 @@ void StdMeshers_NumberOfSegments::SetTableFunction(const std::vector<double>& ta
     }
     prev = par;
   }
-  aCatchSignals.Deactivate();
 
   if( !pos )
     throw SALOME_Exception(LOCALIZED("value of table function is not positive"));
@@ -269,7 +277,7 @@ void StdMeshers_NumberOfSegments::SetTableFunction(const std::vector<double>& ta
  */
 //================================================================================
 
-const std::vector<double>& StdMeshers_NumberOfSegments::GetTableFunction() const
+const vector<double>& StdMeshers_NumberOfSegments::GetTableFunction() const
   throw(SALOME_Exception)
 {
   if (_distrType != DT_TabFunc)
@@ -313,24 +321,18 @@ bool process( const TCollection_AsciiString& str, int convMode,
              bool& non_neg, bool& non_zero,
              bool& singulars, double& sing_point )
 {
-  OSD::SetSignal( true );
-  CASCatch_CatchSignals aCatchSignals;
-  aCatchSignals.Activate();
-
   bool parsed_ok = true;
   Handle( ExprIntrp_GenExp ) myExpr;
-  CASCatch_TRY
-  {
+  try {
+#ifdef NO_CAS_CATCH
+    OCC_CATCH_SIGNALS;
+#endif
     myExpr = ExprIntrp_GenExp::Create();
     myExpr->Process( str.ToCString() );
-  }
-  CASCatch_CATCH(CASCatch_Failure)
-  {
-    aCatchSignals.Deactivate();
-    Handle(CASCatch_Failure) aFail = CASCatch_Failure::Caught();
+  } catch(Standard_Failure) {
+    Handle(Standard_Failure) aFail = Standard_Failure::Caught();
     parsed_ok = false;
   }
-  aCatchSignals.Deactivate();
 
   syntax = false;
   args = false;
@@ -383,7 +385,8 @@ void StdMeshers_NumberOfSegments::SetExpressionFunction(const char* expr)
   throw(SALOME_Exception)
 {
   if (_distrType != DT_ExprFunc)
-    throw SALOME_Exception(LOCALIZED("not an expression function distribution"));
+    _distrType = DT_ExprFunc;
+    //throw SALOME_Exception(LOCALIZED("not an expression function distribution"));
 
   // remove white spaces
   TCollection_AsciiString str((Standard_CString)expr);
@@ -415,7 +418,7 @@ void StdMeshers_NumberOfSegments::SetExpressionFunction(const char* expr)
     return;
   }
   
-  std::string func = expr;
+  string func = expr;
   if( _func != func )
   {
     _func = func;
@@ -446,8 +449,8 @@ const char* StdMeshers_NumberOfSegments::GetExpressionFunction() const
 void StdMeshers_NumberOfSegments::SetConversionMode( int conv )
   throw(SALOME_Exception)
 {
-  if (_distrType != DT_TabFunc && _distrType != DT_ExprFunc)
-    throw SALOME_Exception(LOCALIZED("not a functional distribution"));
+//   if (_distrType != DT_TabFunc && _distrType != DT_ExprFunc)
+//     throw SALOME_Exception(LOCALIZED("not a functional distribution"));
 
   if( conv != _convMode )
   {
@@ -465,8 +468,8 @@ void StdMeshers_NumberOfSegments::SetConversionMode( int conv )
 int StdMeshers_NumberOfSegments::ConversionMode() const
   throw(SALOME_Exception)
 {
-  if (_distrType != DT_TabFunc && _distrType != DT_ExprFunc)
-    throw SALOME_Exception(LOCALIZED("not a functional distribution"));
+//   if (_distrType != DT_TabFunc && _distrType != DT_ExprFunc)
+//     throw SALOME_Exception(LOCALIZED("not a functional distribution"));
   return _convMode;
 }
 
@@ -478,6 +481,7 @@ int StdMeshers_NumberOfSegments::ConversionMode() const
 
 ostream & StdMeshers_NumberOfSegments::SaveTo(ostream & save)
 {
+  int listSize = _edgeIDs.size();
   save << _numberOfSegments << " " << (int)_distrType;
   switch (_distrType)
   {
@@ -500,6 +504,13 @@ ostream & StdMeshers_NumberOfSegments::SaveTo(ostream & save)
 
   if (_distrType == DT_TabFunc || _distrType == DT_ExprFunc)
     save << " " << _convMode;
+
+  if ( _distrType != DT_Regular && listSize > 0 ) {
+    save << " " << listSize;
+    for ( int i = 0; i < listSize; i++ )
+      save << " " << _edgeIDs[i];
+    save << " " << _objEntry;
+  }
   
   return save;
 }
@@ -616,6 +627,18 @@ istream & StdMeshers_NumberOfSegments::LoadFrom(istream & load)
       load.clear(ios::badbit | load.rdstate());
   }
 
+  // load reversed edges IDs
+  int intVal;
+  isOK = (load >> intVal);
+  if ( isOK && _distrType != DT_Regular && intVal > 0 ) {
+    _edgeIDs.reserve( intVal );
+    for (int i = 0; i < _edgeIDs.capacity() && isOK; i++) {
+      isOK = (load >> intVal);
+      if ( isOK ) _edgeIDs.push_back( intVal );
+    }
+    isOK = (load >> _objEntry);
+  }
+
   return load;
 }
 
@@ -640,3 +663,70 @@ istream & operator >>(istream & load, StdMeshers_NumberOfSegments & hyp)
 {
   return hyp.LoadFrom( load );
 }
+
+//================================================================================
+/*!
+ * \brief Initialize number of segments by the mesh built on the geometry
+ * \param theMesh - the built mesh
+ * \param theShape - the geometry of interest
+ * \retval bool - true if parameter values have been successfully defined
+ */
+//================================================================================
+
+bool StdMeshers_NumberOfSegments::SetParametersByMesh(const SMESH_Mesh*   theMesh,
+                                                      const TopoDS_Shape& theShape)
+{
+  if ( !theMesh || theShape.IsNull() )
+    return false;
+
+  _numberOfSegments = 0;
+  _distrType = DT_Regular;
+
+  int nbEdges = 0;
+  TopTools_IndexedMapOfShape edgeMap;
+  TopExp::MapShapes( theShape, TopAbs_EDGE, edgeMap );
+  SMESHDS_Mesh* aMeshDS = const_cast< SMESH_Mesh* >( theMesh )->GetMeshDS();
+  for ( int i = 1; i <= edgeMap.Extent(); ++i )
+  {
+    // get current segment length
+    SMESHDS_SubMesh * eSubMesh = aMeshDS->MeshElements( edgeMap( i ));
+    if ( eSubMesh && eSubMesh->NbElements())
+      _numberOfSegments += eSubMesh->NbElements();
+
+    ++nbEdges;
+  }
+  if ( nbEdges )
+    _numberOfSegments /= nbEdges;
+
+  if (_numberOfSegments == 0) _numberOfSegments = 1;
+
+  return nbEdges;
+}
+//================================================================================
+/*!
+ * \brief Initialize my parameter values by default parameters.
+ *  \retval bool - true if parameter values have been successfully defined
+ */
+//================================================================================
+
+bool StdMeshers_NumberOfSegments::SetParametersByDefaults(const TDefaults&  dflts,
+                                                          const SMESH_Mesh* /*theMesh*/)
+{
+  return (_numberOfSegments = dflts._nbSegments );
+}
+
+//=============================================================================
+/*!
+ *  
+ */
+//=============================================================================
+
+void StdMeshers_NumberOfSegments::SetReversedEdges( std::vector<int>& ids )
+{
+  if ( ids != _edgeIDs ) {
+    _edgeIDs = ids;
+
+    NotifySubMeshesHypothesisModification();
+  }
+}
+