1 // Copyright (C) 2007-2012 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 // File : StdMeshers_QuadToTriaAdaptor.cxx
22 // Created : Wen May 07 16:37:07 2008
23 // Author : Sergey KUUL (skl)
25 #include "StdMeshers_QuadToTriaAdaptor.hxx"
27 #include "SMDS_SetIterator.hxx"
29 #include "SMESH_Algo.hxx"
30 #include "SMESH_MesherHelper.hxx"
32 #include <IntAna_IntConicQuad.hxx>
33 #include <IntAna_Quadric.hxx>
34 #include <TColgp_HArray1OfPnt.hxx>
35 #include <TColgp_HArray1OfVec.hxx>
36 #include <TColgp_HSequenceOfPnt.hxx>
37 #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. * M_PI / 180. );
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
180 //================================================================================
182 * \brief Move medium nodes of merged quadratic pyramids
184 //================================================================================
186 void UpdateQuadraticPyramids(const set<const SMDS_MeshNode*>& commonApex,
187 SMESHDS_Mesh* meshDS)
189 typedef SMDS_StdIterator< const SMDS_MeshElement*, SMDS_ElemIteratorPtr > TStdElemIterator;
190 TStdElemIterator itEnd;
192 // shift of node index to get medium nodes between the 4 base nodes and the apex
193 const int base2MediumShift = 9;
195 set<const SMDS_MeshNode*>::const_iterator nIt = commonApex.begin();
196 for ( ; nIt != commonApex.end(); ++nIt )
198 SMESH_TNodeXYZ apex( *nIt );
200 vector< const SMDS_MeshElement* > pyrams // pyramids sharing the apex node
201 ( TStdElemIterator( apex._node->GetInverseElementIterator( SMDSAbs_Volume )), itEnd );
203 // Select medium nodes to keep and medium nodes to remove
205 typedef map < const SMDS_MeshNode*, const SMDS_MeshNode*, TIDCompare > TN2NMap;
206 TN2NMap base2medium; // to keep
207 vector< const SMDS_MeshNode* > nodesToRemove;
209 for ( unsigned i = 0; i < pyrams.size(); ++i )
210 for ( int baseIndex = 0; baseIndex < PYRAM_APEX; ++baseIndex )
212 SMESH_TNodeXYZ base = pyrams[i]->GetNode( baseIndex );
213 const SMDS_MeshNode* medium = pyrams[i]->GetNode( baseIndex + base2MediumShift );
214 TN2NMap::iterator b2m = base2medium.insert( make_pair( base._node, medium )).first;
215 if ( b2m->second != medium )
217 nodesToRemove.push_back( medium );
221 // move the kept medium node
222 gp_XYZ newXYZ = 0.5 * ( apex + base );
223 meshDS->MoveNode( medium, newXYZ.X(), newXYZ.Y(), newXYZ.Z() );
227 // Within pyramids, replace nodes to remove by nodes to keep
229 for ( unsigned i = 0; i < pyrams.size(); ++i )
231 vector< const SMDS_MeshNode* > nodes( pyrams[i]->begin_nodes(),
232 pyrams[i]->end_nodes() );
233 for ( int baseIndex = 0; baseIndex < PYRAM_APEX; ++baseIndex )
235 const SMDS_MeshNode* base = pyrams[i]->GetNode( baseIndex );
236 nodes[ baseIndex + base2MediumShift ] = base2medium[ base ];
238 meshDS->ChangeElementNodes( pyrams[i], &nodes[0], nodes.size());
241 // Remove the replaced nodes
243 if ( !nodesToRemove.empty() )
245 SMESHDS_SubMesh * sm = meshDS->MeshElements( nodesToRemove[0]->getshapeId() );
246 for ( unsigned i = 0; i < nodesToRemove.size(); ++i )
247 meshDS->RemoveFreeNode( nodesToRemove[i], sm, /*fromGroups=*/false);
254 //================================================================================
256 * \brief Merge the two pyramids (i.e. fuse their apex) and others already merged with them
258 //================================================================================
260 void StdMeshers_QuadToTriaAdaptor::MergePiramids( const SMDS_MeshElement* PrmI,
261 const SMDS_MeshElement* PrmJ,
262 set<const SMDS_MeshNode*> & nodesToMove)
264 const SMDS_MeshNode* Nrem = PrmJ->GetNode(4); // node to remove
265 //int nbJ = Nrem->NbInverseElements( SMDSAbs_Volume );
266 SMESH_TNodeXYZ Pj( Nrem );
268 // an apex node to make common to all merged pyramids
269 SMDS_MeshNode* CommonNode = const_cast<SMDS_MeshNode*>(PrmI->GetNode(4));
270 if ( CommonNode == Nrem ) return; // already merged
271 //int nbI = CommonNode->NbInverseElements( SMDSAbs_Volume );
272 SMESH_TNodeXYZ Pi( CommonNode );
273 gp_XYZ Pnew = /*( nbI*Pi + nbJ*Pj ) / (nbI+nbJ);*/ 0.5 * ( Pi + Pj );
274 CommonNode->setXYZ( Pnew.X(), Pnew.Y(), Pnew.Z() );
276 nodesToMove.insert( CommonNode );
277 nodesToMove.erase ( Nrem );
279 typedef SMDS_StdIterator< const SMDS_MeshElement*, SMDS_ElemIteratorPtr > TStdElemIterator;
280 TStdElemIterator itEnd;
282 // find and remove coincided faces of merged pyramids
283 vector< const SMDS_MeshElement* > inverseElems
284 // copy inverse elements to avoid iteration on changing container
285 ( TStdElemIterator( CommonNode->GetInverseElementIterator(SMDSAbs_Face)), itEnd);
286 for ( unsigned i = 0; i < inverseElems.size(); ++i )
288 const SMDS_MeshElement* FI = inverseElems[i];
289 const SMDS_MeshElement* FJEqual = 0;
290 SMDS_ElemIteratorPtr triItJ = Nrem->GetInverseElementIterator(SMDSAbs_Face);
291 while ( !FJEqual && triItJ->more() )
293 const SMDS_MeshElement* FJ = triItJ->next();
294 if ( EqualTriangles( FJ, FI ))
299 removeTmpElement( FI );
300 removeTmpElement( FJEqual );
301 myRemovedTrias.insert( FI );
302 myRemovedTrias.insert( FJEqual );
306 // set the common apex node to pyramids and triangles merged with J
307 inverseElems.assign( TStdElemIterator( Nrem->GetInverseElementIterator()), itEnd );
308 for ( unsigned i = 0; i < inverseElems.size(); ++i )
310 const SMDS_MeshElement* elem = inverseElems[i];
311 vector< const SMDS_MeshNode* > nodes( elem->begin_nodes(), elem->end_nodes() );
312 nodes[ elem->GetType() == SMDSAbs_Volume ? PYRAM_APEX : TRIA_APEX ] = CommonNode;
313 GetMeshDS()->ChangeElementNodes( elem, &nodes[0], nodes.size());
315 ASSERT( Nrem->NbInverseElements() == 0 );
316 GetMeshDS()->RemoveFreeNode( Nrem,
317 GetMeshDS()->MeshElements( Nrem->getshapeId()),
318 /*fromGroups=*/false);
321 //================================================================================
323 * \brief Merges adjacent pyramids
325 //================================================================================
327 void StdMeshers_QuadToTriaAdaptor::MergeAdjacent(const SMDS_MeshElement* PrmI,
328 set<const SMDS_MeshNode*>& nodesToMove)
330 TIDSortedElemSet adjacentPyrams;
331 bool mergedPyrams = false;
332 for(int k=0; k<4; k++) // loop on 4 base nodes of PrmI
334 const SMDS_MeshNode* n = PrmI->GetNode(k);
335 SMDS_ElemIteratorPtr vIt = n->GetInverseElementIterator( SMDSAbs_Volume );
336 while ( vIt->more() )
338 const SMDS_MeshElement* PrmJ = vIt->next();
339 if ( PrmJ->NbCornerNodes() != 5 || !adjacentPyrams.insert( PrmJ ).second )
341 if ( PrmI != PrmJ && TooCloseAdjacent( PrmI, PrmJ, GetMesh()->HasShapeToMesh() ))
343 MergePiramids( PrmI, PrmJ, nodesToMove );
345 // container of inverse elements can change
346 vIt = n->GetInverseElementIterator( SMDSAbs_Volume );
352 TIDSortedElemSet::iterator prm;
353 for (prm = adjacentPyrams.begin(); prm != adjacentPyrams.end(); ++prm)
354 MergeAdjacent( *prm, nodesToMove );
358 //================================================================================
362 //================================================================================
364 StdMeshers_QuadToTriaAdaptor::StdMeshers_QuadToTriaAdaptor():
369 //================================================================================
373 //================================================================================
375 StdMeshers_QuadToTriaAdaptor::~StdMeshers_QuadToTriaAdaptor()
377 // temporary faces are deleted by ~SMESH_ProxyMesh()
378 if ( myElemSearcher ) delete myElemSearcher;
382 //=======================================================================
383 //function : FindBestPoint
384 //purpose : Return a point P laying on the line (PC,V) so that triangle
385 // (P, P1, P2) to be equilateral as much as possible
386 // V - normal to (P1,P2,PC)
387 //=======================================================================
389 static gp_Pnt FindBestPoint(const gp_Pnt& P1, const gp_Pnt& P2,
390 const gp_Pnt& PC, const gp_Vec& V)
393 const double a = P1.Distance(P2);
394 const double b = P1.Distance(PC);
395 const double c = P2.Distance(PC);
399 // find shift along V in order a to became equal to (b+c)/2
400 const double Vsize = V.Magnitude();
401 if ( fabs( Vsize ) > std::numeric_limits<double>::min() )
403 const double shift = sqrt( a*a + (b*b-c*c)*(b*b-c*c)/16/a/a - (b*b+c*c)/2 );
404 Pbest.ChangeCoord() += shift * V.XYZ() / Vsize;
410 //=======================================================================
411 //function : HasIntersection3
412 //purpose : Auxilare for HasIntersection()
413 // find intersection point between triangle (P1,P2,P3)
414 // and segment [PC,P]
415 //=======================================================================
417 static bool HasIntersection3(const gp_Pnt& P, const gp_Pnt& PC, gp_Pnt& Pint,
418 const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3)
420 //cout<<"HasIntersection3"<<endl;
421 //cout<<" PC("<<PC.X()<<","<<PC.Y()<<","<<PC.Z()<<")"<<endl;
422 //cout<<" P("<<P.X()<<","<<P.Y()<<","<<P.Z()<<")"<<endl;
423 //cout<<" P1("<<P1.X()<<","<<P1.Y()<<","<<P1.Z()<<")"<<endl;
424 //cout<<" P2("<<P2.X()<<","<<P2.Y()<<","<<P2.Z()<<")"<<endl;
425 //cout<<" P3("<<P3.X()<<","<<P3.Y()<<","<<P3.Z()<<")"<<endl;
428 IntAna_Quadric IAQ(gp_Pln(P1,VP1.Crossed(VP2)));
429 IntAna_IntConicQuad IAICQ(gp_Lin(PC,gp_Dir(gp_Vec(PC,P))),IAQ);
431 if( IAICQ.IsInQuadric() )
433 if( IAICQ.NbPoints() == 1 ) {
434 gp_Pnt PIn = IAICQ.Point(1);
435 const double preci = 1.e-10 * P.Distance(PC);
436 // check if this point is internal for segment [PC,P]
438 ( (PC.X()-PIn.X())*(P.X()-PIn.X()) > preci ) ||
439 ( (PC.Y()-PIn.Y())*(P.Y()-PIn.Y()) > preci ) ||
440 ( (PC.Z()-PIn.Z())*(P.Z()-PIn.Z()) > preci );
444 // check if this point is internal for triangle (P1,P2,P3)
448 if( V1.Magnitude()<preci ||
449 V2.Magnitude()<preci ||
450 V3.Magnitude()<preci ) {
454 const double angularTol = 1e-6;
455 gp_Vec VC1 = V1.Crossed(V2);
456 gp_Vec VC2 = V2.Crossed(V3);
457 gp_Vec VC3 = V3.Crossed(V1);
458 if(VC1.Magnitude()<gp::Resolution()) {
459 if(VC2.IsOpposite(VC3,angularTol)) {
463 else if(VC2.Magnitude()<gp::Resolution()) {
464 if(VC1.IsOpposite(VC3,angularTol)) {
468 else if(VC3.Magnitude()<gp::Resolution()) {
469 if(VC1.IsOpposite(VC2,angularTol)) {
474 if( VC1.IsOpposite(VC2,angularTol) || VC1.IsOpposite(VC3,angularTol) ||
475 VC2.IsOpposite(VC3,angularTol) ) {
487 //=======================================================================
488 //function : HasIntersection
489 //purpose : Auxilare for CheckIntersection()
490 //=======================================================================
492 static bool HasIntersection(const gp_Pnt& P, const gp_Pnt& PC, gp_Pnt& Pint,
493 Handle(TColgp_HSequenceOfPnt)& aContour)
495 if(aContour->Length()==3) {
496 return HasIntersection3( P, PC, Pint, aContour->Value(1),
497 aContour->Value(2), aContour->Value(3) );
501 if( (aContour->Value(1).Distance(aContour->Value(2)) > 1.e-6) &&
502 (aContour->Value(1).Distance(aContour->Value(3)) > 1.e-6) &&
503 (aContour->Value(2).Distance(aContour->Value(3)) > 1.e-6) ) {
504 check = HasIntersection3( P, PC, Pint, aContour->Value(1),
505 aContour->Value(2), aContour->Value(3) );
507 if(check) return true;
508 if( (aContour->Value(1).Distance(aContour->Value(4)) > 1.e-6) &&
509 (aContour->Value(1).Distance(aContour->Value(3)) > 1.e-6) &&
510 (aContour->Value(4).Distance(aContour->Value(3)) > 1.e-6) ) {
511 check = HasIntersection3( P, PC, Pint, aContour->Value(1),
512 aContour->Value(3), aContour->Value(4) );
514 if(check) return true;
520 //================================================================================
522 * \brief Checks if a line segment (P,PC) intersects any mesh face.
523 * \param P - first segment end
524 * \param PC - second segment end (it is a gravity center of quadrangle)
525 * \param Pint - (out) intersection point
526 * \param aMesh - mesh
527 * \param aShape - shape to check faces on
528 * \param NotCheckedFace - mesh face not to check
529 * \retval bool - true if there is an intersection
531 //================================================================================
533 bool StdMeshers_QuadToTriaAdaptor::CheckIntersection (const gp_Pnt& P,
537 const TopoDS_Shape& aShape,
538 const SMDS_MeshElement* NotCheckedFace)
540 if ( !myElemSearcher )
541 myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher();
542 SMESH_ElementSearcher* searcher = const_cast<SMESH_ElementSearcher*>(myElemSearcher);
544 //SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
545 //cout<<" CheckIntersection: meshDS->NbFaces() = "<<meshDS->NbFaces()<<endl;
547 double dist = RealLast(); // find intersection closest to the segment
550 gp_Ax1 line( P, gp_Vec(P,PC));
551 vector< const SMDS_MeshElement* > suspectElems;
552 searcher->GetElementsNearLine( line, SMDSAbs_Face, suspectElems);
554 for ( int i = 0; i < suspectElems.size(); ++i )
556 const SMDS_MeshElement* face = suspectElems[i];
557 if ( face == NotCheckedFace ) continue;
558 Handle(TColgp_HSequenceOfPnt) aContour = new TColgp_HSequenceOfPnt;
559 for ( int i = 0; i < face->NbCornerNodes(); ++i )
560 aContour->Append( SMESH_TNodeXYZ( face->GetNode(i) ));
561 if( HasIntersection(P, PC, Pres, aContour) ) {
563 double tmp = PC.Distance(Pres);
573 //================================================================================
575 * \brief Prepare data for the given face
576 * \param PN - coordinates of face nodes
577 * \param VN - cross products of vectors (PC-PN(i)) ^ (PC-PN(i+1))
578 * \param FNodes - face nodes
579 * \param PC - gravity center of nodes
580 * \param VNorm - face normal (sum of VN)
581 * \param volumes - two volumes sharing the given face, the first is in VNorm direction
582 * \retval int - 0 if given face is not quad,
583 * 1 if given face is quad,
584 * 2 if given face is degenerate quad (two nodes are coincided)
586 //================================================================================
588 int StdMeshers_QuadToTriaAdaptor::Preparation(const SMDS_MeshElement* face,
589 Handle(TColgp_HArray1OfPnt)& PN,
590 Handle(TColgp_HArray1OfVec)& VN,
591 vector<const SMDS_MeshNode*>& FNodes,
594 const SMDS_MeshElement** volumes)
596 if( face->NbCornerNodes() != 4 )
602 gp_XYZ xyzC(0., 0., 0.);
603 for ( i = 0; i < 4; ++i )
605 gp_XYZ p = SMESH_TNodeXYZ( FNodes[i] = face->GetNode(i) );
606 PN->SetValue( i+1, p );
617 if( PN->Value(i).Distance(PN->Value(j)) < 1.e-6 )
622 //int deg_num = IsDegenarate(PN);
626 //cout<<"find degeneration"<<endl;
628 gp_Pnt Pdeg = PN->Value(i);
630 list< const SMDS_MeshNode* >::iterator itdg = myDegNodes.begin();
631 const SMDS_MeshNode* DegNode = 0;
632 for(; itdg!=myDegNodes.end(); itdg++) {
633 const SMDS_MeshNode* N = (*itdg);
634 gp_Pnt Ptmp(N->X(),N->Y(),N->Z());
635 if(Pdeg.Distance(Ptmp)<1.e-6) {
637 //DegNode = const_cast<SMDS_MeshNode*>(N);
642 DegNode = FNodes[i-1];
643 myDegNodes.push_back(DegNode);
646 FNodes[i-1] = DegNode;
649 PN->SetValue(i,PN->Value(i+1));
650 FNodes[i-1] = FNodes[i];
655 PN->SetValue(nbp+1,PN->Value(1));
656 FNodes[nbp] = FNodes[0];
657 // find normal direction
658 gp_Vec V1(PC,PN->Value(nbp));
659 gp_Vec V2(PC,PN->Value(1));
660 VNorm = V1.Crossed(V2);
661 VN->SetValue(nbp,VNorm);
662 for(i=1; i<nbp; i++) {
663 V1 = gp_Vec(PC,PN->Value(i));
664 V2 = gp_Vec(PC,PN->Value(i+1));
665 gp_Vec Vtmp = V1.Crossed(V2);
666 VN->SetValue(i,Vtmp);
670 // find volumes sharing the face
673 volumes[0] = volumes[1] = 0;
674 SMDS_ElemIteratorPtr vIt = FNodes[0]->GetInverseElementIterator( SMDSAbs_Volume );
675 while ( vIt->more() )
677 const SMDS_MeshElement* vol = vIt->next();
678 bool volSharesAllNodes = true;
679 for ( int i = 1; i < face->NbNodes() && volSharesAllNodes; ++i )
680 volSharesAllNodes = ( vol->GetNodeIndex( FNodes[i] ) >= 0 );
681 if ( volSharesAllNodes )
682 volumes[ volumes[0] ? 1 : 0 ] = vol;
683 // we could additionally check that vol has all FNodes in its one face using SMDS_VolumeTool
685 // define volume position relating to the face normal
689 SMDS_ElemIteratorPtr nodeIt = volumes[0]->nodesIterator();
691 volGC = accumulate( TXyzIterator(nodeIt), TXyzIterator(), volGC ) / volumes[0]->NbNodes();
693 if ( VNorm * gp_Vec( PC, volGC ) < 0 )
694 swap( volumes[0], volumes[1] );
698 //cout<<" VNorm("<<VNorm.X()<<","<<VNorm.Y()<<","<<VNorm.Z()<<")"<<endl;
699 return hasdeg ? DEGEN_QUAD : QUAD;
703 //=======================================================================
706 //=======================================================================
708 bool StdMeshers_QuadToTriaAdaptor::Compute(SMESH_Mesh& aMesh,
709 const TopoDS_Shape& aShape,
710 SMESH_ProxyMesh* aProxyMesh)
712 SMESH_ProxyMesh::setMesh( aMesh );
714 if ( aShape.ShapeType() != TopAbs_SOLID &&
715 aShape.ShapeType() != TopAbs_SHELL )
720 vector<const SMDS_MeshElement*> myPyramids;
722 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
723 SMESH_MesherHelper helper(aMesh);
724 helper.IsQuadraticSubMesh(aShape);
725 helper.SetElementsOnShape( true );
727 if ( myElemSearcher ) delete myElemSearcher;
729 myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher( aProxyMesh->GetFaces(aShape));
731 myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher();
733 const SMESHDS_SubMesh * aSubMeshDSFace;
734 Handle(TColgp_HArray1OfPnt) PN = new TColgp_HArray1OfPnt(1,5);
735 Handle(TColgp_HArray1OfVec) VN = new TColgp_HArray1OfVec(1,4);
736 vector<const SMDS_MeshNode*> FNodes(5);
740 for (TopExp_Explorer exp(aShape,TopAbs_FACE);exp.More();exp.Next())
742 const TopoDS_Shape& aShapeFace = exp.Current();
744 aSubMeshDSFace = aProxyMesh->GetSubMesh( aShapeFace );
746 aSubMeshDSFace = meshDS->MeshElements( aShapeFace );
748 vector<const SMDS_MeshElement*> trias, quads;
749 bool hasNewTrias = false;
751 if ( aSubMeshDSFace )
754 if ( helper.NbAncestors( aShapeFace, aMesh, aShape.ShapeType() ) > 1 )
755 isRev = SMESH_Algo::IsReversedSubMesh( TopoDS::Face(aShapeFace), meshDS );
757 SMDS_ElemIteratorPtr iteratorElem = aSubMeshDSFace->GetElements();
758 while ( iteratorElem->more() ) // loop on elements on a geometrical face
760 const SMDS_MeshElement* face = iteratorElem->next();
761 // preparation step to get face info
762 int stat = Preparation(face, PN, VN, FNodes, PC, VNorm);
767 trias.push_back( face );
773 // add triangles to result map
774 SMDS_MeshFace* NewFace;
776 NewFace = meshDS->AddFace( FNodes[0], FNodes[1], FNodes[2] );
778 NewFace = meshDS->AddFace( FNodes[0], FNodes[2], FNodes[1] );
779 storeTmpElement( NewFace );
780 trias.push_back ( NewFace );
781 quads.push_back( face );
788 if(!isRev) VNorm.Reverse();
789 double xc = 0., yc = 0., zc = 0.;
794 Pbest = FindBestPoint(PN->Value(i), PN->Value(i+1), PC, VN->Value(i).Reversed());
796 Pbest = FindBestPoint(PN->Value(i), PN->Value(i+1), PC, VN->Value(i));
801 gp_Pnt PCbest(xc/4., yc/4., zc/4.);
804 double height = PCbest.Distance(PC);
806 // create new PCbest using a bit shift along VNorm
807 PCbest = PC.XYZ() + VNorm.XYZ() * 0.001;
810 // check possible intersection with other faces
812 bool check = CheckIntersection(PCbest, PC, Pint, aMesh, aShape, face);
814 //cout<<"--PC("<<PC.X()<<","<<PC.Y()<<","<<PC.Z()<<")"<<endl;
815 //cout<<" PCbest("<<PCbest.X()<<","<<PCbest.Y()<<","<<PCbest.Z()<<")"<<endl;
816 double dist = PC.Distance(Pint)/3.;
817 gp_Dir aDir(gp_Vec(PC,PCbest));
818 PCbest = PC.XYZ() + aDir.XYZ() * dist;
821 gp_Vec VB(PC,PCbest);
822 gp_Pnt PCbestTmp = PC.XYZ() + VB.XYZ() * 3.0;
823 check = CheckIntersection(PCbestTmp, PC, Pint, aMesh, aShape, face);
825 double dist = PC.Distance(Pint)/3.;
827 gp_Dir aDir(gp_Vec(PC,PCbest));
828 PCbest = PC.XYZ() + aDir.XYZ() * dist;
833 // create node for PCbest
834 SMDS_MeshNode* NewNode = helper.AddNode( PCbest.X(), PCbest.Y(), PCbest.Z() );
836 // add triangles to result map
839 trias.push_back ( meshDS->AddFace( NewNode, FNodes[i], FNodes[i+1] ));
840 storeTmpElement( trias.back() );
843 if ( isRev ) swap( FNodes[1], FNodes[3]);
844 SMDS_MeshVolume* aPyram =
845 helper.AddVolume( FNodes[0], FNodes[1], FNodes[2], FNodes[3], NewNode );
846 myPyramids.push_back(aPyram);
848 quads.push_back( face );
855 } // end loop on elements on a face submesh
857 bool sourceSubMeshIsProxy = false;
860 // move proxy sub-mesh from other proxy mesh to this
861 sourceSubMeshIsProxy = takeProxySubMesh( aShapeFace, aProxyMesh );
862 // move also tmp elements added in mesh
863 takeTmpElemsInMesh( aProxyMesh );
867 SMESH_ProxyMesh::SubMesh* prxSubMesh = getProxySubMesh( aShapeFace );
868 prxSubMesh->ChangeElements( trias.begin(), trias.end() );
870 // delete tmp quadrangles removed from aProxyMesh
871 if ( sourceSubMeshIsProxy )
873 for ( unsigned i = 0; i < quads.size(); ++i )
874 removeTmpElement( quads[i] );
876 delete myElemSearcher;
878 SMESH_MeshEditor(&aMesh).GetElementSearcher( aProxyMesh->GetFaces(aShape));
882 } // end for(TopExp_Explorer exp(aShape,TopAbs_FACE);exp.More();exp.Next()) {
884 return Compute2ndPart(aMesh, myPyramids);
887 //================================================================================
889 * \brief Computes pyramids in mesh with no shape
891 //================================================================================
893 bool StdMeshers_QuadToTriaAdaptor::Compute(SMESH_Mesh& aMesh)
895 SMESH_ProxyMesh::setMesh( aMesh );
896 SMESH_ProxyMesh::_allowedTypes.push_back( SMDSEntity_Triangle );
897 SMESH_ProxyMesh::_allowedTypes.push_back( SMDSEntity_Quad_Triangle );
898 if ( aMesh.NbQuadrangles() < 1 )
901 vector<const SMDS_MeshElement*> myPyramids;
902 SMESH_MesherHelper helper(aMesh);
903 helper.IsQuadraticSubMesh(aMesh.GetShapeToMesh());
904 helper.SetElementsOnShape( true );
906 if ( !myElemSearcher )
907 myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher();
908 SMESH_ElementSearcher* searcher = const_cast<SMESH_ElementSearcher*>(myElemSearcher);
910 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
911 SMESH_ProxyMesh::SubMesh* prxSubMesh = getProxySubMesh();
913 SMDS_FaceIteratorPtr fIt = meshDS->facesIterator(/*idInceasingOrder=*/true);
916 const SMDS_MeshElement* face = fIt->next();
917 if ( !face ) continue;
918 // retrieve needed information about a face
919 Handle(TColgp_HArray1OfPnt) PN = new TColgp_HArray1OfPnt(1,5);
920 Handle(TColgp_HArray1OfVec) VN = new TColgp_HArray1OfVec(1,4);
921 vector<const SMDS_MeshNode*> FNodes(5);
924 const SMDS_MeshElement* volumes[2];
925 int what = Preparation(face, PN, VN, FNodes, PC, VNorm, volumes);
926 if ( what == NOT_QUAD )
928 if ( volumes[0] && volumes[1] )
929 continue; // face is shared by two volumes - no space for a pyramid
931 if ( what == DEGEN_QUAD )
934 // add a triangle to the proxy mesh
935 SMDS_MeshFace* NewFace;
938 double tmp = PN->Value(1).Distance(PN->Value(2)) + PN->Value(2).Distance(PN->Value(3));
939 // far points in VNorm direction
940 gp_Pnt Ptmp1 = PC.XYZ() + VNorm.XYZ() * tmp * 1.e6;
941 gp_Pnt Ptmp2 = PC.XYZ() - VNorm.XYZ() * tmp * 1.e6;
942 // check intersection for Ptmp1 and Ptmp2
946 double dist1 = RealLast();
947 double dist2 = RealLast();
950 gp_Ax1 line( PC, VNorm );
951 vector< const SMDS_MeshElement* > suspectElems;
952 searcher->GetElementsNearLine( line, SMDSAbs_Face, suspectElems);
954 for ( int iF = 0; iF < suspectElems.size(); ++iF ) {
955 const SMDS_MeshElement* F = suspectElems[iF];
956 if(F==face) continue;
957 Handle(TColgp_HSequenceOfPnt) aContour = new TColgp_HSequenceOfPnt;
958 for ( int i = 0; i < 4; ++i )
959 aContour->Append( SMESH_TNodeXYZ( F->GetNode(i) ));
961 if( !volumes[0] && HasIntersection(Ptmp1, PC, PPP, aContour) ) {
963 double tmp = PC.Distance(PPP);
969 if( !volumes[1] && HasIntersection(Ptmp2, PC, PPP, aContour) ) {
971 double tmp = PC.Distance(PPP);
979 if( IsOK1 && !IsOK2 ) {
980 // using existed direction
982 else if( !IsOK1 && IsOK2 ) {
983 // using opposite direction
986 else { // IsOK1 && IsOK2
987 double tmp1 = PC.Distance(Pres1);
988 double tmp2 = PC.Distance(Pres2);
990 // using existed direction
993 // using opposite direction
998 NewFace = meshDS->AddFace( FNodes[0], FNodes[1], FNodes[2] );
1000 NewFace = meshDS->AddFace( FNodes[0], FNodes[2], FNodes[1] );
1001 storeTmpElement( NewFace );
1002 prxSubMesh->AddElement( NewFace );
1006 // Case of non-degenerated quadrangle
1008 // Find pyramid peak
1010 gp_XYZ PCbest(0., 0., 0.); // pyramid peak
1013 gp_Pnt Pbest = FindBestPoint(PN->Value(i), PN->Value(i+1), PC, VN->Value(i));
1014 PCbest += Pbest.XYZ();
1018 double height = PC.Distance(PCbest); // pyramid height to precise
1019 if ( height < 1.e-6 ) {
1020 // create new PCbest using a bit shift along VNorm
1021 PCbest = PC.XYZ() + VNorm.XYZ() * 0.001;
1022 height = PC.Distance(PCbest);
1023 if ( height < std::numeric_limits<double>::min() )
1024 return false; // batterfly element
1027 // Restrict pyramid height by intersection with other faces
1028 gp_Vec tmpDir(PC,PCbest); tmpDir.Normalize();
1029 double tmp = PN->Value(1).Distance(PN->Value(3)) + PN->Value(2).Distance(PN->Value(4));
1030 // far points: in (PC, PCbest) direction and vice-versa
1031 gp_Pnt farPnt[2] = { PC.XYZ() + tmpDir.XYZ() * tmp * 1.e6,
1032 PC.XYZ() - tmpDir.XYZ() * tmp * 1.e6 };
1033 // check intersection for farPnt1 and farPnt2
1034 bool intersected[2] = { false, false };
1035 double dist [2] = { RealLast(), RealLast() };
1038 gp_Ax1 line( PC, tmpDir );
1039 vector< const SMDS_MeshElement* > suspectElems;
1040 searcher->GetElementsNearLine( line, SMDSAbs_Face, suspectElems);
1042 for ( int iF = 0; iF < suspectElems.size(); ++iF )
1044 const SMDS_MeshElement* F = suspectElems[iF];
1045 if(F==face) continue;
1046 Handle(TColgp_HSequenceOfPnt) aContour = new TColgp_HSequenceOfPnt;
1047 int nbN = F->NbNodes() / ( F->IsQuadratic() ? 2 : 1 );
1048 for ( i = 0; i < nbN; ++i )
1049 aContour->Append( SMESH_TNodeXYZ( F->GetNode(i) ));
1051 for ( int isRev = 0; isRev < 2; ++isRev )
1053 if( !volumes[isRev] && HasIntersection(farPnt[isRev], PC, intP, aContour) ) {
1054 intersected[isRev] = true;
1055 double d = PC.Distance( intP );
1056 if( d < dist[isRev] )
1058 intPnt[isRev] = intP;
1065 // Create one or two pyramids
1067 for ( int isRev = 0; isRev < 2; ++isRev )
1069 if( !intersected[isRev] ) continue;
1070 double pyramidH = Min( height, PC.Distance(intPnt[isRev])/3.);
1071 PCbest = PC.XYZ() + tmpDir.XYZ() * (isRev ? -pyramidH : pyramidH);
1073 // create node for PCbest
1074 SMDS_MeshNode* NewNode = helper.AddNode( PCbest.X(), PCbest.Y(), PCbest.Z() );
1076 // add triangles to result map
1077 for(i=0; i<4; i++) {
1078 SMDS_MeshFace* NewFace;
1080 NewFace = meshDS->AddFace( NewNode, FNodes[i], FNodes[i+1] );
1082 NewFace = meshDS->AddFace( NewNode, FNodes[i+1], FNodes[i] );
1083 storeTmpElement( NewFace );
1084 prxSubMesh->AddElement( NewFace );
1087 SMDS_MeshVolume* aPyram;
1089 aPyram = helper.AddVolume( FNodes[0], FNodes[1], FNodes[2], FNodes[3], NewNode );
1091 aPyram = helper.AddVolume( FNodes[0], FNodes[3], FNodes[2], FNodes[1], NewNode );
1092 myPyramids.push_back(aPyram);
1094 } // end loop on all faces
1096 return Compute2ndPart(aMesh, myPyramids);
1099 //================================================================================
1101 * \brief Update created pyramids and faces to avoid their intersection
1103 //================================================================================
1105 bool StdMeshers_QuadToTriaAdaptor::Compute2ndPart(SMESH_Mesh& aMesh,
1106 const vector<const SMDS_MeshElement*>& myPyramids)
1108 if(myPyramids.empty())
1111 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
1112 int i, j, k, myShapeID = myPyramids[0]->GetNode(4)->getshapeId();
1114 if ( myElemSearcher ) delete myElemSearcher;
1115 myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher();
1116 SMESH_ElementSearcher* searcher = const_cast<SMESH_ElementSearcher*>(myElemSearcher);
1118 set<const SMDS_MeshNode*> nodesToMove;
1120 // check adjacent pyramids
1122 for ( i = 0; i < myPyramids.size(); ++i )
1124 const SMDS_MeshElement* PrmI = myPyramids[i];
1125 MergeAdjacent( PrmI, nodesToMove );
1128 // iterate on all pyramids
1129 for ( i = 0; i < myPyramids.size(); ++i )
1131 const SMDS_MeshElement* PrmI = myPyramids[i];
1133 // compare PrmI with all the rest pyramids
1135 // collect adjacent pyramids and nodes coordinates of PrmI
1136 set<const SMDS_MeshElement*> checkedPyrams;
1137 vector<gp_Pnt> PsI(5);
1138 for(k=0; k<5; k++) // loop on 4 base nodes of PrmI
1140 const SMDS_MeshNode* n = PrmI->GetNode(k);
1141 PsI[k] = SMESH_TNodeXYZ( n );
1142 SMDS_ElemIteratorPtr vIt = n->GetInverseElementIterator( SMDSAbs_Volume );
1143 while ( vIt->more() )
1145 const SMDS_MeshElement* PrmJ = vIt->next();
1146 if ( SMESH_Algo::GetCommonNodes( PrmI, PrmJ ).size() > 1 )
1147 checkedPyrams.insert( PrmJ );
1151 // check intersection with distant pyramids
1152 for(k=0; k<4; k++) // loop on 4 base nodes of PrmI
1154 gp_Vec Vtmp(PsI[k],PsI[4]);
1155 gp_Ax1 line( PsI[k], Vtmp );
1156 vector< const SMDS_MeshElement* > suspectPyrams;
1157 searcher->GetElementsNearLine( line, SMDSAbs_Volume, suspectPyrams);
1159 for ( j = 0; j < suspectPyrams.size(); ++j )
1161 const SMDS_MeshElement* PrmJ = suspectPyrams[j];
1162 if ( PrmJ == PrmI || PrmJ->NbCornerNodes() != 5 )
1164 if ( myShapeID != PrmJ->GetNode(4)->getshapeId())
1165 continue; // pyramid from other SOLID
1166 if ( PrmI->GetNode(4) == PrmJ->GetNode(4) )
1167 continue; // pyramids PrmI and PrmJ already merged
1168 if ( !checkedPyrams.insert( PrmJ ).second )
1169 continue; // already checked
1171 TXyzIterator xyzIt( PrmJ->nodesIterator() );
1172 vector<gp_Pnt> PsJ( xyzIt, TXyzIterator() );
1176 for(k=0; k<4 && !hasInt; k++) {
1177 gp_Vec Vtmp(PsI[k],PsI[4]);
1178 gp_Pnt Pshift = PsI[k].XYZ() + Vtmp.XYZ() * 0.01; // base node moved a bit to apex
1180 ( HasIntersection3( Pshift, PsI[4], Pint, PsJ[0], PsJ[1], PsJ[4]) ||
1181 HasIntersection3( Pshift, PsI[4], Pint, PsJ[1], PsJ[2], PsJ[4]) ||
1182 HasIntersection3( Pshift, PsI[4], Pint, PsJ[2], PsJ[3], PsJ[4]) ||
1183 HasIntersection3( Pshift, PsI[4], Pint, PsJ[3], PsJ[0], PsJ[4]) );
1185 for(k=0; k<4 && !hasInt; k++) {
1186 gp_Vec Vtmp(PsJ[k],PsJ[4]);
1187 gp_Pnt Pshift = PsJ[k].XYZ() + Vtmp.XYZ() * 0.01;
1189 ( HasIntersection3( Pshift, PsJ[4], Pint, PsI[0], PsI[1], PsI[4]) ||
1190 HasIntersection3( Pshift, PsJ[4], Pint, PsI[1], PsI[2], PsI[4]) ||
1191 HasIntersection3( Pshift, PsJ[4], Pint, PsI[2], PsI[3], PsI[4]) ||
1192 HasIntersection3( Pshift, PsJ[4], Pint, PsI[3], PsI[0], PsI[4]) );
1197 // count common nodes of base faces of two pyramids
1200 nbc += int ( PrmI->GetNodeIndex( PrmJ->GetNode(k) ) >= 0 );
1203 continue; // pyrams have a common base face
1207 // Merge the two pyramids and others already merged with them
1208 MergePiramids( PrmI, PrmJ, nodesToMove );
1212 // decrease height of pyramids
1213 gp_XYZ PCi(0,0,0), PCj(0,0,0);
1214 for(k=0; k<4; k++) {
1215 PCi += PsI[k].XYZ();
1216 PCj += PsJ[k].XYZ();
1219 gp_Vec VN1(PCi,PsI[4]);
1220 gp_Vec VN2(PCj,PsJ[4]);
1221 gp_Vec VI1(PCi,Pint);
1222 gp_Vec VI2(PCj,Pint);
1223 double ang1 = fabs(VN1.Angle(VI1));
1224 double ang2 = fabs(VN2.Angle(VI2));
1225 double coef1 = 0.5 - (( ang1 < M_PI/3. ) ? cos(ang1)*0.25 : 0 );
1226 double coef2 = 0.5 - (( ang2 < M_PI/3. ) ? cos(ang2)*0.25 : 0 ); // cos(ang2) ?
1227 // double coef2 = 0.5;
1229 // coef2 -= cos(ang1)*0.25;
1233 SMDS_MeshNode* aNode1 = const_cast<SMDS_MeshNode*>(PrmI->GetNode(4));
1234 aNode1->setXYZ( PCi.X()+VN1.X(), PCi.Y()+VN1.Y(), PCi.Z()+VN1.Z() );
1235 SMDS_MeshNode* aNode2 = const_cast<SMDS_MeshNode*>(PrmJ->GetNode(4));
1236 aNode2->setXYZ( PCj.X()+VN2.X(), PCj.Y()+VN2.Y(), PCj.Z()+VN2.Z() );
1237 nodesToMove.insert( aNode1 );
1238 nodesToMove.insert( aNode2 );
1240 // fix intersections that could appear after apex movement
1241 MergeAdjacent( PrmI, nodesToMove );
1242 MergeAdjacent( PrmJ, nodesToMove );
1245 } // loop on suspectPyrams
1246 } // loop on 4 base nodes of PrmI
1248 } // loop on all pyramids
1250 if( !nodesToMove.empty() && !meshDS->IsEmbeddedMode() )
1252 set<const SMDS_MeshNode*>::iterator n = nodesToMove.begin();
1253 for ( ; n != nodesToMove.end(); ++n )
1254 meshDS->MoveNode( *n, (*n)->X(), (*n)->Y(), (*n)->Z() );
1257 // move medium nodes of merged quadratic pyramids
1258 if ( myPyramids[0]->IsQuadratic() )
1259 UpdateQuadraticPyramids( nodesToMove, GetMeshDS() );
1261 // erase removed triangles from the proxy mesh
1262 if ( !myRemovedTrias.empty() )
1264 for ( int i = 0; i <= meshDS->MaxShapeIndex(); ++i )
1265 if ( SMESH_ProxyMesh::SubMesh* sm = findProxySubMesh(i))
1267 vector<const SMDS_MeshElement *> faces;
1268 faces.reserve( sm->NbElements() );
1269 SMDS_ElemIteratorPtr fIt = sm->GetElements();
1270 while ( fIt->more() )
1272 const SMDS_MeshElement* tria = fIt->next();
1273 set<const SMDS_MeshElement*>::iterator rmTria = myRemovedTrias.find( tria );
1274 if ( rmTria != myRemovedTrias.end() )
1275 myRemovedTrias.erase( rmTria );
1277 faces.push_back( tria );
1279 sm->ChangeElements( faces.begin(), faces.end() );
1285 delete myElemSearcher;