Salome HOME
Regression of 21397: EDF SMESH: a quadrangle face mesh can't be projected to a cylinder
authoreap <eap@opencascade.com>
Thu, 9 Apr 2015 18:07:54 +0000 (21:07 +0300)
committereap <eap@opencascade.com>
Thu, 9 Apr 2015 18:07:54 +0000 (21:07 +0300)
+ 52675: Viscous Layers construction fails with certain hypotheses
    Fix normal to be well visible from all faces -> thickness reaches 0.0007

src/StdMeshers/StdMeshers_ProjectionUtils.cxx
src/StdMeshers/StdMeshers_Projection_2D.cxx
src/StdMeshers/StdMeshers_ViscousLayers.cxx

index 85ecfac3a8650408b82005747fb417b41e25b0a6..443ddc43ce37222e646b605d7d7338e362c366d7 100644 (file)
@@ -1410,7 +1410,7 @@ int StdMeshers_ProjectionUtils::FindFaceAssociation(const TopoDS_Face&    face1,
       }
       if ( TopExp::FirstVertex( *edgeIt ).IsSame( TopExp::LastVertex( *edgeIt )) &&
            SMESH_Algo::isDegenerated( *edgeIt )) {
-        --edgeIt; // skip a degenerated edge (www.salome-platform.org/forum/forum_11/173031193)
+        --edgeIt; // skip a degenerated edge (test 3D_mesh_Projection_00/A3)
       }
       if ( !VV1[1].IsSame( TopExp::FirstVertex( *edgeIt, true ))) {
         CONT_BAD_RESULT("GetOrderedEdges() failed");
index 2db02b2ce110fcd767728e8109252e6da13e4bcb..853875c5d2d300c895ce88d4cb94b6884c6f4a26 100644 (file)
@@ -56,6 +56,7 @@
 #include <TopExp_Explorer.hxx>
 #include <TopTools_DataMapIteratorOfDataMapOfShapeShape.hxx>
 #include <TopTools_ListIteratorOfListOfShape.hxx>
+#include <TopTools_MapOfShape.hxx>
 #include <TopoDS.hxx>
 #include <gp_Ax2.hxx>
 #include <gp_Ax3.hxx>
@@ -1183,6 +1184,8 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
         // and target faces in the mapper. Thus we select srcV1 so that
         // GetOrderedEdges() to return EDGEs in a needed order
         TopoDS_Face tgtFaceBis = tgtFace;
+        TopTools_MapOfShape checkedVMap( tgtEdges.size() );
+        checkedVMap.Add ( srcV1 );
         for ( vSrcExp.Next(); vSrcExp.More(); )
         {
           tgtFaceBis.Reverse();
@@ -1191,7 +1194,7 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
           bool ok = true;
           list< TopoDS_Edge >::iterator edgeS = srcEdges.begin(), edgeT = tgtEdges.begin();
           for ( ; edgeS != srcEdges.end() && ok ; ++edgeS, ++edgeT )
-            ok = edgeS->IsSame( shape2ShapeMap( *edgeT ));
+            ok = edgeT->IsSame( shape2ShapeMap( *edgeS, /*isSrc=*/true ));
           if ( ok )
             break; // FOUND!
 
@@ -1199,6 +1202,8 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
           if ( reverse )
           {
             vSrcExp.Next();
+            while ( vSrcExp.More() && !checkedVMap.Add( vSrcExp.Current() ))
+              vSrcExp.Next();
           }
           else
           {
@@ -1209,6 +1214,38 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
           }
         }
       }
+      // for the case: project to a closed face from a non-closed face w/o vertex assoc;
+      // avoid projecting to a seam from two EDGEs with different nb nodes on them
+      // ( test mesh_Projection_2D_01/B1 )
+      if ( !_sourceHypo->HasVertexAssociation() &&
+           nbEdgesInWires.front() > 2 &&
+           helper.IsRealSeam( tgtEdges.front() ))
+      {
+        TopoDS_Shape srcEdge1 = shape2ShapeMap( tgtEdges.front() );
+        list< TopoDS_Edge >::iterator srcEdge2 =
+          std::find( srcEdges.begin(), srcEdges.end(), srcEdge1);
+        list< TopoDS_Edge >::iterator srcEdge3 =
+          std::find( srcEdges.begin(), srcEdges.end(), srcEdge1.Reversed());
+        if ( srcEdge2 == srcEdges.end() || srcEdge3 == srcEdges.end() ) // srcEdge1 is not a seam
+        {
+          // find srcEdge2 which also will be projected to tgtEdges.front()
+          for ( srcEdge2 = srcEdges.begin(); srcEdge2 != srcEdges.end(); ++srcEdge2 )
+            if ( !srcEdge1.IsSame( *srcEdge2 ) &&
+                 tgtEdges.front().IsSame( shape2ShapeMap( *srcEdge2, /*isSrc=*/true )))
+              break;
+          // compare nb nodes on srcEdge1 and srcEdge2
+          if ( srcEdge2 != srcEdges.end() )
+          {
+            int nbN1 = 0, nbN2 = 0;
+            if ( SMESHDS_SubMesh* sm = srcMesh->GetMeshDS()->MeshElements( srcEdge1 ))
+              nbN1 = sm->NbNodes();
+            if ( SMESHDS_SubMesh* sm = srcMesh->GetMeshDS()->MeshElements( *srcEdge2 ))
+              nbN2 = sm->NbNodes();
+            if ( nbN1 != nbN2 )
+              srcV1 = helper.IthVertex( 1, srcEdges.front() );
+          }
+        }
+      }
     }
     else if ( nbEdgesInWires.front() == 1 )
     {
index 61df856fee365049cd71932ddf99b2e2a1fdfea0..8167f864a2b4098759020098080ea1fccfe34517 100644 (file)
@@ -3239,74 +3239,87 @@ gp_XYZ _ViscousBuilder::getWeigthedNormal( const SMDS_MeshNode*             n,
   }
 
   // exclude equal normals
-  //int nbUniqNorms = nbFaces;
-  for ( int i = 0; i < nbFaces; ++i )
+  int nbUniqNorms = nbFaces;
+  for ( int i = 0; i < nbFaces; ++i ) {
     for ( int j = i+1; j < nbFaces; ++j )
       if ( fId2Normal[i].second.IsEqual( fId2Normal[j].second, 0.1 ))
       {
         fId2Normal[i].second.SetCoord( 0,0,0 );
-        //--nbUniqNorms;
+        --nbUniqNorms;
         break;
       }
-  //if ( nbUniqNorms < 3 )
-  {
-    for ( int i = 0; i < nbFaces; ++i )
-      resNorm += fId2Normal[i].second;
-    return resNorm;
   }
-
-  double angles[30];
   for ( int i = 0; i < nbFaces; ++i )
-  {
-    const TopoDS_Face& F = fId2Normal[i].first;
+    resNorm += fId2Normal[i].second;
 
-    // look for two EDGEs shared by F and other FACEs within fId2Normal
-    TopoDS_Edge ee[2];
-    int nbE = 0;
-    PShapeIteratorPtr eIt = SMESH_MesherHelper::GetAncestors( V, *_mesh, TopAbs_EDGE );
-    while ( const TopoDS_Shape* E = eIt->next() )
-    {
-      if ( !SMESH_MesherHelper::IsSubShape( *E, F ))
-        continue;
-      bool isSharedEdge = false;
-      for ( int j = 0; j < nbFaces && !isSharedEdge; ++j )
-      {
-        if ( i == j ) continue;
-        const TopoDS_Shape& otherF = fId2Normal[j].first;
-        isSharedEdge = SMESH_MesherHelper::IsSubShape( *E, otherF );
-      }
-      if ( !isSharedEdge )
-        continue;
-      ee[ nbE ] = TopoDS::Edge( *E );
-      ee[ nbE ].Orientation( SMESH_MesherHelper::GetSubShapeOri( F, *E ));
-      if ( ++nbE == 2 )
-        break;
-    }
-
-    // get an angle between the two EDGEs
-    angles[i] = 0;
-    if ( nbE < 1 ) continue;
-    if ( nbE == 1 )
-    {
-      ee[ 1 ] == ee[ 0 ];
-    }
-    else
+  // assure that resNorm is visible by every FACE (IPAL0052675)
+  if ( nbUniqNorms > 3 )
+  {
+    bool change = false;
+    for ( int nbAttempts = 0; nbAttempts < nbFaces; ++nbAttempts)
     {
-      if ( !V.IsSame( SMESH_MesherHelper::IthVertex( 0, ee[ 1 ] )))
-        std::swap( ee[0], ee[1] );
+      for ( int i = 0; i < nbFaces; ++i )
+        if ( resNorm * fId2Normal[i].second < 0.5 )
+        {
+          resNorm += fId2Normal[i].second;
+          change = true;
+        }
+      if ( !change ) break;
     }
-    angles[i] = SMESH_MesherHelper::GetAngle( ee[0], ee[1], F, TopoDS::Vertex( V ));
   }
 
-  // compute a weighted normal
-  double sumAngle = 0;
-  for ( int i = 0; i < nbFaces; ++i )
-  {
-    angles[i] = ( angles[i] > 2*M_PI )  ?  0  :  M_PI - angles[i];
-    sumAngle += angles[i];
-  }
-  for ( int i = 0; i < nbFaces; ++i )
-    resNorm += angles[i] / sumAngle * fId2Normal[i].second;
+  // double angles[30];
+  // for ( int i = 0; i < nbFaces; ++i )
+  // {
+  //   const TopoDS_Face& F = fId2Normal[i].first;
+
+  //   // look for two EDGEs shared by F and other FACEs within fId2Normal
+  //   TopoDS_Edge ee[2];
+  //   int nbE = 0;
+  //   PShapeIteratorPtr eIt = SMESH_MesherHelper::GetAncestors( V, *_mesh, TopAbs_EDGE );
+  //   while ( const TopoDS_Shape* E = eIt->next() )
+  //   {
+  //     if ( !SMESH_MesherHelper::IsSubShape( *E, F ))
+  //       continue;
+  //     bool isSharedEdge = false;
+  //     for ( int j = 0; j < nbFaces && !isSharedEdge; ++j )
+  //     {
+  //       if ( i == j ) continue;
+  //       const TopoDS_Shape& otherF = fId2Normal[j].first;
+  //       isSharedEdge = SMESH_MesherHelper::IsSubShape( *E, otherF );
+  //     }
+  //     if ( !isSharedEdge )
+  //       continue;
+  //     ee[ nbE ] = TopoDS::Edge( *E );
+  //     ee[ nbE ].Orientation( SMESH_MesherHelper::GetSubShapeOri( F, *E ));
+  //     if ( ++nbE == 2 )
+  //       break;
+  //   }
+
+  //   // get an angle between the two EDGEs
+  //   angles[i] = 0;
+  //   if ( nbE < 1 ) continue;
+  //   if ( nbE == 1 )
+  //   {
+  //     ee[ 1 ] == ee[ 0 ];
+  //   }
+  //   else
+  //   {
+  //     if ( !V.IsSame( SMESH_MesherHelper::IthVertex( 0, ee[ 1 ] )))
+  //       std::swap( ee[0], ee[1] );
+  //   }
+  //   angles[i] = SMESH_MesherHelper::GetAngle( ee[0], ee[1], F, TopoDS::Vertex( V ));
+  // }
+
+  // // compute a weighted normal
+  // double sumAngle = 0;
+  // for ( int i = 0; i < nbFaces; ++i )
+  // {
+  //   angles[i] = ( angles[i] > 2*M_PI )  ?  0  :  M_PI - angles[i];
+  //   sumAngle += angles[i];
+  // }
+  // for ( int i = 0; i < nbFaces; ++i )
+  //   resNorm += angles[i] / sumAngle * fId2Normal[i].second;
 
   return resNorm;
 }