X-Git-Url: http://git.salome-platform.org/gitweb/?a=blobdiff_plain;f=src%2FStdMeshers%2FStdMeshers_QuadToTriaAdaptor.cxx;h=54c3b7609f79c418c945cc1454a0095b68c07c0a;hb=4cd2499bddcd3da3ec8900fe825bc98669b789b5;hp=047cee6f62fd97b7ee5eacd32afb59d4d2d549e5;hpb=d8a28e48bbb291563352ba6b1e284e335bc8302a;p=modules%2Fsmesh.git diff --git a/src/StdMeshers/StdMeshers_QuadToTriaAdaptor.cxx b/src/StdMeshers/StdMeshers_QuadToTriaAdaptor.cxx index 047cee6f6..54c3b7609 100644 --- a/src/StdMeshers/StdMeshers_QuadToTriaAdaptor.cxx +++ b/src/StdMeshers/StdMeshers_QuadToTriaAdaptor.cxx @@ -1,34 +1,35 @@ -// Copyright (C) 2007-2010 CEA/DEN, EDF R&D, OPEN CASCADE +// Copyright (C) 2007-2012 CEA/DEN, EDF R&D, OPEN CASCADE // -// This library is free software; you can redistribute it and/or -// modify it under the terms of the GNU Lesser General Public -// License as published by the Free Software Foundation; either -// version 2.1 of the License. +// This library is free software; you can redistribute it and/or +// modify it under the terms of the GNU Lesser General Public +// License as published by the Free Software Foundation; either +// version 2.1 of the License. // -// This library is distributed in the hope that it will be useful, -// but WITHOUT ANY WARRANTY; without even the implied warranty of -// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU -// Lesser General Public License for more details. +// This library is distributed in the hope that it will be useful, +// but WITHOUT ANY WARRANTY; without even the implied warranty of +// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU +// Lesser General Public License for more details. // -// You should have received a copy of the GNU Lesser General Public -// License along with this library; if not, write to the Free Software -// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA +// You should have received a copy of the GNU Lesser General Public +// License along with this library; if not, write to the Free Software +// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA // -// See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com +// See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com // -// SMESH SMESH : implementaion of SMESH idl descriptions // File : StdMeshers_QuadToTriaAdaptor.cxx // Module : SMESH // Created : Wen May 07 16:37:07 2008 // Author : Sergey KUUL (skl) -// + #include "StdMeshers_QuadToTriaAdaptor.hxx" #include "SMDS_SetIterator.hxx" #include "SMESH_Algo.hxx" #include "SMESH_MesherHelper.hxx" +#include "SMESH_Group.hxx" +#include "SMESHDS_GroupBase.hxx" #include #include @@ -39,8 +40,11 @@ #include #include #include +#include "utilities.h" +#include #include +#include using namespace std; @@ -111,15 +115,14 @@ namespace // Check angle between normals double angle = nI.Angle( nJ ); - bool tooClose = ( angle < 15 * PI180 ); + bool tooClose = ( angle < 15. * M_PI / 180. ); // Check if pyramids collide - bool isOutI, isOutJ; if ( !tooClose && baI * baJ > 0 ) { // find out if nI points outside of PrmI or inside int dInd = baseNodesIndI[1] - baseNodesIndI[0]; - isOutI = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0; + bool isOutI = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0; // find out sign of projection of nJ to baI double proj = baI * nJ; @@ -133,12 +136,15 @@ namespace // check order of baseNodes within pyramids, it must be opposite int dInd; dInd = baseNodesIndI[1] - baseNodesIndI[0]; - isOutI = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0; + bool isOutI = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0; dInd = baseNodesIndJ[1] - baseNodesIndJ[0]; - isOutJ = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0; + bool isOutJ = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0; if ( isOutJ == isOutI ) return false; // other domain + // direct both normals outside pyramid + ( isOutI ? nJ : nI ).Reverse(); + // check absence of a face separating domains between pyramids TIDSortedElemSet emptySet, avoidSet; int i1, i2; @@ -153,6 +159,9 @@ namespace while ( otherNodeInd == i1 || otherNodeInd == i2 ) otherNodeInd++; const SMDS_MeshNode* otherFaceNode = f->GetNode( otherNodeInd ); + if ( otherFaceNode == nApexI || otherFaceNode == nApexJ ) + continue; // f is a temporary triangle + // check if f is a base face of either of pyramids if ( f->NbCornerNodes() == 4 && ( PrmI->GetNodeIndex( otherFaceNode ) >= 0 || @@ -161,7 +170,6 @@ namespace // check projections of face direction (baOFN) to triange normals (nI and nJ) gp_Vec baOFN( base1, SMESH_TNodeXYZ( otherFaceNode )); - ( isOutI ? nJ : nI ).Reverse(); if ( nI * baOFN > 0 && nJ * baOFN > 0 ) { tooClose = false; // f is between pyramids @@ -172,6 +180,79 @@ namespace return tooClose; } + + //================================================================================ + /*! + * \brief Move medium nodes of merged quadratic pyramids + */ + //================================================================================ + + void UpdateQuadraticPyramids(const set& commonApex, + SMESHDS_Mesh* meshDS) + { + typedef SMDS_StdIterator< const SMDS_MeshElement*, SMDS_ElemIteratorPtr > TStdElemIterator; + TStdElemIterator itEnd; + + // shift of node index to get medium nodes between the 4 base nodes and the apex + const int base2MediumShift = 9; + + set::const_iterator nIt = commonApex.begin(); + for ( ; nIt != commonApex.end(); ++nIt ) + { + SMESH_TNodeXYZ apex( *nIt ); + + vector< const SMDS_MeshElement* > pyrams // pyramids sharing the apex node + ( TStdElemIterator( apex._node->GetInverseElementIterator( SMDSAbs_Volume )), itEnd ); + + // Select medium nodes to keep and medium nodes to remove + + typedef map < const SMDS_MeshNode*, const SMDS_MeshNode*, TIDCompare > TN2NMap; + TN2NMap base2medium; // to keep + vector< const SMDS_MeshNode* > nodesToRemove; + + for ( unsigned i = 0; i < pyrams.size(); ++i ) + for ( int baseIndex = 0; baseIndex < PYRAM_APEX; ++baseIndex ) + { + SMESH_TNodeXYZ base = pyrams[i]->GetNode( baseIndex ); + const SMDS_MeshNode* medium = pyrams[i]->GetNode( baseIndex + base2MediumShift ); + TN2NMap::iterator b2m = base2medium.insert( make_pair( base._node, medium )).first; + if ( b2m->second != medium ) + { + nodesToRemove.push_back( medium ); + } + else + { + // move the kept medium node + gp_XYZ newXYZ = 0.5 * ( apex + base ); + meshDS->MoveNode( medium, newXYZ.X(), newXYZ.Y(), newXYZ.Z() ); + } + } + + // Within pyramids, replace nodes to remove by nodes to keep + + for ( unsigned i = 0; i < pyrams.size(); ++i ) + { + vector< const SMDS_MeshNode* > nodes( pyrams[i]->begin_nodes(), + pyrams[i]->end_nodes() ); + for ( int baseIndex = 0; baseIndex < PYRAM_APEX; ++baseIndex ) + { + const SMDS_MeshNode* base = pyrams[i]->GetNode( baseIndex ); + nodes[ baseIndex + base2MediumShift ] = base2medium[ base ]; + } + meshDS->ChangeElementNodes( pyrams[i], &nodes[0], nodes.size()); + } + + // Remove the replaced nodes + + if ( !nodesToRemove.empty() ) + { + SMESHDS_SubMesh * sm = meshDS->MeshElements( nodesToRemove[0]->getshapeId() ); + for ( unsigned i = 0; i < nodesToRemove.size(); ++i ) + meshDS->RemoveFreeNode( nodesToRemove[i], sm, /*fromGroups=*/false); + } + } + } + } //================================================================================ @@ -204,7 +285,7 @@ void StdMeshers_QuadToTriaAdaptor::MergePiramids( const SMDS_MeshElement* Pr // find and remove coincided faces of merged pyramids vector< const SMDS_MeshElement* > inverseElems - // copy inverse elements to avoid iteration on changing container + // copy inverse elements to avoid iteration on changing container ( TStdElemIterator( CommonNode->GetInverseElementIterator(SMDSAbs_Face)), itEnd); for ( unsigned i = 0; i < inverseElems.size(); ++i ) { @@ -302,37 +383,41 @@ StdMeshers_QuadToTriaAdaptor::~StdMeshers_QuadToTriaAdaptor() myElemSearcher=0; } - //======================================================================= //function : FindBestPoint //purpose : Return a point P laying on the line (PC,V) so that triangle // (P, P1, P2) to be equilateral as much as possible // V - normal to (P1,P2,PC) //======================================================================= + static gp_Pnt FindBestPoint(const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& PC, const gp_Vec& V) { - double a = P1.Distance(P2); - double b = P1.Distance(PC); - double c = P2.Distance(PC); + gp_Pnt Pbest = PC; + const double a = P1.Distance(P2); + const double b = P1.Distance(PC); + const double c = P2.Distance(PC); if( a < (b+c)/2 ) - return PC; + return Pbest; else { // find shift along V in order a to became equal to (b+c)/2 - double shift = sqrt( a*a + (b*b-c*c)*(b*b-c*c)/16/a/a - (b*b+c*c)/2 ); - gp_Dir aDir(V); - gp_Pnt Pbest = PC.XYZ() + aDir.XYZ() * shift; - return Pbest; + const double Vsize = V.Magnitude(); + if ( fabs( Vsize ) > std::numeric_limits::min() ) + { + const double shift = sqrt( a*a + (b*b-c*c)*(b*b-c*c)/16/a/a - (b*b+c*c)/2 ); + Pbest.ChangeCoord() += shift * V.XYZ() / Vsize; + } } + return Pbest; } - //======================================================================= //function : HasIntersection3 //purpose : Auxilare for HasIntersection() // find intersection point between triangle (P1,P2,P3) // and segment [PC,P] //======================================================================= + static bool HasIntersection3(const gp_Pnt& P, const gp_Pnt& PC, gp_Pnt& Pint, const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3) { @@ -403,7 +488,6 @@ static bool HasIntersection3(const gp_Pnt& P, const gp_Pnt& PC, gp_Pnt& Pint, return false; } - //======================================================================= //function : HasIntersection //purpose : Auxilare for CheckIntersection() @@ -470,13 +554,13 @@ bool StdMeshers_QuadToTriaAdaptor::CheckIntersection (const gp_Pnt& P, gp_Ax1 line( P, gp_Vec(P,PC)); vector< const SMDS_MeshElement* > suspectElems; searcher->GetElementsNearLine( line, SMDSAbs_Face, suspectElems); - + for ( int i = 0; i < suspectElems.size(); ++i ) { const SMDS_MeshElement* face = suspectElems[i]; if ( face == NotCheckedFace ) continue; Handle(TColgp_HSequenceOfPnt) aContour = new TColgp_HSequenceOfPnt; - for ( int i = 0; i < face->NbCornerNodes(); ++i ) + for ( int i = 0; i < face->NbCornerNodes(); ++i ) aContour->Append( SMESH_TNodeXYZ( face->GetNode(i) )); if( HasIntersection(P, PC, Pres, aContour) ) { res = true; @@ -622,12 +706,12 @@ int StdMeshers_QuadToTriaAdaptor::Preparation(const SMDS_MeshElement* face //======================================================================= //function : Compute -//purpose : +//purpose : //======================================================================= -bool StdMeshers_QuadToTriaAdaptor::Compute(SMESH_Mesh& aMesh, - const TopoDS_Shape& aShape, - SMESH_ProxyMesh* aProxyMesh) +bool StdMeshers_QuadToTriaAdaptor::Compute(SMESH_Mesh& aMesh, + const TopoDS_Shape& aShape, + SMESH_ProxyMesh* aProxyMesh) { SMESH_ProxyMesh::setMesh( aMesh ); @@ -635,6 +719,8 @@ bool StdMeshers_QuadToTriaAdaptor::Compute(SMESH_Mesh& aMesh, aShape.ShapeType() != TopAbs_SHELL ) return false; + myShape = aShape; + vector myPyramids; SMESHDS_Mesh * meshDS = aMesh.GetMeshDS(); @@ -670,7 +756,7 @@ bool StdMeshers_QuadToTriaAdaptor::Compute(SMESH_Mesh& aMesh, { bool isRev = false; if ( helper.NbAncestors( aShapeFace, aMesh, aShape.ShapeType() ) > 1 ) - isRev = SMESH_Algo::IsReversedSubMesh( TopoDS::Face(aShapeFace), meshDS ); + isRev = helper.IsReversedSubMesh( TopoDS::Face(aShapeFace) ); SMDS_ElemIteratorPtr iteratorElem = aSubMeshDSFace->GetElements(); while ( iteratorElem->more() ) // loop on elements on a geometrical face @@ -816,6 +902,36 @@ bool StdMeshers_QuadToTriaAdaptor::Compute(SMESH_Mesh& aMesh) if ( aMesh.NbQuadrangles() < 1 ) return false; + // find if there is a group of faces identified as skin faces, with normal going outside the volume + std::string groupName = "skinFaces"; + SMESHDS_GroupBase* groupDS = 0; + SMESH_Mesh::GroupIteratorPtr groupIt = aMesh.GetGroups(); + while ( groupIt->more() ) + { + groupDS = 0; + SMESH_Group * group = groupIt->next(); + if ( !group ) continue; + groupDS = group->GetGroupDS(); + if ( !groupDS || groupDS->IsEmpty() ) + { + groupDS = 0; + continue; + } + if (groupDS->GetType() != SMDSAbs_Face) + { + groupDS = 0; + continue; + } + std::string grpName = group->GetName(); + if (grpName == groupName) + { + MESSAGE("group skinFaces provided"); + break; + } + else + groupDS = 0; + } + vector myPyramids; SMESH_MesherHelper helper(aMesh); helper.IsQuadraticSubMesh(aMesh.GetShapeToMesh()); @@ -829,7 +945,7 @@ bool StdMeshers_QuadToTriaAdaptor::Compute(SMESH_Mesh& aMesh) SMESH_ProxyMesh::SubMesh* prxSubMesh = getProxySubMesh(); SMDS_FaceIteratorPtr fIt = meshDS->facesIterator(/*idInceasingOrder=*/true); - while( fIt->more()) + while( fIt->more()) { const SMDS_MeshElement* face = fIt->next(); if ( !face ) continue; @@ -921,7 +1037,7 @@ bool StdMeshers_QuadToTriaAdaptor::Compute(SMESH_Mesh& aMesh) continue; } - // Case of non-degenerated quadrangle + // Case of non-degenerated quadrangle // Find pyramid peak @@ -934,10 +1050,12 @@ bool StdMeshers_QuadToTriaAdaptor::Compute(SMESH_Mesh& aMesh) PCbest /= 4; double height = PC.Distance(PCbest); // pyramid height to precise - if(height<1.e-6) { + if ( height < 1.e-6 ) { // create new PCbest using a bit shift along VNorm PCbest = PC.XYZ() + VNorm.XYZ() * 0.001; height = PC.Distance(PCbest); + if ( height < std::numeric_limits::min() ) + return false; // batterfly element } // Restrict pyramid height by intersection with other faces @@ -978,6 +1096,10 @@ bool StdMeshers_QuadToTriaAdaptor::Compute(SMESH_Mesh& aMesh) } } + // if the face belong to the group of skinFaces, do not build a pyramid outside + if (groupDS && groupDS->Contains(face)) + intersected[0] = false; + // Create one or two pyramids for ( int isRev = 0; isRev < 2; ++isRev ) @@ -1027,8 +1149,8 @@ bool StdMeshers_QuadToTriaAdaptor::Compute2ndPart(SMESH_Mesh& SMESHDS_Mesh * meshDS = aMesh.GetMeshDS(); int i, j, k, myShapeID = myPyramids[0]->GetNode(4)->getshapeId(); - if ( !myElemSearcher ) - myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher(); + if ( myElemSearcher ) delete myElemSearcher; + myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher(); SMESH_ElementSearcher* searcher = const_cast(myElemSearcher); set nodesToMove; @@ -1088,11 +1210,11 @@ bool StdMeshers_QuadToTriaAdaptor::Compute2ndPart(SMESH_Mesh& vector PsJ( xyzIt, TXyzIterator() ); gp_Pnt Pint; - bool hasInt; + bool hasInt=false; for(k=0; k<4 && !hasInt; k++) { gp_Vec Vtmp(PsI[k],PsI[4]); gp_Pnt Pshift = PsI[k].XYZ() + Vtmp.XYZ() * 0.01; // base node moved a bit to apex - hasInt = + hasInt = ( HasIntersection3( Pshift, PsI[4], Pint, PsJ[0], PsJ[1], PsJ[4]) || HasIntersection3( Pshift, PsI[4], Pint, PsJ[1], PsJ[2], PsJ[4]) || HasIntersection3( Pshift, PsI[4], Pint, PsJ[2], PsJ[3], PsJ[4]) || @@ -1101,7 +1223,7 @@ bool StdMeshers_QuadToTriaAdaptor::Compute2ndPart(SMESH_Mesh& for(k=0; k<4 && !hasInt; k++) { gp_Vec Vtmp(PsJ[k],PsJ[4]); gp_Pnt Pshift = PsJ[k].XYZ() + Vtmp.XYZ() * 0.01; - hasInt = + hasInt = ( HasIntersection3( Pshift, PsJ[4], Pint, PsI[0], PsI[1], PsI[4]) || HasIntersection3( Pshift, PsJ[4], Pint, PsI[1], PsI[2], PsI[4]) || HasIntersection3( Pshift, PsJ[4], Pint, PsI[2], PsI[3], PsI[4]) || @@ -1131,15 +1253,15 @@ bool StdMeshers_QuadToTriaAdaptor::Compute2ndPart(SMESH_Mesh& PCi += PsI[k].XYZ(); PCj += PsJ[k].XYZ(); } - PCi /= 4; PCj /= 4; + PCi /= 4; PCj /= 4; gp_Vec VN1(PCi,PsI[4]); gp_Vec VN2(PCj,PsJ[4]); gp_Vec VI1(PCi,Pint); gp_Vec VI2(PCj,Pint); double ang1 = fabs(VN1.Angle(VI1)); double ang2 = fabs(VN2.Angle(VI2)); - double coef1 = 0.5 - (( ang1MoveNode( *n, (*n)->X(), (*n)->Y(), (*n)->Z() ); } + // move medium nodes of merged quadratic pyramids + if ( myPyramids[0]->IsQuadratic() ) + UpdateQuadraticPyramids( nodesToMove, GetMeshDS() ); + // erase removed triangles from the proxy mesh if ( !myRemovedTrias.empty() ) { @@ -1199,15 +1325,3 @@ bool StdMeshers_QuadToTriaAdaptor::Compute2ndPart(SMESH_Mesh& return true; } - -//================================================================================ -/*! - * \brief Return list of created triangles for given face - */ -//================================================================================ - -// const list* StdMeshers_QuadToTriaAdaptor::GetTriangles (const SMDS_MeshElement* aQuad) -// { -// TQuad2Trias::iterator it = myResMap.find(aQuad); -// return ( it != myResMap.end() ? & it->second : 0 ); -// }