Additional fix for issue 16186.

Quadrangle split into 2 triangles by smallest diagonal.
This commit is contained in:
rnv 2008-09-17 13:47:00 +00:00
parent 6eb4c26173
commit 06e979fa3e
2 changed files with 44 additions and 17 deletions

View File

@ -329,10 +329,7 @@ bool StdMeshers_Quadrangle_2D::Compute (SMESH_Mesh& aMesh,
meshDS->SetMeshElementOnShape(face, geomFaceID); meshDS->SetMeshElementOnShape(face, geomFaceID);
} }
else { else {
SMDS_MeshFace* face = myTool->AddFace(a, b, c); SplitQuad(meshDS, geomFaceID, a, b, c, d);
meshDS->SetMeshElementOnShape(face, geomFaceID);
face = myTool->AddFace(a, c, d);
meshDS->SetMeshElementOnShape(face, geomFaceID);
} }
// if node d is not at position g - make additional triangles // 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); meshDS->SetMeshElementOnShape(face, geomFaceID);
} }
else { else {
SMDS_MeshFace* face = myTool->AddFace(a, b, c); SplitQuad(meshDS, geomFaceID, a, b, c, d);
meshDS->SetMeshElementOnShape(face, geomFaceID);
face = myTool->AddFace(a, c, d);
meshDS->SetMeshElementOnShape(face, geomFaceID);
} }
if (near + 1 < g) { // if d not is at g - make additional triangles 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); meshDS->SetMeshElementOnShape(face, geomFaceID);
} }
else { else {
SMDS_MeshFace* face = myTool->AddFace(a, b, c); SplitQuad(meshDS, geomFaceID, a, b, c, d);
meshDS->SetMeshElementOnShape(face, geomFaceID);
face = myTool->AddFace(a, c, d);
meshDS->SetMeshElementOnShape(face, geomFaceID);
} }
if (near - 1 > g) { // if d not is at g - make additional triangles 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); meshDS->SetMeshElementOnShape(face, geomFaceID);
} }
else { else {
SMDS_MeshFace* face = myTool->AddFace(a, b, c); SplitQuad(meshDS, geomFaceID, a, b, c, d);
meshDS->SetMeshElementOnShape(face, geomFaceID);
face = myTool->AddFace(a, c, d);
meshDS->SetMeshElementOnShape(face, geomFaceID);
} }
if (near + 1 < g) { // if d not is at g - make additional triangles 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; 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 );
}
}

View File

@ -38,9 +38,9 @@
class SMESH_Mesh; class SMESH_Mesh;
class SMESH_MesherHelper; class SMESH_MesherHelper;
class StdMeshers_FaceSide; class StdMeshers_FaceSide;
class SMDS_MeshNode;
struct uvPtStruct; struct uvPtStruct;
//class SMDS_MeshNode;
enum TSideID { BOTTOM_SIDE=0, RIGHT_SIDE, TOP_SIDE, LEFT_SIDE, NB_SIDES }; enum TSideID { BOTTOM_SIDE=0, RIGHT_SIDE, TOP_SIDE, LEFT_SIDE, NB_SIDES };
@ -79,10 +79,18 @@ protected:
const TopoDS_Shape& aShape, const TopoDS_Shape& aShape,
FaceQuadStruct*& quad); 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 * Special function for creation only quandrangle faces
*/ */
bool ComputeQuadPref(SMESH_Mesh& aMesh, bool ComputeQuadPref(SMESH_Mesh& aMesh,
const TopoDS_Shape& aShape, const TopoDS_Shape& aShape,
FaceQuadStruct* quad); FaceQuadStruct* quad);