1 // Copyright (C) 2007-2010 CEA/DEN, EDF R&D, OPEN CASCADE
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License.
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
11 // Lesser General Public License for more details.
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
20 // SMESH SMESH : implementaion of SMESH idl descriptions
21 // File : StdMeshers_QuadToTriaAdaptor.cxx
23 // Created : Wen May 07 16:37:07 2008
24 // Author : Sergey KUUL (skl)
26 #include "StdMeshers_QuadToTriaAdaptor.hxx"
28 #include "SMDS_SetIterator.hxx"
30 #include "SMESH_Algo.hxx"
31 #include "SMESH_MesherHelper.hxx"
33 #include <IntAna_IntConicQuad.hxx>
34 #include <IntAna_Quadric.hxx>
35 #include <TColgp_HArray1OfPnt.hxx>
36 #include <TColgp_HArray1OfVec.hxx>
37 #include <TColgp_HSequenceOfPnt.hxx>
38 #include <TopExp_Explorer.hxx>
47 enum EQuadNature { NOT_QUAD, QUAD, DEGEN_QUAD, PYRAM_APEX = 4, TRIA_APEX = 0 };
49 // std-like iterator used to get coordinates of nodes of mesh element
50 typedef SMDS_StdIterator< SMESH_TNodeXYZ, SMDS_ElemIteratorPtr > TXyzIterator;
54 //================================================================================
56 * \brief Return true if two nodes of triangles are equal
58 //================================================================================
60 bool EqualTriangles(const SMDS_MeshElement* F1,const SMDS_MeshElement* F2)
63 ( F1->GetNode(1)==F2->GetNode(2) && F1->GetNode(2)==F2->GetNode(1) ) ||
64 ( F1->GetNode(1)==F2->GetNode(1) && F1->GetNode(2)==F2->GetNode(2) );
66 //================================================================================
68 * \brief Return true if two adjacent pyramids are too close one to another
69 * so that a tetrahedron to built between them would have too poor quality
71 //================================================================================
73 bool TooCloseAdjacent( const SMDS_MeshElement* PrmI,
74 const SMDS_MeshElement* PrmJ,
77 const SMDS_MeshNode* nApexI = PrmI->GetNode(4);
78 const SMDS_MeshNode* nApexJ = PrmJ->GetNode(4);
79 if ( nApexI == nApexJ ||
80 nApexI->getshapeId() != nApexJ->getshapeId() )
83 // Find two common base nodes and their indices within PrmI and PrmJ
84 const SMDS_MeshNode* baseNodes[2] = { 0,0 };
85 int baseNodesIndI[2], baseNodesIndJ[2];
86 for ( int i = 0; i < 4 ; ++i )
88 int j = PrmJ->GetNodeIndex( PrmI->GetNode(i));
91 int ind = baseNodes[0] ? 1:0;
92 if ( baseNodes[ ind ])
93 return false; // pyramids with a common base face
94 baseNodes [ ind ] = PrmI->GetNode(i);
95 baseNodesIndI[ ind ] = i;
96 baseNodesIndJ[ ind ] = j;
99 if ( !baseNodes[1] ) return false; // not adjacent
101 // Get normals of triangles sharing baseNodes
102 gp_XYZ apexI = SMESH_TNodeXYZ( nApexI );
103 gp_XYZ apexJ = SMESH_TNodeXYZ( nApexJ );
104 gp_XYZ base1 = SMESH_TNodeXYZ( baseNodes[0]);
105 gp_XYZ base2 = SMESH_TNodeXYZ( baseNodes[1]);
106 gp_Vec baseVec( base1, base2 );
107 gp_Vec baI( base1, apexI );
108 gp_Vec baJ( base1, apexJ );
109 gp_Vec nI = baseVec.Crossed( baI );
110 gp_Vec nJ = baseVec.Crossed( baJ );
112 // Check angle between normals
113 double angle = nI.Angle( nJ );
114 bool tooClose = ( angle < 15 * PI180 );
116 // Check if pyramids collide
117 if ( !tooClose && baI * baJ > 0 )
119 // find out if nI points outside of PrmI or inside
120 int dInd = baseNodesIndI[1] - baseNodesIndI[0];
121 bool isOutI = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0;
123 // find out sign of projection of nJ to baI
124 double proj = baI * nJ;
126 tooClose = isOutI ? proj > 0 : proj < 0;
129 // Check if PrmI and PrmJ are in same domain
130 if ( tooClose && !hasShape )
132 // check order of baseNodes within pyramids, it must be opposite
134 dInd = baseNodesIndI[1] - baseNodesIndI[0];
135 bool isOutI = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0;
136 dInd = baseNodesIndJ[1] - baseNodesIndJ[0];
137 bool isOutJ = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0;
138 if ( isOutJ == isOutI )
139 return false; // other domain
141 // direct both normals outside pyramid
142 ( isOutI ? nJ : nI ).Reverse();
144 // check absence of a face separating domains between pyramids
145 TIDSortedElemSet emptySet, avoidSet;
147 while ( const SMDS_MeshElement* f =
148 SMESH_MeshEditor::FindFaceInSet( baseNodes[0], baseNodes[1],
149 emptySet, avoidSet, &i1, &i2 ))
151 avoidSet.insert( f );
153 // face node other than baseNodes
154 int otherNodeInd = 0;
155 while ( otherNodeInd == i1 || otherNodeInd == i2 ) otherNodeInd++;
156 const SMDS_MeshNode* otherFaceNode = f->GetNode( otherNodeInd );
158 if ( otherFaceNode == nApexI || otherFaceNode == nApexJ )
159 continue; // f is a temporary triangle
161 // check if f is a base face of either of pyramids
162 if ( f->NbCornerNodes() == 4 &&
163 ( PrmI->GetNodeIndex( otherFaceNode ) >= 0 ||
164 PrmJ->GetNodeIndex( otherFaceNode ) >= 0 ))
165 continue; // f is a base quadrangle
167 // check projections of face direction (baOFN) to triange normals (nI and nJ)
168 gp_Vec baOFN( base1, SMESH_TNodeXYZ( otherFaceNode ));
169 if ( nI * baOFN > 0 && nJ * baOFN > 0 )
171 tooClose = false; // f is between pyramids
181 //================================================================================
183 * \brief Merge the two pyramids (i.e. fuse their apex) and others already merged with them
185 //================================================================================
187 void StdMeshers_QuadToTriaAdaptor::MergePiramids( const SMDS_MeshElement* PrmI,
188 const SMDS_MeshElement* PrmJ,
189 set<const SMDS_MeshNode*> & nodesToMove)
191 const SMDS_MeshNode* Nrem = PrmJ->GetNode(4); // node to remove
192 //int nbJ = Nrem->NbInverseElements( SMDSAbs_Volume );
193 SMESH_TNodeXYZ Pj( Nrem );
195 // an apex node to make common to all merged pyramids
196 SMDS_MeshNode* CommonNode = const_cast<SMDS_MeshNode*>(PrmI->GetNode(4));
197 if ( CommonNode == Nrem ) return; // already merged
198 //int nbI = CommonNode->NbInverseElements( SMDSAbs_Volume );
199 SMESH_TNodeXYZ Pi( CommonNode );
200 gp_XYZ Pnew = /*( nbI*Pi + nbJ*Pj ) / (nbI+nbJ);*/ 0.5 * ( Pi + Pj );
201 CommonNode->setXYZ( Pnew.X(), Pnew.Y(), Pnew.Z() );
203 nodesToMove.insert( CommonNode );
204 nodesToMove.erase ( Nrem );
206 typedef SMDS_StdIterator< const SMDS_MeshElement*, SMDS_ElemIteratorPtr > TStdElemIterator;
207 TStdElemIterator itEnd;
209 // find and remove coincided faces of merged pyramids
210 vector< const SMDS_MeshElement* > inverseElems
211 // copy inverse elements to avoid iteration on changing container
212 ( TStdElemIterator( CommonNode->GetInverseElementIterator(SMDSAbs_Face)), itEnd);
213 for ( unsigned i = 0; i < inverseElems.size(); ++i )
215 const SMDS_MeshElement* FI = inverseElems[i];
216 const SMDS_MeshElement* FJEqual = 0;
217 SMDS_ElemIteratorPtr triItJ = Nrem->GetInverseElementIterator(SMDSAbs_Face);
218 while ( !FJEqual && triItJ->more() )
220 const SMDS_MeshElement* FJ = triItJ->next();
221 if ( EqualTriangles( FJ, FI ))
226 removeTmpElement( FI );
227 removeTmpElement( FJEqual );
228 myRemovedTrias.insert( FI );
229 myRemovedTrias.insert( FJEqual );
233 // set the common apex node to pyramids and triangles merged with J
234 inverseElems.assign( TStdElemIterator( Nrem->GetInverseElementIterator()), itEnd );
235 for ( unsigned i = 0; i < inverseElems.size(); ++i )
237 const SMDS_MeshElement* elem = inverseElems[i];
238 vector< const SMDS_MeshNode* > nodes( elem->begin_nodes(), elem->end_nodes() );
239 nodes[ elem->GetType() == SMDSAbs_Volume ? PYRAM_APEX : TRIA_APEX ] = CommonNode;
240 GetMeshDS()->ChangeElementNodes( elem, &nodes[0], nodes.size());
242 ASSERT( Nrem->NbInverseElements() == 0 );
243 GetMeshDS()->RemoveFreeNode( Nrem,
244 GetMeshDS()->MeshElements( Nrem->getshapeId()),
245 /*fromGroups=*/false);
248 //================================================================================
250 * \brief Merges adjacent pyramids
252 //================================================================================
254 void StdMeshers_QuadToTriaAdaptor::MergeAdjacent(const SMDS_MeshElement* PrmI,
255 set<const SMDS_MeshNode*>& nodesToMove)
257 TIDSortedElemSet adjacentPyrams;
258 bool mergedPyrams = false;
259 for(int k=0; k<4; k++) // loop on 4 base nodes of PrmI
261 const SMDS_MeshNode* n = PrmI->GetNode(k);
262 SMDS_ElemIteratorPtr vIt = n->GetInverseElementIterator( SMDSAbs_Volume );
263 while ( vIt->more() )
265 const SMDS_MeshElement* PrmJ = vIt->next();
266 if ( PrmJ->NbCornerNodes() != 5 || !adjacentPyrams.insert( PrmJ ).second )
268 if ( PrmI != PrmJ && TooCloseAdjacent( PrmI, PrmJ, GetMesh()->HasShapeToMesh() ))
270 MergePiramids( PrmI, PrmJ, nodesToMove );
272 // container of inverse elements can change
273 vIt = n->GetInverseElementIterator( SMDSAbs_Volume );
279 TIDSortedElemSet::iterator prm;
280 for (prm = adjacentPyrams.begin(); prm != adjacentPyrams.end(); ++prm)
281 MergeAdjacent( *prm, nodesToMove );
285 //================================================================================
289 //================================================================================
291 StdMeshers_QuadToTriaAdaptor::StdMeshers_QuadToTriaAdaptor():
296 //================================================================================
300 //================================================================================
302 StdMeshers_QuadToTriaAdaptor::~StdMeshers_QuadToTriaAdaptor()
304 // temporary faces are deleted by ~SMESH_ProxyMesh()
305 if ( myElemSearcher ) delete myElemSearcher;
310 //=======================================================================
311 //function : FindBestPoint
312 //purpose : Return a point P laying on the line (PC,V) so that triangle
313 // (P, P1, P2) to be equilateral as much as possible
314 // V - normal to (P1,P2,PC)
315 //=======================================================================
316 static gp_Pnt FindBestPoint(const gp_Pnt& P1, const gp_Pnt& P2,
317 const gp_Pnt& PC, const gp_Vec& V)
319 double a = P1.Distance(P2);
320 double b = P1.Distance(PC);
321 double c = P2.Distance(PC);
325 // find shift along V in order a to became equal to (b+c)/2
326 double shift = sqrt( a*a + (b*b-c*c)*(b*b-c*c)/16/a/a - (b*b+c*c)/2 );
328 gp_Pnt Pbest = PC.XYZ() + aDir.XYZ() * shift;
334 //=======================================================================
335 //function : HasIntersection3
336 //purpose : Auxilare for HasIntersection()
337 // find intersection point between triangle (P1,P2,P3)
338 // and segment [PC,P]
339 //=======================================================================
340 static bool HasIntersection3(const gp_Pnt& P, const gp_Pnt& PC, gp_Pnt& Pint,
341 const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3)
343 //cout<<"HasIntersection3"<<endl;
344 //cout<<" PC("<<PC.X()<<","<<PC.Y()<<","<<PC.Z()<<")"<<endl;
345 //cout<<" P("<<P.X()<<","<<P.Y()<<","<<P.Z()<<")"<<endl;
346 //cout<<" P1("<<P1.X()<<","<<P1.Y()<<","<<P1.Z()<<")"<<endl;
347 //cout<<" P2("<<P2.X()<<","<<P2.Y()<<","<<P2.Z()<<")"<<endl;
348 //cout<<" P3("<<P3.X()<<","<<P3.Y()<<","<<P3.Z()<<")"<<endl;
351 IntAna_Quadric IAQ(gp_Pln(P1,VP1.Crossed(VP2)));
352 IntAna_IntConicQuad IAICQ(gp_Lin(PC,gp_Dir(gp_Vec(PC,P))),IAQ);
354 if( IAICQ.IsInQuadric() )
356 if( IAICQ.NbPoints() == 1 ) {
357 gp_Pnt PIn = IAICQ.Point(1);
358 const double preci = 1.e-10 * P.Distance(PC);
359 // check if this point is internal for segment [PC,P]
361 ( (PC.X()-PIn.X())*(P.X()-PIn.X()) > preci ) ||
362 ( (PC.Y()-PIn.Y())*(P.Y()-PIn.Y()) > preci ) ||
363 ( (PC.Z()-PIn.Z())*(P.Z()-PIn.Z()) > preci );
367 // check if this point is internal for triangle (P1,P2,P3)
371 if( V1.Magnitude()<preci ||
372 V2.Magnitude()<preci ||
373 V3.Magnitude()<preci ) {
377 const double angularTol = 1e-6;
378 gp_Vec VC1 = V1.Crossed(V2);
379 gp_Vec VC2 = V2.Crossed(V3);
380 gp_Vec VC3 = V3.Crossed(V1);
381 if(VC1.Magnitude()<gp::Resolution()) {
382 if(VC2.IsOpposite(VC3,angularTol)) {
386 else if(VC2.Magnitude()<gp::Resolution()) {
387 if(VC1.IsOpposite(VC3,angularTol)) {
391 else if(VC3.Magnitude()<gp::Resolution()) {
392 if(VC1.IsOpposite(VC2,angularTol)) {
397 if( VC1.IsOpposite(VC2,angularTol) || VC1.IsOpposite(VC3,angularTol) ||
398 VC2.IsOpposite(VC3,angularTol) ) {
411 //=======================================================================
412 //function : HasIntersection
413 //purpose : Auxilare for CheckIntersection()
414 //=======================================================================
416 static bool HasIntersection(const gp_Pnt& P, const gp_Pnt& PC, gp_Pnt& Pint,
417 Handle(TColgp_HSequenceOfPnt)& aContour)
419 if(aContour->Length()==3) {
420 return HasIntersection3( P, PC, Pint, aContour->Value(1),
421 aContour->Value(2), aContour->Value(3) );
425 if( (aContour->Value(1).Distance(aContour->Value(2)) > 1.e-6) &&
426 (aContour->Value(1).Distance(aContour->Value(3)) > 1.e-6) &&
427 (aContour->Value(2).Distance(aContour->Value(3)) > 1.e-6) ) {
428 check = HasIntersection3( P, PC, Pint, aContour->Value(1),
429 aContour->Value(2), aContour->Value(3) );
431 if(check) return true;
432 if( (aContour->Value(1).Distance(aContour->Value(4)) > 1.e-6) &&
433 (aContour->Value(1).Distance(aContour->Value(3)) > 1.e-6) &&
434 (aContour->Value(4).Distance(aContour->Value(3)) > 1.e-6) ) {
435 check = HasIntersection3( P, PC, Pint, aContour->Value(1),
436 aContour->Value(3), aContour->Value(4) );
438 if(check) return true;
444 //================================================================================
446 * \brief Checks if a line segment (P,PC) intersects any mesh face.
447 * \param P - first segment end
448 * \param PC - second segment end (it is a gravity center of quadrangle)
449 * \param Pint - (out) intersection point
450 * \param aMesh - mesh
451 * \param aShape - shape to check faces on
452 * \param NotCheckedFace - mesh face not to check
453 * \retval bool - true if there is an intersection
455 //================================================================================
457 bool StdMeshers_QuadToTriaAdaptor::CheckIntersection (const gp_Pnt& P,
461 const TopoDS_Shape& aShape,
462 const SMDS_MeshElement* NotCheckedFace)
464 if ( !myElemSearcher )
465 myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher();
466 SMESH_ElementSearcher* searcher = const_cast<SMESH_ElementSearcher*>(myElemSearcher);
468 //SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
469 //cout<<" CheckIntersection: meshDS->NbFaces() = "<<meshDS->NbFaces()<<endl;
471 double dist = RealLast(); // find intersection closest to the segment
474 gp_Ax1 line( P, gp_Vec(P,PC));
475 vector< const SMDS_MeshElement* > suspectElems;
476 searcher->GetElementsNearLine( line, SMDSAbs_Face, suspectElems);
478 for ( int i = 0; i < suspectElems.size(); ++i )
480 const SMDS_MeshElement* face = suspectElems[i];
481 if ( face == NotCheckedFace ) continue;
482 Handle(TColgp_HSequenceOfPnt) aContour = new TColgp_HSequenceOfPnt;
483 for ( int i = 0; i < face->NbCornerNodes(); ++i )
484 aContour->Append( SMESH_TNodeXYZ( face->GetNode(i) ));
485 if( HasIntersection(P, PC, Pres, aContour) ) {
487 double tmp = PC.Distance(Pres);
497 //================================================================================
499 * \brief Prepare data for the given face
500 * \param PN - coordinates of face nodes
501 * \param VN - cross products of vectors (PC-PN(i)) ^ (PC-PN(i+1))
502 * \param FNodes - face nodes
503 * \param PC - gravity center of nodes
504 * \param VNorm - face normal (sum of VN)
505 * \param volumes - two volumes sharing the given face, the first is in VNorm direction
506 * \retval int - 0 if given face is not quad,
507 * 1 if given face is quad,
508 * 2 if given face is degenerate quad (two nodes are coincided)
510 //================================================================================
512 int StdMeshers_QuadToTriaAdaptor::Preparation(const SMDS_MeshElement* face,
513 Handle(TColgp_HArray1OfPnt)& PN,
514 Handle(TColgp_HArray1OfVec)& VN,
515 vector<const SMDS_MeshNode*>& FNodes,
518 const SMDS_MeshElement** volumes)
520 if( face->NbCornerNodes() != 4 )
526 gp_XYZ xyzC(0., 0., 0.);
527 for ( i = 0; i < 4; ++i )
529 gp_XYZ p = SMESH_TNodeXYZ( FNodes[i] = face->GetNode(i) );
530 PN->SetValue( i+1, p );
541 if( PN->Value(i).Distance(PN->Value(j)) < 1.e-6 )
546 //int deg_num = IsDegenarate(PN);
550 //cout<<"find degeneration"<<endl;
552 gp_Pnt Pdeg = PN->Value(i);
554 list< const SMDS_MeshNode* >::iterator itdg = myDegNodes.begin();
555 const SMDS_MeshNode* DegNode = 0;
556 for(; itdg!=myDegNodes.end(); itdg++) {
557 const SMDS_MeshNode* N = (*itdg);
558 gp_Pnt Ptmp(N->X(),N->Y(),N->Z());
559 if(Pdeg.Distance(Ptmp)<1.e-6) {
561 //DegNode = const_cast<SMDS_MeshNode*>(N);
566 DegNode = FNodes[i-1];
567 myDegNodes.push_back(DegNode);
570 FNodes[i-1] = DegNode;
573 PN->SetValue(i,PN->Value(i+1));
574 FNodes[i-1] = FNodes[i];
579 PN->SetValue(nbp+1,PN->Value(1));
580 FNodes[nbp] = FNodes[0];
581 // find normal direction
582 gp_Vec V1(PC,PN->Value(nbp));
583 gp_Vec V2(PC,PN->Value(1));
584 VNorm = V1.Crossed(V2);
585 VN->SetValue(nbp,VNorm);
586 for(i=1; i<nbp; i++) {
587 V1 = gp_Vec(PC,PN->Value(i));
588 V2 = gp_Vec(PC,PN->Value(i+1));
589 gp_Vec Vtmp = V1.Crossed(V2);
590 VN->SetValue(i,Vtmp);
594 // find volumes sharing the face
597 volumes[0] = volumes[1] = 0;
598 SMDS_ElemIteratorPtr vIt = FNodes[0]->GetInverseElementIterator( SMDSAbs_Volume );
599 while ( vIt->more() )
601 const SMDS_MeshElement* vol = vIt->next();
602 bool volSharesAllNodes = true;
603 for ( int i = 1; i < face->NbNodes() && volSharesAllNodes; ++i )
604 volSharesAllNodes = ( vol->GetNodeIndex( FNodes[i] ) >= 0 );
605 if ( volSharesAllNodes )
606 volumes[ volumes[0] ? 1 : 0 ] = vol;
607 // we could additionally check that vol has all FNodes in its one face using SMDS_VolumeTool
609 // define volume position relating to the face normal
613 SMDS_ElemIteratorPtr nodeIt = volumes[0]->nodesIterator();
615 volGC = accumulate( TXyzIterator(nodeIt), TXyzIterator(), volGC ) / volumes[0]->NbNodes();
617 if ( VNorm * gp_Vec( PC, volGC ) < 0 )
618 swap( volumes[0], volumes[1] );
622 //cout<<" VNorm("<<VNorm.X()<<","<<VNorm.Y()<<","<<VNorm.Z()<<")"<<endl;
623 return hasdeg ? DEGEN_QUAD : QUAD;
627 //=======================================================================
630 //=======================================================================
632 bool StdMeshers_QuadToTriaAdaptor::Compute(SMESH_Mesh& aMesh,
633 const TopoDS_Shape& aShape,
634 SMESH_ProxyMesh* aProxyMesh)
636 SMESH_ProxyMesh::setMesh( aMesh );
638 if ( aShape.ShapeType() != TopAbs_SOLID &&
639 aShape.ShapeType() != TopAbs_SHELL )
642 vector<const SMDS_MeshElement*> myPyramids;
644 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
645 SMESH_MesherHelper helper(aMesh);
646 helper.IsQuadraticSubMesh(aShape);
647 helper.SetElementsOnShape( true );
649 if ( myElemSearcher ) delete myElemSearcher;
651 myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher( aProxyMesh->GetFaces(aShape));
653 myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher();
655 const SMESHDS_SubMesh * aSubMeshDSFace;
656 Handle(TColgp_HArray1OfPnt) PN = new TColgp_HArray1OfPnt(1,5);
657 Handle(TColgp_HArray1OfVec) VN = new TColgp_HArray1OfVec(1,4);
658 vector<const SMDS_MeshNode*> FNodes(5);
662 for (TopExp_Explorer exp(aShape,TopAbs_FACE);exp.More();exp.Next())
664 const TopoDS_Shape& aShapeFace = exp.Current();
666 aSubMeshDSFace = aProxyMesh->GetSubMesh( aShapeFace );
668 aSubMeshDSFace = meshDS->MeshElements( aShapeFace );
670 vector<const SMDS_MeshElement*> trias, quads;
671 bool hasNewTrias = false;
673 if ( aSubMeshDSFace )
676 if ( helper.NbAncestors( aShapeFace, aMesh, aShape.ShapeType() ) > 1 )
677 isRev = SMESH_Algo::IsReversedSubMesh( TopoDS::Face(aShapeFace), meshDS );
679 SMDS_ElemIteratorPtr iteratorElem = aSubMeshDSFace->GetElements();
680 while ( iteratorElem->more() ) // loop on elements on a geometrical face
682 const SMDS_MeshElement* face = iteratorElem->next();
683 // preparation step to get face info
684 int stat = Preparation(face, PN, VN, FNodes, PC, VNorm);
689 trias.push_back( face );
695 // add triangles to result map
696 SMDS_MeshFace* NewFace;
698 NewFace = meshDS->AddFace( FNodes[0], FNodes[1], FNodes[2] );
700 NewFace = meshDS->AddFace( FNodes[0], FNodes[2], FNodes[1] );
701 storeTmpElement( NewFace );
702 trias.push_back ( NewFace );
703 quads.push_back( face );
710 if(!isRev) VNorm.Reverse();
711 double xc = 0., yc = 0., zc = 0.;
716 Pbest = FindBestPoint(PN->Value(i), PN->Value(i+1), PC, VN->Value(i).Reversed());
718 Pbest = FindBestPoint(PN->Value(i), PN->Value(i+1), PC, VN->Value(i));
723 gp_Pnt PCbest(xc/4., yc/4., zc/4.);
726 double height = PCbest.Distance(PC);
728 // create new PCbest using a bit shift along VNorm
729 PCbest = PC.XYZ() + VNorm.XYZ() * 0.001;
732 // check possible intersection with other faces
734 bool check = CheckIntersection(PCbest, PC, Pint, aMesh, aShape, face);
736 //cout<<"--PC("<<PC.X()<<","<<PC.Y()<<","<<PC.Z()<<")"<<endl;
737 //cout<<" PCbest("<<PCbest.X()<<","<<PCbest.Y()<<","<<PCbest.Z()<<")"<<endl;
738 double dist = PC.Distance(Pint)/3.;
739 gp_Dir aDir(gp_Vec(PC,PCbest));
740 PCbest = PC.XYZ() + aDir.XYZ() * dist;
743 gp_Vec VB(PC,PCbest);
744 gp_Pnt PCbestTmp = PC.XYZ() + VB.XYZ() * 3.0;
745 check = CheckIntersection(PCbestTmp, PC, Pint, aMesh, aShape, face);
747 double dist = PC.Distance(Pint)/3.;
749 gp_Dir aDir(gp_Vec(PC,PCbest));
750 PCbest = PC.XYZ() + aDir.XYZ() * dist;
755 // create node for PCbest
756 SMDS_MeshNode* NewNode = helper.AddNode( PCbest.X(), PCbest.Y(), PCbest.Z() );
758 // add triangles to result map
761 trias.push_back ( meshDS->AddFace( NewNode, FNodes[i], FNodes[i+1] ));
762 storeTmpElement( trias.back() );
765 if ( isRev ) swap( FNodes[1], FNodes[3]);
766 SMDS_MeshVolume* aPyram =
767 helper.AddVolume( FNodes[0], FNodes[1], FNodes[2], FNodes[3], NewNode );
768 myPyramids.push_back(aPyram);
770 quads.push_back( face );
777 } // end loop on elements on a face submesh
779 bool sourceSubMeshIsProxy = false;
782 // move proxy sub-mesh from other proxy mesh to this
783 sourceSubMeshIsProxy = takeProxySubMesh( aShapeFace, aProxyMesh );
784 // move also tmp elements added in mesh
785 takeTmpElemsInMesh( aProxyMesh );
789 SMESH_ProxyMesh::SubMesh* prxSubMesh = getProxySubMesh( aShapeFace );
790 prxSubMesh->ChangeElements( trias.begin(), trias.end() );
792 // delete tmp quadrangles removed from aProxyMesh
793 if ( sourceSubMeshIsProxy )
795 for ( unsigned i = 0; i < quads.size(); ++i )
796 removeTmpElement( quads[i] );
798 delete myElemSearcher;
800 SMESH_MeshEditor(&aMesh).GetElementSearcher( aProxyMesh->GetFaces(aShape));
804 } // end for(TopExp_Explorer exp(aShape,TopAbs_FACE);exp.More();exp.Next()) {
806 return Compute2ndPart(aMesh, myPyramids);
809 //================================================================================
811 * \brief Computes pyramids in mesh with no shape
813 //================================================================================
815 bool StdMeshers_QuadToTriaAdaptor::Compute(SMESH_Mesh& aMesh)
817 SMESH_ProxyMesh::setMesh( aMesh );
818 SMESH_ProxyMesh::_allowedTypes.push_back( SMDSEntity_Triangle );
819 SMESH_ProxyMesh::_allowedTypes.push_back( SMDSEntity_Quad_Triangle );
820 if ( aMesh.NbQuadrangles() < 1 )
823 vector<const SMDS_MeshElement*> myPyramids;
824 SMESH_MesherHelper helper(aMesh);
825 helper.IsQuadraticSubMesh(aMesh.GetShapeToMesh());
826 helper.SetElementsOnShape( true );
828 if ( !myElemSearcher )
829 myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher();
830 SMESH_ElementSearcher* searcher = const_cast<SMESH_ElementSearcher*>(myElemSearcher);
832 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
833 SMESH_ProxyMesh::SubMesh* prxSubMesh = getProxySubMesh();
835 SMDS_FaceIteratorPtr fIt = meshDS->facesIterator(/*idInceasingOrder=*/true);
838 const SMDS_MeshElement* face = fIt->next();
839 if ( !face ) continue;
840 // retrieve needed information about a face
841 Handle(TColgp_HArray1OfPnt) PN = new TColgp_HArray1OfPnt(1,5);
842 Handle(TColgp_HArray1OfVec) VN = new TColgp_HArray1OfVec(1,4);
843 vector<const SMDS_MeshNode*> FNodes(5);
846 const SMDS_MeshElement* volumes[2];
847 int what = Preparation(face, PN, VN, FNodes, PC, VNorm, volumes);
848 if ( what == NOT_QUAD )
850 if ( volumes[0] && volumes[1] )
851 continue; // face is shared by two volumes - no space for a pyramid
853 if ( what == DEGEN_QUAD )
856 // add a triangle to the proxy mesh
857 SMDS_MeshFace* NewFace;
860 double tmp = PN->Value(1).Distance(PN->Value(2)) + PN->Value(2).Distance(PN->Value(3));
861 // far points in VNorm direction
862 gp_Pnt Ptmp1 = PC.XYZ() + VNorm.XYZ() * tmp * 1.e6;
863 gp_Pnt Ptmp2 = PC.XYZ() - VNorm.XYZ() * tmp * 1.e6;
864 // check intersection for Ptmp1 and Ptmp2
868 double dist1 = RealLast();
869 double dist2 = RealLast();
872 gp_Ax1 line( PC, VNorm );
873 vector< const SMDS_MeshElement* > suspectElems;
874 searcher->GetElementsNearLine( line, SMDSAbs_Face, suspectElems);
876 for ( int iF = 0; iF < suspectElems.size(); ++iF ) {
877 const SMDS_MeshElement* F = suspectElems[iF];
878 if(F==face) continue;
879 Handle(TColgp_HSequenceOfPnt) aContour = new TColgp_HSequenceOfPnt;
880 for ( int i = 0; i < 4; ++i )
881 aContour->Append( SMESH_TNodeXYZ( F->GetNode(i) ));
883 if( !volumes[0] && HasIntersection(Ptmp1, PC, PPP, aContour) ) {
885 double tmp = PC.Distance(PPP);
891 if( !volumes[1] && HasIntersection(Ptmp2, PC, PPP, aContour) ) {
893 double tmp = PC.Distance(PPP);
901 if( IsOK1 && !IsOK2 ) {
902 // using existed direction
904 else if( !IsOK1 && IsOK2 ) {
905 // using opposite direction
908 else { // IsOK1 && IsOK2
909 double tmp1 = PC.Distance(Pres1);
910 double tmp2 = PC.Distance(Pres2);
912 // using existed direction
915 // using opposite direction
920 NewFace = meshDS->AddFace( FNodes[0], FNodes[1], FNodes[2] );
922 NewFace = meshDS->AddFace( FNodes[0], FNodes[2], FNodes[1] );
923 storeTmpElement( NewFace );
924 prxSubMesh->AddElement( NewFace );
928 // Case of non-degenerated quadrangle
932 gp_XYZ PCbest(0., 0., 0.); // pyramid peak
935 gp_Pnt Pbest = FindBestPoint(PN->Value(i), PN->Value(i+1), PC, VN->Value(i));
936 PCbest += Pbest.XYZ();
940 double height = PC.Distance(PCbest); // pyramid height to precise
942 // create new PCbest using a bit shift along VNorm
943 PCbest = PC.XYZ() + VNorm.XYZ() * 0.001;
944 height = PC.Distance(PCbest);
947 // Restrict pyramid height by intersection with other faces
948 gp_Vec tmpDir(PC,PCbest); tmpDir.Normalize();
949 double tmp = PN->Value(1).Distance(PN->Value(3)) + PN->Value(2).Distance(PN->Value(4));
950 // far points: in (PC, PCbest) direction and vice-versa
951 gp_Pnt farPnt[2] = { PC.XYZ() + tmpDir.XYZ() * tmp * 1.e6,
952 PC.XYZ() - tmpDir.XYZ() * tmp * 1.e6 };
953 // check intersection for farPnt1 and farPnt2
954 bool intersected[2] = { false, false };
955 double dist [2] = { RealLast(), RealLast() };
958 gp_Ax1 line( PC, tmpDir );
959 vector< const SMDS_MeshElement* > suspectElems;
960 searcher->GetElementsNearLine( line, SMDSAbs_Face, suspectElems);
962 for ( int iF = 0; iF < suspectElems.size(); ++iF )
964 const SMDS_MeshElement* F = suspectElems[iF];
965 if(F==face) continue;
966 Handle(TColgp_HSequenceOfPnt) aContour = new TColgp_HSequenceOfPnt;
967 int nbN = F->NbNodes() / ( F->IsQuadratic() ? 2 : 1 );
968 for ( i = 0; i < nbN; ++i )
969 aContour->Append( SMESH_TNodeXYZ( F->GetNode(i) ));
971 for ( int isRev = 0; isRev < 2; ++isRev )
973 if( !volumes[isRev] && HasIntersection(farPnt[isRev], PC, intP, aContour) ) {
974 intersected[isRev] = true;
975 double d = PC.Distance( intP );
976 if( d < dist[isRev] )
978 intPnt[isRev] = intP;
985 // Create one or two pyramids
987 for ( int isRev = 0; isRev < 2; ++isRev )
989 if( !intersected[isRev] ) continue;
990 double pyramidH = Min( height, PC.Distance(intPnt[isRev])/3.);
991 PCbest = PC.XYZ() + tmpDir.XYZ() * (isRev ? -pyramidH : pyramidH);
993 // create node for PCbest
994 SMDS_MeshNode* NewNode = helper.AddNode( PCbest.X(), PCbest.Y(), PCbest.Z() );
996 // add triangles to result map
998 SMDS_MeshFace* NewFace;
1000 NewFace = meshDS->AddFace( NewNode, FNodes[i], FNodes[i+1] );
1002 NewFace = meshDS->AddFace( NewNode, FNodes[i+1], FNodes[i] );
1003 storeTmpElement( NewFace );
1004 prxSubMesh->AddElement( NewFace );
1007 SMDS_MeshVolume* aPyram;
1009 aPyram = helper.AddVolume( FNodes[0], FNodes[1], FNodes[2], FNodes[3], NewNode );
1011 aPyram = helper.AddVolume( FNodes[0], FNodes[3], FNodes[2], FNodes[1], NewNode );
1012 myPyramids.push_back(aPyram);
1014 } // end loop on all faces
1016 return Compute2ndPart(aMesh, myPyramids);
1019 //================================================================================
1021 * \brief Update created pyramids and faces to avoid their intersection
1023 //================================================================================
1025 bool StdMeshers_QuadToTriaAdaptor::Compute2ndPart(SMESH_Mesh& aMesh,
1026 const vector<const SMDS_MeshElement*>& myPyramids)
1028 if(myPyramids.empty())
1031 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
1032 int i, j, k, myShapeID = myPyramids[0]->GetNode(4)->getshapeId();
1034 if ( !myElemSearcher )
1035 myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher();
1036 SMESH_ElementSearcher* searcher = const_cast<SMESH_ElementSearcher*>(myElemSearcher);
1038 set<const SMDS_MeshNode*> nodesToMove;
1040 // check adjacent pyramids
1042 for ( i = 0; i < myPyramids.size(); ++i )
1044 const SMDS_MeshElement* PrmI = myPyramids[i];
1045 MergeAdjacent( PrmI, nodesToMove );
1048 // iterate on all pyramids
1049 for ( i = 0; i < myPyramids.size(); ++i )
1051 const SMDS_MeshElement* PrmI = myPyramids[i];
1053 // compare PrmI with all the rest pyramids
1055 // collect adjacent pyramids and nodes coordinates of PrmI
1056 set<const SMDS_MeshElement*> checkedPyrams;
1057 vector<gp_Pnt> PsI(5);
1058 for(k=0; k<5; k++) // loop on 4 base nodes of PrmI
1060 const SMDS_MeshNode* n = PrmI->GetNode(k);
1061 PsI[k] = SMESH_TNodeXYZ( n );
1062 SMDS_ElemIteratorPtr vIt = n->GetInverseElementIterator( SMDSAbs_Volume );
1063 while ( vIt->more() )
1065 const SMDS_MeshElement* PrmJ = vIt->next();
1066 if ( SMESH_Algo::GetCommonNodes( PrmI, PrmJ ).size() > 1 )
1067 checkedPyrams.insert( PrmJ );
1071 // check intersection with distant pyramids
1072 for(k=0; k<4; k++) // loop on 4 base nodes of PrmI
1074 gp_Vec Vtmp(PsI[k],PsI[4]);
1075 gp_Ax1 line( PsI[k], Vtmp );
1076 vector< const SMDS_MeshElement* > suspectPyrams;
1077 searcher->GetElementsNearLine( line, SMDSAbs_Volume, suspectPyrams);
1079 for ( j = 0; j < suspectPyrams.size(); ++j )
1081 const SMDS_MeshElement* PrmJ = suspectPyrams[j];
1082 if ( PrmJ == PrmI || PrmJ->NbCornerNodes() != 5 )
1084 if ( myShapeID != PrmJ->GetNode(4)->getshapeId())
1085 continue; // pyramid from other SOLID
1086 if ( PrmI->GetNode(4) == PrmJ->GetNode(4) )
1087 continue; // pyramids PrmI and PrmJ already merged
1088 if ( !checkedPyrams.insert( PrmJ ).second )
1089 continue; // already checked
1091 TXyzIterator xyzIt( PrmJ->nodesIterator() );
1092 vector<gp_Pnt> PsJ( xyzIt, TXyzIterator() );
1096 for(k=0; k<4 && !hasInt; k++) {
1097 gp_Vec Vtmp(PsI[k],PsI[4]);
1098 gp_Pnt Pshift = PsI[k].XYZ() + Vtmp.XYZ() * 0.01; // base node moved a bit to apex
1100 ( HasIntersection3( Pshift, PsI[4], Pint, PsJ[0], PsJ[1], PsJ[4]) ||
1101 HasIntersection3( Pshift, PsI[4], Pint, PsJ[1], PsJ[2], PsJ[4]) ||
1102 HasIntersection3( Pshift, PsI[4], Pint, PsJ[2], PsJ[3], PsJ[4]) ||
1103 HasIntersection3( Pshift, PsI[4], Pint, PsJ[3], PsJ[0], PsJ[4]) );
1105 for(k=0; k<4 && !hasInt; k++) {
1106 gp_Vec Vtmp(PsJ[k],PsJ[4]);
1107 gp_Pnt Pshift = PsJ[k].XYZ() + Vtmp.XYZ() * 0.01;
1109 ( HasIntersection3( Pshift, PsJ[4], Pint, PsI[0], PsI[1], PsI[4]) ||
1110 HasIntersection3( Pshift, PsJ[4], Pint, PsI[1], PsI[2], PsI[4]) ||
1111 HasIntersection3( Pshift, PsJ[4], Pint, PsI[2], PsI[3], PsI[4]) ||
1112 HasIntersection3( Pshift, PsJ[4], Pint, PsI[3], PsI[0], PsI[4]) );
1117 // count common nodes of base faces of two pyramids
1120 nbc += int ( PrmI->GetNodeIndex( PrmJ->GetNode(k) ) >= 0 );
1123 continue; // pyrams have a common base face
1127 // Merge the two pyramids and others already merged with them
1128 MergePiramids( PrmI, PrmJ, nodesToMove );
1132 // decrease height of pyramids
1133 gp_XYZ PCi(0,0,0), PCj(0,0,0);
1134 for(k=0; k<4; k++) {
1135 PCi += PsI[k].XYZ();
1136 PCj += PsJ[k].XYZ();
1139 gp_Vec VN1(PCi,PsI[4]);
1140 gp_Vec VN2(PCj,PsJ[4]);
1141 gp_Vec VI1(PCi,Pint);
1142 gp_Vec VI2(PCj,Pint);
1143 double ang1 = fabs(VN1.Angle(VI1));
1144 double ang2 = fabs(VN2.Angle(VI2));
1145 double coef1 = 0.5 - (( ang1<PI/3 ) ? cos(ang1)*0.25 : 0 );
1146 double coef2 = 0.5 - (( ang2<PI/3 ) ? cos(ang2)*0.25 : 0 ); // cos(ang2) ?
1147 // double coef2 = 0.5;
1149 // coef2 -= cos(ang1)*0.25;
1153 SMDS_MeshNode* aNode1 = const_cast<SMDS_MeshNode*>(PrmI->GetNode(4));
1154 aNode1->setXYZ( PCi.X()+VN1.X(), PCi.Y()+VN1.Y(), PCi.Z()+VN1.Z() );
1155 SMDS_MeshNode* aNode2 = const_cast<SMDS_MeshNode*>(PrmJ->GetNode(4));
1156 aNode2->setXYZ( PCj.X()+VN2.X(), PCj.Y()+VN2.Y(), PCj.Z()+VN2.Z() );
1157 nodesToMove.insert( aNode1 );
1158 nodesToMove.insert( aNode2 );
1160 // fix intersections that could appear after apex movement
1161 MergeAdjacent( PrmI, nodesToMove );
1162 MergeAdjacent( PrmJ, nodesToMove );
1165 } // loop on suspectPyrams
1166 } // loop on 4 base nodes of PrmI
1168 } // loop on all pyramids
1170 if( !nodesToMove.empty() && !meshDS->IsEmbeddedMode() )
1172 set<const SMDS_MeshNode*>::iterator n = nodesToMove.begin();
1173 for ( ; n != nodesToMove.end(); ++n )
1174 meshDS->MoveNode( *n, (*n)->X(), (*n)->Y(), (*n)->Z() );
1177 // erase removed triangles from the proxy mesh
1178 if ( !myRemovedTrias.empty() )
1180 for ( int i = 0; i <= meshDS->MaxShapeIndex(); ++i )
1181 if ( SMESH_ProxyMesh::SubMesh* sm = findProxySubMesh(i))
1183 vector<const SMDS_MeshElement *> faces;
1184 faces.reserve( sm->NbElements() );
1185 SMDS_ElemIteratorPtr fIt = sm->GetElements();
1186 while ( fIt->more() )
1188 const SMDS_MeshElement* tria = fIt->next();
1189 set<const SMDS_MeshElement*>::iterator rmTria = myRemovedTrias.find( tria );
1190 if ( rmTria != myRemovedTrias.end() )
1191 myRemovedTrias.erase( rmTria );
1193 faces.push_back( tria );
1195 sm->ChangeElements( faces.begin(), faces.end() );
1201 delete myElemSearcher;
1207 //================================================================================
1209 * \brief Return list of created triangles for given face
1211 //================================================================================
1213 // const list<const SMDS_MeshFace* >* StdMeshers_QuadToTriaAdaptor::GetTriangles (const SMDS_MeshElement* aQuad)
1215 // TQuad2Trias::iterator it = myResMap.find(aQuad);
1216 // return ( it != myResMap.end() ? & it->second : 0 );