From: nge Date: Thu, 23 Jul 2009 15:12:56 +0000 (+0000) Subject: 0020330: Pb with ghs3d as a submesh X-Git-Url: http://git.salome-platform.org/gitweb/?a=commitdiff_plain;h=97342297fe40ab536e787d3d396d946fd6b7cd90;p=plugins%2Fghs3dplugin.git 0020330: Pb with ghs3d as a submesh 0020042: EDF 864 SMESH: Mesh of holes (GHS3D/BLSurf) --- diff --git a/src/GHS3DPlugin_GHS3D.cxx b/src/GHS3DPlugin_GHS3D.cxx index b972040..23d5949 100644 --- a/src/GHS3DPlugin_GHS3D.cxx +++ b/src/GHS3DPlugin_GHS3D.cxx @@ -162,9 +162,14 @@ static TopoDS_Shape findShape(const SMDS_MeshNode *aNode[], gp_XYZ aPnt(0,0,0); int j, iShape, nbNode = 4; - for ( j=0; jX(), aNode[j]->Y(), aNode[j]->Z() ); - aPnt /= nbNode; + for ( j=0; jX(), aNode[j]->Y(), aNode[j]->Z() ); + if ( aNode[j]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE ) { + aPnt = p; + break; + } + aPnt += p; + } BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion()); if (state) *state = SC.State(); @@ -624,11 +629,23 @@ static int findShapeID(SMESH_Mesh& mesh, if ( toMeshHoles ) return meshDS->ShapeToIndex( solid1 ); - // - are we at a hole boundary face? + // - Are we at a hole boundary face? if ( shells(1).IsSame( BRepTools::OuterShell( solid1 )) ) - return meshDS->ShapeToIndex( solid1 ); // - no + { // - No, but maybe a hole is bound by two shapes? Does shells(1) touches another shell? + bool touch = false; + TopExp_Explorer eExp( shells(1), TopAbs_EDGE ); + // check if any edge of shells(1) belongs to another shell + for ( ; eExp.More() && !touch; eExp.Next() ) { + ansIt = mesh.GetAncestors( eExp.Current() ); + for ( ; ansIt.More() && !touch; ansIt.Next() ) { + if ( ansIt.Value().ShapeType() == TopAbs_SHELL ) + touch = ( !ansIt.Value().IsSame( shells(1) )); + } + } + if (!touch) + return meshDS->ShapeToIndex( solid1 ); + } } - // find orientation of geom face within the first solid TopExp_Explorer fExp( solid1, TopAbs_FACE ); for ( ; fExp.More(); fExp.Next() ) @@ -639,24 +656,8 @@ static int findShapeID(SMESH_Mesh& mesh, if ( !fExp.More() ) return invalidID; // face not found - // find UV of node1 on geomFace - SMESH_MesherHelper helper( mesh ); - gp_XY uv = helper.GetNodeUV( geomFace, node1 ); - - // check that uv is correct + // normale to triangle gp_Pnt node1Pnt ( node1->X(), node1->Y(), node1->Z() ); - double tol = BRep_Tool::Tolerance( geomFace ); - BRepAdaptor_Surface surface( geomFace ); - if ( node1Pnt.Distance( surface.Value( uv.X(), uv.Y() )) > 2 * tol ) { - // project node1 onto geomFace to get right UV - GeomAPI_ProjectPointOnSurf projector( node1Pnt, surface.Surface().Surface() ); - if ( !projector.IsDone() || projector.NbPoints() < 1 ) - return invalidID; - Quantity_Parameter U,V; - projector.LowerDistanceParameters(U,V); - uv = gp_XY( U,V ); - } - // normale to face at node1 gp_Pnt node2Pnt ( node2->X(), node2->Y(), node2->Z() ); gp_Pnt node3Pnt ( node3->X(), node3->Y(), node3->Z() ); gp_Vec vec12( node1Pnt, node2Pnt ); @@ -665,6 +666,28 @@ static int findShapeID(SMESH_Mesh& mesh, if ( meshNormal.SquareMagnitude() < DBL_MIN ) return invalidID; + // find UV of node1 on geomFace + SMESH_MesherHelper helper( mesh ); helper.SetSubShape( geomFace ); + const SMDS_MeshNode* nNotOnSeamEdge = 0; + if ( helper.IsSeamShape( node1->GetPosition()->GetShapeId() )) + if ( helper.IsSeamShape( node2->GetPosition()->GetShapeId() )) + nNotOnSeamEdge = node3; + else + nNotOnSeamEdge = node2; + gp_XY uv = helper.GetNodeUV( geomFace, node1, nNotOnSeamEdge ); + // check that uv is correct + BRepAdaptor_Surface surface( geomFace ); + double tol = BRep_Tool::Tolerance( geomFace ); + if ( node1Pnt.Distance( surface.Value( uv.X(), uv.Y() )) > 2 * tol ) { + GeomAPI_ProjectPointOnSurf projector( node1Pnt, surface.Surface().Surface(), 2 * tol ); + if ( !projector.IsDone() || projector.NbPoints() < 1 || projector.LowerDistance() > 2 * tol) + return invalidID; + Quantity_Parameter U,V; + projector.LowerDistanceParameters(U,V); + if ( node1Pnt.Distance( surface.Value( U, V )) > 2 * tol ) + return invalidID; + uv.SetCoord( U,V ); + } // normale to geomFace at UV gp_Vec du, dv; surface.D1( uv.X(), uv.Y(), node1Pnt, du, dv ); @@ -790,8 +813,19 @@ static bool readResultFile(const int fileOpen, try { OCC_CATCH_SIGNALS; tabID[i] = findShapeID( theMesh, n1, n2, n3, toMeshHoles ); + // -- 0020330: Pb with ghs3d as a submesh + // check that found shape is to be meshed + if ( tabID[i] > 0 ) { + const TopoDS_Shape& foundShape = theMeshDS->IndexToShape( tabID[i] ); + bool isToBeMeshed = false; + for ( int iS = 0; !isToBeMeshed && iS < nbShape; ++iS ) + isToBeMeshed = foundShape.IsSame( tabShape[ iS ]); + if ( !isToBeMeshed ) + tabID[i] = HOLE_ID; + } + // END -- 0020330: Pb with ghs3d as a submesh #ifdef _DEBUG_ - std::cout << i+1 << " subdomain: findShapeID() returns " << tabID[i] << std::endl; + cout << i+1 << " subdomain: findShapeID() returns " << tabID[i] << endl; #endif } catch ( Standard_Failure ) { } catch (...) {} @@ -1011,7 +1045,12 @@ bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh, bool Ok(false); SMESHDS_Mesh* meshDS = theMesh.GetMeshDS(); - _nbShape = countShape( meshDS, TopAbs_SOLID ); // we count the number of shapes + // we count the number of shapes + // _nbShape = countShape( meshDS, TopAbs_SOLID ); -- 0020330: Pb with ghs3d as a submesh + _nbShape = 0; + TopExp_Explorer expBox ( theShape, TopAbs_SOLID ); + for ( ; expBox.More(); expBox.Next() ) + _nbShape++; // create bounding box for every shape inside the compound @@ -1024,8 +1063,8 @@ bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh, tabBox[i] = new double[6]; Standard_Real Xmin, Ymin, Zmin, Xmax, Ymax, Zmax; - TopExp_Explorer expBox (meshDS->ShapeToMesh(), TopAbs_SOLID); - for (; expBox.More(); expBox.Next()) { + // TopExp_Explorer expBox (meshDS->ShapeToMesh(), TopAbs_SOLID); -- 0020330:...ghs3d as a submesh + for (expBox.ReInit(); expBox.More(); expBox.Next()) { tabShape[iShape] = expBox.Current(); Bnd_Box BoundingBox; BRepBndLib::Add(expBox.Current(), BoundingBox);