]> SALOME platform Git repositories - modules/smesh.git/blobdiff - src/SMESH/SMESH_Algo.cxx
Salome HOME
PR: synchro V7_main tag mergefrom_V6_main_28Feb13
[modules/smesh.git] / src / SMESH / SMESH_Algo.cxx
index f1d3fad659d984b9f8912471f98850eb5455503d..9f3ed80532907460f7c3078a881b8ecf339a59bf 100644 (file)
@@ -1,4 +1,4 @@
-// Copyright (C) 2007-2011  CEA/DEN, EDF R&D, OPEN CASCADE
+// Copyright (C) 2007-2012  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
@@ -18,6 +18,7 @@
 // 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   : SMESH_Algo.cxx
@@ -38,6 +39,7 @@
 #include "SMESH_HypoFilter.hxx"
 #include "SMESH_Mesh.hxx"
 #include "SMESH_TypeDefs.hxx"
+#include "SMESH_subMesh.hxx"
 
 #include <Basics_OCCTVersion.hxx>
 
@@ -47,7 +49,9 @@
 #include <GCPnts_AbscissaPoint.hxx>
 #include <GeomAdaptor_Curve.hxx>
 #include <Geom_Surface.hxx>
+#include <LDOMParser.hxx>
 #include <TopExp.hxx>
+#include <TopExp_Explorer.hxx>
 #include <TopLoc_Location.hxx>
 #include <TopTools_ListIteratorOfListOfShape.hxx>
 #include <TopTools_ListOfShape.hxx>
@@ -55,6 +59,7 @@
 #include <TopoDS_Edge.hxx>
 #include <TopoDS_Face.hxx>
 #include <TopoDS_Vertex.hxx>
+#include <TopoDS_Wire.hxx>
 #include <gp_Pnt.hxx>
 #include <gp_Pnt2d.hxx>
 #include <gp_Vec.hxx>
 
 using namespace std;
 
+//================================================================================
+/*!
+ * \brief Returns \a true if two algorithms (described by \a this and the given
+ *        algo data) are compatible by their output and input types of elements.
+ */
+//================================================================================
+
+bool SMESH_Algo::Features::IsCompatible( const SMESH_Algo::Features& algo2 ) const
+{
+  if ( _dim > algo2._dim ) return algo2.IsCompatible( *this );
+  // algo2 is of highter dimension
+  if ( _outElemTypes.empty() || algo2._inElemTypes.empty() )
+    return false;
+  bool compatible = true;
+  set<SMDSAbs_GeometryType>::const_iterator myOutType = _outElemTypes.begin();
+  for ( ; myOutType != _outElemTypes.end() && compatible; ++myOutType )
+    compatible = algo2._inElemTypes.count( *myOutType );
+  return compatible;
+}
+
+//================================================================================
+/*!
+ * \brief Return Data of the algorithm
+ */
+//================================================================================
+
+const SMESH_Algo::Features& SMESH_Algo::GetFeatures( const std::string& algoType )
+{
+  static map< string, SMESH_Algo::Features > theFeaturesByName;
+  if ( theFeaturesByName.empty() )
+  {
+    // Read Plugin.xml files
+    vector< string > xmlPaths = SMESH_Gen::GetPluginXMLPaths();
+    LDOMParser xmlParser;
+    for ( size_t iXML = 0; iXML < xmlPaths.size(); ++iXML )
+    {
+      bool error = xmlParser.parse( xmlPaths[iXML].c_str() );
+      if ( error )
+      {
+        TCollection_AsciiString data;
+        INFOS( xmlParser.GetError(data) );
+        continue;
+      }
+      // <algorithm type="Regular_1D"
+      //            ...
+      //            input="EDGE"
+      //            output="QUAD,TRIA">
+      //
+      LDOM_Document xmlDoc = xmlParser.getDocument();
+      LDOM_NodeList algoNodeList = xmlDoc.getElementsByTagName( "algorithm" );
+      for ( int i = 0; i < algoNodeList.getLength(); ++i )
+      {
+        LDOM_Node     algoNode           = algoNodeList.item( i );
+        LDOM_Element& algoElem           = (LDOM_Element&) algoNode;
+        TCollection_AsciiString algoType = algoElem.getAttribute("type");
+        TCollection_AsciiString input    = algoElem.getAttribute("input");
+        TCollection_AsciiString output   = algoElem.getAttribute("output");
+        TCollection_AsciiString dim      = algoElem.getAttribute("dim");
+        TCollection_AsciiString label    = algoElem.getAttribute("label-id");
+        if ( algoType.IsEmpty() ) continue;
+
+        Features & data = theFeaturesByName[ algoType.ToCString() ];
+        data._dim   = dim.IntegerValue();
+        data._label = label.ToCString();
+        for ( int isInput = 0; isInput < 2; ++isInput )
+        {
+          TCollection_AsciiString&   typeStr = isInput ? input : output;
+          set<SMDSAbs_GeometryType>& typeSet = isInput ? data._inElemTypes : data._outElemTypes;
+          int beg = 1, end;
+          while ( beg <= typeStr.Length() )
+          {
+            while ( beg < typeStr.Length() && !isalpha( typeStr.Value( beg ) ))
+              ++beg;
+            end = beg;
+            while ( end < typeStr.Length() && isalpha( typeStr.Value( end + 1 ) ))
+              ++end;
+            if ( end > beg )
+            {
+              TCollection_AsciiString typeName = typeStr.SubString( beg, end );
+              if      ( typeName == "EDGE" ) typeSet.insert( SMDSGeom_EDGE );
+              else if ( typeName == "TRIA" ) typeSet.insert( SMDSGeom_TRIANGLE );
+              else if ( typeName == "QUAD" ) typeSet.insert( SMDSGeom_QUADRANGLE );
+            }
+            beg = end + 1;
+          }
+        }
+      }
+    }
+  }
+  return theFeaturesByName[ algoType ];
+}
+
 //=============================================================================
 /*!
  *  
@@ -80,9 +177,11 @@ SMESH_Algo::SMESH_Algo (int hypId, int studyId, SMESH_Gen * gen)
 {
   gen->_mapAlgo[hypId] = this;
 
-  _onlyUnaryInput = _requireDescretBoundary = _requireShape = true;
+  _onlyUnaryInput = _requireDiscreteBoundary = _requireShape = true;
   _quadraticMesh = _supportSubmeshes = false;
   _error = COMPERR_OK;
+  for ( int i = 0; i < 4; ++i )
+    _neededLowerHyps[ i ] = false;
 }
 
 //=============================================================================
@@ -95,6 +194,41 @@ SMESH_Algo::~SMESH_Algo()
 {
 }
 
+//=============================================================================
+/*!
+ *  
+ */
+//=============================================================================
+
+SMESH_0D_Algo::SMESH_0D_Algo(int hypId, int studyId, SMESH_Gen* gen)
+  : SMESH_Algo(hypId, studyId, gen)
+{
+  _shapeType = (1 << TopAbs_VERTEX);
+  _type = ALGO_0D;
+  gen->_map0D_Algo[hypId] = this;
+}
+SMESH_1D_Algo::SMESH_1D_Algo(int hypId, int studyId, SMESH_Gen* gen)
+  : SMESH_Algo(hypId, studyId, gen)
+{
+  _shapeType = (1 << TopAbs_EDGE);
+  _type = ALGO_1D;
+  gen->_map1D_Algo[hypId] = this;
+}
+SMESH_2D_Algo::SMESH_2D_Algo(int hypId, int studyId, SMESH_Gen* gen)
+  : SMESH_Algo(hypId, studyId, gen)
+{
+  _shapeType = (1 << TopAbs_FACE);
+  _type = ALGO_2D;
+  gen->_map2D_Algo[hypId] = this;
+}
+SMESH_3D_Algo::SMESH_3D_Algo(int hypId, int studyId, SMESH_Gen* gen)
+  : SMESH_Algo(hypId, studyId, gen)
+{
+  _shapeType = (1 << TopAbs_SOLID);
+  _type = ALGO_3D;
+  gen->_map3D_Algo[hypId] = this;
+}
+
 //=============================================================================
 /*!
  * Usually an algoritm has nothing to save
@@ -211,113 +345,13 @@ bool SMESH_Algo::FaceNormal(const SMDS_MeshElement* F, gp_XYZ& normal, bool norm
   return ok;
 }
 
-//================================================================================
-/*!
- * \brief Find out elements orientation on a geometrical face
- * \param theFace - The face correctly oriented in the shape being meshed
- * \param theMeshDS - The mesh data structure
- * \retval bool - true if the face normal and the normal of first element
- *                in the correspoding submesh point in different directions
+/*
+ * Moved to SMESH_MesherHelper
  */
-//================================================================================
-
-bool SMESH_Algo::IsReversedSubMesh (const TopoDS_Face&  theFace,
-                                    SMESHDS_Mesh*       theMeshDS)
-{
-  if ( theFace.IsNull() || !theMeshDS )
-    return false;
-
-  // find out orientation of a meshed face
-  int faceID = theMeshDS->ShapeToIndex( theFace );
-  TopoDS_Shape aMeshedFace = theMeshDS->IndexToShape( faceID );
-  bool isReversed = ( theFace.Orientation() != aMeshedFace.Orientation() );
-
-  const SMESHDS_SubMesh * aSubMeshDSFace = theMeshDS->MeshElements( faceID );
-  if ( !aSubMeshDSFace )
-    return isReversed;
-
-  // find element with node located on face and get its normal
-  const SMDS_FacePosition* facePos = 0;
-  int vertexID = 0;
-  gp_Pnt nPnt[3];
-  gp_Vec Ne;
-  bool normalOK = false;
-  SMDS_ElemIteratorPtr iteratorElem = aSubMeshDSFace->GetElements();
-  while ( iteratorElem->more() ) // loop on elements on theFace
-  {
-    const SMDS_MeshElement* elem = iteratorElem->next();
-    if ( elem && elem->NbNodes() > 2 ) {
-      SMDS_ElemIteratorPtr nodesIt = elem->nodesIterator();
-      const SMDS_FacePosition* fPos = 0;
-      int i = 0, vID = 0;
-      while ( nodesIt->more() ) { // loop on nodes
-        const SMDS_MeshNode* node
-          = static_cast<const SMDS_MeshNode *>(nodesIt->next());
-        if ( i == 3 ) i = 2;
-        nPnt[ i++ ].SetCoord( node->X(), node->Y(), node->Z() );
-        // check position
-        const SMDS_PositionPtr& pos = node->GetPosition();
-        if ( !pos ) continue;
-        if ( pos->GetTypeOfPosition() == SMDS_TOP_FACE ) {
-          fPos = dynamic_cast< const SMDS_FacePosition* >( pos );
-        }
-        else if ( pos->GetTypeOfPosition() == SMDS_TOP_VERTEX ) {
-          vID = node->getshapeId();
-        }
-      }
-      if ( fPos || ( !normalOK && vID )) {
-        // compute normal
-        gp_Vec v01( nPnt[0], nPnt[1] ), v02( nPnt[0], nPnt[2] );
-        if ( v01.SquareMagnitude() > RealSmall() &&
-             v02.SquareMagnitude() > RealSmall() )
-        {
-          Ne = v01 ^ v02;
-          normalOK = ( Ne.SquareMagnitude() > RealSmall() );
-        }
-        // we need position on theFace or at least on vertex
-        if ( normalOK ) {
-          vertexID = vID;
-          if ((facePos = fPos))
-            break;
-        }
-      }
-    }
-  }
-  if ( !normalOK )
-    return isReversed;
-
-  // node position on face
-  double u,v;
-  if ( facePos ) {
-    u = facePos->GetUParameter();
-    v = facePos->GetVParameter();
-  }
-  else if ( vertexID ) {
-    TopoDS_Shape V = theMeshDS->IndexToShape( vertexID );
-    if ( V.IsNull() || V.ShapeType() != TopAbs_VERTEX )
-      return isReversed;
-    gp_Pnt2d uv = BRep_Tool::Parameters( TopoDS::Vertex( V ), theFace );
-    u = uv.X();
-    v = uv.Y();
-  }
-  else
-  {
-    return isReversed;
-  }
-
-  // face normal at node position
-  TopLoc_Location loc;
-  Handle(Geom_Surface) surf = BRep_Tool::Surface( theFace, loc );
-  if ( surf.IsNull() || surf->Continuity() < GeomAbs_C1 ) return isReversed;
-  gp_Vec d1u, d1v;
-  surf->D1( u, v, nPnt[0], d1u, d1v );
-  gp_Vec Nf = (d1u ^ d1v).Transformed( loc );
-
-  if ( theFace.Orientation() == TopAbs_REVERSED )
-    Nf.Reverse();
-
-  return Ne * Nf < 0.;
-}
+// bool SMESH_Algo::IsReversedSubMesh (const TopoDS_Face&  theFace,
+//                                     SMESHDS_Mesh*       theMeshDS)
+// {
+// }
 
 //================================================================================
 /*!
@@ -500,13 +534,19 @@ GeomAbs_Shape SMESH_Algo::Continuity(TopoDS_Edge E1,
 {
   //E1.Orientation(TopAbs_FORWARD), E2.Orientation(TopAbs_FORWARD); // avoid pb with internal edges
   if (E1.Orientation() > TopAbs_REVERSED) // INTERNAL
-      E1.Orientation( TopAbs_FORWARD );
+    E1.Orientation( TopAbs_FORWARD );
   if (E2.Orientation() > TopAbs_REVERSED) // INTERNAL
-      E2.Orientation( TopAbs_FORWARD );
-  TopoDS_Vertex V = TopExp::LastVertex (E1, true);
-  if ( !V.IsSame( TopExp::FirstVertex(E2, true )))
-    if ( !TopExp::CommonVertex( E1, E2, V ))
-      return GeomAbs_C0;
+    E2.Orientation( TopAbs_FORWARD );
+
+  TopoDS_Vertex V, VV1[2], VV2[2];
+  TopExp::Vertices( E1, VV1[0], VV1[1], true );
+  TopExp::Vertices( E2, VV2[0], VV2[1], true );
+  if      ( VV1[1].IsSame( VV2[0] ))  { V = VV1[1]; }
+  else if ( VV1[0].IsSame( VV2[1] ))  { V = VV1[0]; }
+  else if ( VV1[1].IsSame( VV2[1] ))  { V = VV1[1]; E1.Reverse(); }
+  else if ( VV1[0].IsSame( VV2[0] ))  { V = VV1[0]; E1.Reverse(); }
+  else { return GeomAbs_C0; }
+
   Standard_Real u1 = BRep_Tool::Parameter( V, E1 );
   Standard_Real u2 = BRep_Tool::Parameter( V, E2 );
   BRepAdaptor_Curve C1( E1 ), C2( E2 );
@@ -681,11 +721,17 @@ bool SMESH_Algo::Compute(SMESH_Mesh & /*aMesh*/, SMESH_MesherHelper* /*aHelper*/
   return error( COMPERR_BAD_INPUT_MESH, "Mesh built on shape expected");
 }
 
-#ifdef WITH_SMESH_CANCEL_COMPUTE
+//=======================================================================
+//function : CancelCompute
+//purpose  : Sets _computeCanceled to true. It's usage depends on
+//  *        implementation of a particular mesher.
+//=======================================================================
+
 void SMESH_Algo::CancelCompute()
 {
+  _computeCanceled = true;
+  _error = COMPERR_CANCELED;
 }
-#endif
 
 //================================================================================
 /*!
@@ -747,6 +793,8 @@ void SMESH_Algo::InitComputeError()
     if ( (*elem)->GetID() < 1 )
       delete *elem;
   _badInputElements.clear();
+
+  _computeCanceled = false;
 }
 
 //================================================================================
@@ -762,3 +810,61 @@ void SMESH_Algo::addBadInputElement(const SMDS_MeshElement* elem)
   if ( elem )
     _badInputElements.push_back( elem );
 }
+
+//=======================================================================
+//function : addBadInputElements
+//purpose  : store a bad input elements or nodes preventing computation
+//=======================================================================
+
+void SMESH_Algo::addBadInputElements(const SMESHDS_SubMesh* sm,
+                                     const bool             addNodes)
+{
+  if ( sm )
+  {
+    if ( addNodes )
+    {
+      SMDS_NodeIteratorPtr nIt = sm->GetNodes();
+      while ( nIt->more() ) addBadInputElement( nIt->next() );
+    }
+    else
+    {
+      SMDS_ElemIteratorPtr eIt = sm->GetElements();
+      while ( eIt->more() ) addBadInputElement( eIt->next() );
+    }
+  }
+}
+
+//=============================================================================
+/*!
+ *  
+ */
+//=============================================================================
+
+// int SMESH_Algo::NumberOfWires(const TopoDS_Shape& S)
+// {
+//   int i = 0;
+//   for (TopExp_Explorer exp(S,TopAbs_WIRE); exp.More(); exp.Next())
+//     i++;
+//   return i;
+// }
+
+//=============================================================================
+/*!
+ *  
+ */
+//=============================================================================
+
+int SMESH_Algo::NumberOfPoints(SMESH_Mesh& aMesh, const TopoDS_Wire& W)
+{
+  int nbPoints = 0;
+  for (TopExp_Explorer exp(W,TopAbs_EDGE); exp.More(); exp.Next()) {
+    const TopoDS_Edge& E = TopoDS::Edge(exp.Current());
+    int nb = aMesh.GetSubMesh(E)->GetSubMeshDS()->NbNodes();
+    if(_quadraticMesh)
+      nb = nb/2;
+    nbPoints += nb + 1; // internal points plus 1 vertex of 2 (last point ?)
+  }
+  return nbPoints;
+}
+
+