]> SALOME platform Git repositories - modules/smesh.git/commitdiff
Salome HOME
Additional fix for issue 16186
authorrnv <rnv@opencascade.com>
Wed, 17 Sep 2008 14:21:07 +0000 (14:21 +0000)
committerrnv <rnv@opencascade.com>
Wed, 17 Sep 2008 14:21:07 +0000 (14:21 +0000)
Now quadrangle split into 2 triangles by smallest diagonal.

src/StdMeshers/StdMeshers_Quadrangle_2D.cxx
src/StdMeshers/StdMeshers_Quadrangle_2D.hxx

index aa6b281654842906c46aeae75e3f4e7c07319dd7..8dd7c4a5814e005a187a508c19586b122dab516d 100644 (file)
@@ -329,10 +329,7 @@ bool StdMeshers_Quadrangle_2D::Compute (SMESH_Mesh& aMesh,
           meshDS->SetMeshElementOnShape(face, geomFaceID);
         }
         else {
-          SMDS_MeshFace* face = myTool->AddFace(a, b, c);
-          meshDS->SetMeshElementOnShape(face, geomFaceID);
-          face = myTool->AddFace(a, c, d);
-          meshDS->SetMeshElementOnShape(face, geomFaceID);
+          SplitQuad(meshDS, geomFaceID, a, b, c, d);
         }
 
         // if node d is not at position g - make additional triangles
@@ -421,10 +418,7 @@ bool StdMeshers_Quadrangle_2D::Compute (SMESH_Mesh& aMesh,
             meshDS->SetMeshElementOnShape(face, geomFaceID);
           }
           else {
-            SMDS_MeshFace* face = myTool->AddFace(a, b, c);
-            meshDS->SetMeshElementOnShape(face, geomFaceID);
-            face = myTool->AddFace(a, c, d);
-            meshDS->SetMeshElementOnShape(face, geomFaceID);
+            SplitQuad(meshDS, geomFaceID, a, b, c, d);
           }
 
           if (near + 1 < g) { // if d not is at g - make additional triangles
@@ -499,10 +493,7 @@ bool StdMeshers_Quadrangle_2D::Compute (SMESH_Mesh& aMesh,
           meshDS->SetMeshElementOnShape(face, geomFaceID);
         }
         else {
-          SMDS_MeshFace* face = myTool->AddFace(a, b, c);
-          meshDS->SetMeshElementOnShape(face, geomFaceID);
-          face = myTool->AddFace(a, c, d);
-          meshDS->SetMeshElementOnShape(face, geomFaceID);
+          SplitQuad(meshDS, geomFaceID, a, b, c, d);
         }
 
         if (near - 1 > g) { // if d not is at g - make additional triangles
@@ -573,10 +564,7 @@ bool StdMeshers_Quadrangle_2D::Compute (SMESH_Mesh& aMesh,
             meshDS->SetMeshElementOnShape(face, geomFaceID);
           }
           else {
-            SMDS_MeshFace* face = myTool->AddFace(a, b, c);
-            meshDS->SetMeshElementOnShape(face, geomFaceID);
-            face = myTool->AddFace(a, c, d);
-            meshDS->SetMeshElementOnShape(face, geomFaceID);
+            SplitQuad(meshDS, geomFaceID, a, b, c, d);
           }
 
           if (near + 1 < g) { // if d not is at g - make additional triangles
@@ -1524,3 +1512,34 @@ bool StdMeshers_Quadrangle_2D::ComputeQuadPref (SMESH_Mesh &        aMesh,
   return isOk;
 }
 
+//=============================================================================
+/*! Split quadrangle in to 2 triangles by smallest diagonal
+ *   
+ */
+//=============================================================================
+void StdMeshers_Quadrangle_2D::SplitQuad(SMESHDS_Mesh *theMeshDS,
+                                    int theFaceID,
+                                    const SMDS_MeshNode* theNode1,
+                                    const SMDS_MeshNode* theNode2,
+                                    const SMDS_MeshNode* theNode3,
+                                    const SMDS_MeshNode* theNode4)
+{
+  gp_Pnt a(theNode1->X(),theNode1->Y(),theNode1->Z());
+  gp_Pnt b(theNode2->X(),theNode2->Y(),theNode2->Z());
+  gp_Pnt c(theNode3->X(),theNode3->Y(),theNode3->Z());
+  gp_Pnt d(theNode4->X(),theNode4->Y(),theNode4->Z());
+  SMDS_MeshFace* face;
+  if(a.Distance(c) > b.Distance(d)){
+    face = myTool->AddFace(theNode2, theNode4, theNode1);
+    theMeshDS->SetMeshElementOnShape(face, theFaceID );
+    face = myTool->AddFace(theNode2, theNode3, theNode4);
+    theMeshDS->SetMeshElementOnShape(face, theFaceID );
+
+  }
+  else{
+    face = myTool->AddFace(theNode1, theNode2 ,theNode3);
+    theMeshDS->SetMeshElementOnShape(face, theFaceID );
+    face = myTool->AddFace(theNode1, theNode3, theNode4);
+    theMeshDS->SetMeshElementOnShape(face, theFaceID );
+  }
+}
index cf9d78461b70939f53b1c8a12c0c8fdf340f0e06..18ae6aae9024f7e10f3af501b6e089b27b3232c1 100644 (file)
 class SMESH_Mesh;
 class SMESH_MesherHelper;
 class StdMeshers_FaceSide;
+class SMDS_MeshNode;
 struct uvPtStruct;
 
-//class SMDS_MeshNode;
-
 enum TSideID { BOTTOM_SIDE=0, RIGHT_SIDE, TOP_SIDE, LEFT_SIDE, NB_SIDES };
 
 typedef uvPtStruct UVPtStruct;
@@ -78,6 +77,13 @@ protected:
   bool SetNormalizedGrid(SMESH_Mesh& aMesh,
                         const TopoDS_Shape& aShape,
                         FaceQuadStruct*& quad);
+  
+  void SplitQuad(SMESHDS_Mesh *theMeshDS,
+                 const int theFaceID,
+                 const SMDS_MeshNode* theNode1,
+                 const SMDS_MeshNode* theNode2,
+                 const SMDS_MeshNode* theNode3,
+                 const SMDS_MeshNode* theNode4);
 
   /**
    * Special function for creation only quandrangle faces