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

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

index 7984444c920c1472276eb1319d70098c65f9b6d1..19755db871f9686fbdc6a672b2fd246f8e9680bc 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 2b426188c103cbfab6a42761dd252f90e0957e18..0877ca3a68f4bbdde23ea8b2cb2d9cf1121b57a4 100644 (file)
@@ -38,9 +38,9 @@
 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 };
 
@@ -78,11 +78,19 @@ 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
    */
   bool ComputeQuadPref(SMESH_Mesh& aMesh,
+                       
                        const TopoDS_Shape& aShape,
                        FaceQuadStruct* quad);