Salome HOME
PAL13615 (EDF PAL 315/31 GEOM SMESH : meshing of a "5 edges quadrangle")
authoreap <eap@opencascade.com>
Thu, 22 Feb 2007 07:25:37 +0000 (07:25 +0000)
committereap <eap@opencascade.com>
Thu, 22 Feb 2007 07:25:37 +0000 (07:25 +0000)
   take into account medium nodes if not all edges have equadratic elements

src/StdMeshers/StdMeshers_CompositeSegment_1D.cxx
src/StdMeshers/StdMeshers_FaceSide.cxx
src/StdMeshers/StdMeshers_FaceSide.hxx
src/StdMeshers/StdMeshers_MEFISTO_2D.cxx
src/StdMeshers/StdMeshers_Quadrangle_2D.cxx

index 4540ef2a41f421fd634124fc81a35de087afe1da..0f0a2411be071f1f1c6051aebf2aeea3f01c5393 100644 (file)
@@ -292,7 +292,7 @@ StdMeshers_CompositeSegment_1D::GetFaceSide(SMESH_Mesh&        aMesh,
       eNext = nextC1Edge( eNext, aMesh, forward );
     }
   }
       eNext = nextC1Edge( eNext, aMesh, forward );
     }
   }
-  return new StdMeshers_FaceSide( aFace, edges, &aMesh, true );
+  return new StdMeshers_FaceSide( aFace, edges, &aMesh, true, false );
 }
 
 //=============================================================================
 }
 
 //=============================================================================
index b83bfb4b0685f16427c415c3d40d1942ec6a58e7..081d2a79f14391059b06f933b1514ad037b43644 100644 (file)
 StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face& theFace,
                                          const TopoDS_Edge& theEdge,
                                          SMESH_Mesh*        theMesh,
 StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face& theFace,
                                          const TopoDS_Edge& theEdge,
                                          SMESH_Mesh*        theMesh,
-                                         const bool         theIsForward)
+                                         const bool         theIsForward,
+                                         const bool         theIgnoreMediumNodes)
 {
   list<TopoDS_Edge> edges(1,theEdge);
 {
   list<TopoDS_Edge> edges(1,theEdge);
-  *this = StdMeshers_FaceSide( theFace, edges, theMesh, theIsForward );
+  *this = StdMeshers_FaceSide( theFace, edges, theMesh, theIsForward, theIgnoreMediumNodes );
 }
 
 //================================================================================
 }
 
 //================================================================================
@@ -78,7 +79,8 @@ StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face& theFace,
 StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face& theFace,
                                          list<TopoDS_Edge>& theEdges,
                                          SMESH_Mesh*        theMesh,
 StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face& theFace,
                                          list<TopoDS_Edge>& theEdges,
                                          SMESH_Mesh*        theMesh,
-                                         const bool         theIsForward)
+                                         const bool         theIsForward,
+                                         const bool         theIgnoreMediumNodes)
 {
   int nbEdges = theEdges.size();
   myEdge.resize( nbEdges );
 {
   int nbEdges = theEdges.size();
   myEdge.resize( nbEdges );
@@ -90,6 +92,7 @@ StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face& theFace,
   myNbPonits = myNbSegments = 0;
   myMesh = theMesh;
   myMissingVertexNodes = false;
   myNbPonits = myNbSegments = 0;
   myMesh = theMesh;
   myMissingVertexNodes = false;
+  myIgnoreMediumNodes = theIgnoreMediumNodes;
   if ( nbEdges == 0 ) return;
 
   SMESHDS_Mesh* meshDS = theMesh->GetMeshDS();
   if ( nbEdges == 0 ) return;
 
   SMESHDS_Mesh* meshDS = theMesh->GetMeshDS();
@@ -115,9 +118,11 @@ StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face& theFace,
 
     if ( SMESHDS_SubMesh* sm = meshDS->MeshElements( *edge )) {
       int nbN = sm->NbNodes();
 
     if ( SMESHDS_SubMesh* sm = meshDS->MeshElements( *edge )) {
       int nbN = sm->NbNodes();
-      SMDS_ElemIteratorPtr elemIt = sm->GetElements();
-      if ( elemIt->more() && elemIt->next()->IsQuadratic() )
-        nbN -= sm->NbElements();
+      if ( theIgnoreMediumNodes ) {
+        SMDS_ElemIteratorPtr elemIt = sm->GetElements();
+        if ( elemIt->more() && elemIt->next()->IsQuadratic() )
+          nbN -= sm->NbElements();
+      }
       myNbPonits += nbN;
       myNbSegments += sm->NbElements();
     }
       myNbPonits += nbN;
       myNbSegments += sm->NbElements();
     }
@@ -199,7 +204,7 @@ const vector<UVPtStruct>& StdMeshers_FaceSide::GetUVPtStruct(bool   isXConst,
       double paramSize = myLast[i] - myFirst[i], r = myNormPar[i] - prevNormPar;
       while ( nItr->more() ) {
         const SMDS_MeshNode* node = nItr->next();
       double paramSize = myLast[i] - myFirst[i], r = myNormPar[i] - prevNormPar;
       while ( nItr->more() ) {
         const SMDS_MeshNode* node = nItr->next();
-        if ( SMESH_MeshEditor::IsMedium( node, SMDSAbs_Edge ))
+        if ( myIgnoreMediumNodes && SMESH_MeshEditor::IsMedium( node, SMDSAbs_Edge ))
           continue;
         const SMDS_EdgePosition* epos =
           static_cast<const SMDS_EdgePosition*>(node->GetPosition().get());
           continue;
         const SMDS_EdgePosition* epos =
           static_cast<const SMDS_EdgePosition*>(node->GetPosition().get());
index 4320bccdc1cbf7fd961ba6ae730f9bbdcd4e88e1..c1fdd595838253ebbca84d0236f21886a51904fb 100644 (file)
@@ -78,14 +78,16 @@ public:
   StdMeshers_FaceSide(const TopoDS_Face& theFace,
                       const TopoDS_Edge& theEdge,
                       SMESH_Mesh*        theMesh,
   StdMeshers_FaceSide(const TopoDS_Face& theFace,
                       const TopoDS_Edge& theEdge,
                       SMESH_Mesh*        theMesh,
-                      const bool         theIsForward);
+                      const bool         theIsForward,
+                      const bool         theIgnoreMediumNodes);
   /*!
    * \brief Wrap several edges. Edges must be properly ordered and oriented.
    */
   StdMeshers_FaceSide(const TopoDS_Face& theFace,
                       list<TopoDS_Edge>& theEdges,
                       SMESH_Mesh*        theMesh,
   /*!
    * \brief Wrap several edges. Edges must be properly ordered and oriented.
    */
   StdMeshers_FaceSide(const TopoDS_Face& theFace,
                       list<TopoDS_Edge>& theEdges,
                       SMESH_Mesh*        theMesh,
-                      const bool         theIsForward);
+                      const bool         theIsForward,
+                      const bool         theIgnoreMediumNodes);
   /*!
    * \brief Change orientation of side geometry
    */
   /*!
    * \brief Change orientation of side geometry
    */
@@ -180,13 +182,12 @@ protected:
   vector<uvPtStruct>           myPoints, myFalsePoints;
   vector<TopoDS_Edge>          myEdge;
   vector<Handle(Geom2d_Curve)> myC2d;
   vector<uvPtStruct>           myPoints, myFalsePoints;
   vector<TopoDS_Edge>          myEdge;
   vector<Handle(Geom2d_Curve)> myC2d;
-  vector<double>               myFirst;
-  vector<double>               myLast;
+  vector<double>               myFirst, myLast;
   vector<double>               myNormPar;
   double                       myLength;
   int                          myNbPonits, myNbSegments;
   SMESH_Mesh*                  myMesh;
   vector<double>               myNormPar;
   double                       myLength;
   int                          myNbPonits, myNbSegments;
   SMESH_Mesh*                  myMesh;
-  bool                         myMissingVertexNodes;
+  bool                         myMissingVertexNodes, myIgnoreMediumNodes;
 };
 
 
 };
 
 
index 942f5fa74d0a476cfce382aa33edc8be9009e0f2..db356207ea6411c033efddb3443f73a32c41c1a1 100644 (file)
@@ -178,6 +178,7 @@ bool StdMeshers_MEFISTO_2D::Compute(SMESH_Mesh & aMesh, const TopoDS_Shape & aSh
   MESSAGE("StdMeshers_MEFISTO_2D::Compute");
 
   TopoDS_Face F = TopoDS::Face(aShape.Oriented(TopAbs_FORWARD));
   MESSAGE("StdMeshers_MEFISTO_2D::Compute");
 
   TopoDS_Face F = TopoDS::Face(aShape.Oriented(TopAbs_FORWARD));
+  const bool ignoreMediumNodes = _quadraticMesh;
 
   // get all edges of a face
   TopoDS_Vertex V1;
 
   // get all edges of a face
   TopoDS_Vertex V1;
@@ -208,7 +209,8 @@ bool StdMeshers_MEFISTO_2D::Compute(SMESH_Mesh & aMesh, const TopoDS_Shape & aSh
         return false;
       }
     }
         return false;
       }
     }
-    StdMeshers_FaceSide* wire = new StdMeshers_FaceSide( F, wireEdges, &aMesh, true );
+    StdMeshers_FaceSide* wire = new StdMeshers_FaceSide( F, wireEdges, &aMesh,
+                                                         true, ignoreMediumNodes);
     wires[ iW ] = StdMeshers_FaceSidePtr( wire );
     if (_hypLengthFromEdges && wire->NbSegments() )
       _edgeLength += wire->Length() / wire->NbSegments();
     wires[ iW ] = StdMeshers_FaceSidePtr( wire );
     if (_hypLengthFromEdges && wire->NbSegments() )
       _edgeLength += wire->Length() / wire->NbSegments();
index 52a35e8f097d3fc844048ff9741fd82ea57708b3..79e56a799f674df1791513d6a574eea573150848 100644 (file)
@@ -565,6 +565,7 @@ FaceQuadStruct* StdMeshers_Quadrangle_2D::CheckNbEdges(SMESH_Mesh &         aMes
   Unexpect aCatch(SalomeException);
 
   const TopoDS_Face & F = TopoDS::Face(aShape);
   Unexpect aCatch(SalomeException);
 
   const TopoDS_Face & F = TopoDS::Face(aShape);
+  const bool ignoreMediumNodes = _quadraticMesh;
 
   // verify 1 wire only, with 4 edges
   TopoDS_Vertex V;
 
   // verify 1 wire only, with 4 edges
   TopoDS_Vertex V;
@@ -584,9 +585,10 @@ FaceQuadStruct* StdMeshers_Quadrangle_2D::CheckNbEdges(SMESH_Mesh &         aMes
   list< TopoDS_Edge >::iterator edgeIt = edges.begin();
   if ( nbEdgesInWire.front() == 4 ) { // exactly 4 edges
     for ( ; edgeIt != edges.end(); ++edgeIt, nbSides++ )
   list< TopoDS_Edge >::iterator edgeIt = edges.begin();
   if ( nbEdgesInWire.front() == 4 ) { // exactly 4 edges
     for ( ; edgeIt != edges.end(); ++edgeIt, nbSides++ )
-      quad->side[nbSides] = new StdMeshers_FaceSide(F,*edgeIt,&aMesh,nbSides<TOP_SIDE);
+      quad->side[nbSides] = new StdMeshers_FaceSide(F, *edgeIt, &aMesh,
+                                                    nbSides<TOP_SIDE, ignoreMediumNodes);
   }
   }
-  else if ( nbEdgesInWire.front() > 4 ) { // more than 4 edges - try to unite
+  else if ( nbEdgesInWire.front() > 4 ) { // more than 4 edges - try to unite some
     list< TopoDS_Edge > sideEdges;
     while ( edgeIt != edges.end()) {
       sideEdges.clear();
     list< TopoDS_Edge > sideEdges;
     while ( edgeIt != edges.end()) {
       sideEdges.clear();
@@ -598,7 +600,8 @@ FaceQuadStruct* StdMeshers_Quadrangle_2D::CheckNbEdges(SMESH_Mesh &         aMes
         if ( sameSide )
           sideEdges.push_back( *edgeIt++ );
       }
         if ( sameSide )
           sideEdges.push_back( *edgeIt++ );
       }
-      quad->side[nbSides] = new StdMeshers_FaceSide(F,sideEdges,&aMesh,nbSides<TOP_SIDE);
+      quad->side[nbSides] = new StdMeshers_FaceSide(F, sideEdges, &aMesh,
+                                                    nbSides<TOP_SIDE, ignoreMediumNodes);
       ++nbSides;
     }
   }
       ++nbSides;
     }
   }