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 };
49 // std-like iterator used to get coordinates of nodes of mesh element
50 typedef SMDS_StdIterator< SMESH_MeshEditor::TNodeXYZ, SMDS_ElemIteratorPtr > TXyzIterator;
54 const int Q2TAbs_TmpTriangle = SMDSEntity_Last + 1;
56 //================================================================================
58 * \brief Temporary face. It's key feature is that it adds itself to an apex node
59 * as inverse element, so that tmp triangles of a piramid can be easily found
61 //================================================================================
63 class STDMESHERS_EXPORT Q2TAdaptor_Triangle : public SMDS_MeshFace
65 const SMDS_MeshNode* _nodes[3];
67 Q2TAdaptor_Triangle(const SMDS_MeshNode* apexNode,
68 const SMDS_MeshNode* node2,
69 const SMDS_MeshNode* node3)
71 _nodes[0]=0; ChangeApex(apexNode);
75 ~Q2TAdaptor_Triangle() { MarkAsRemoved(); }
76 void ChangeApex(const SMDS_MeshNode* node)
80 const_cast<SMDS_MeshNode*>(node)->AddInverseElement(this);
85 const_cast<SMDS_MeshNode*>(_nodes[0])->RemoveInverseElement(this), _nodes[0] = 0;
87 bool IsRemoved() const { return !_nodes[0]; }
88 virtual int NbNodes() const { return 3; }
89 virtual const SMDS_MeshNode* GetNode(const int ind) const { return _nodes[ind]; }
90 virtual SMDSAbs_EntityType GetEntityType() const
92 return SMDSAbs_EntityType( Q2TAbs_TmpTriangle );
94 virtual SMDS_ElemIteratorPtr elementsIterator(SMDSAbs_ElementType type) const
96 if ( type == SMDSAbs_Node )
97 return SMDS_ElemIteratorPtr( new SMDS_NodeArrayElemIterator( _nodes, _nodes+3 ));
98 throw SALOME_Exception(LOCALIZED("Not implemented"));
102 //================================================================================
104 * \brief Return true if two nodes of triangles are equal
106 //================================================================================
108 bool EqualTriangles(const SMDS_MeshElement* F1,const SMDS_MeshElement* F2)
111 ( F1->GetNode(1)==F2->GetNode(2) && F1->GetNode(2)==F2->GetNode(1) ) ||
112 ( F1->GetNode(1)==F2->GetNode(1) && F1->GetNode(2)==F2->GetNode(2) );
115 //================================================================================
117 * \brief Merge the two pyramids (i.e. fuse their apex) and others already merged with them
119 //================================================================================
121 void MergePiramids( const SMDS_MeshElement* PrmI,
122 const SMDS_MeshElement* PrmJ,
123 SMESHDS_Mesh* meshDS,
124 set<const SMDS_MeshNode*> & nodesToMove)
126 const SMDS_MeshNode* Nrem = PrmJ->GetNode(4); // node to remove
127 int nbJ = Nrem->NbInverseElements( SMDSAbs_Volume );
128 SMESH_MeshEditor::TNodeXYZ Pj( Nrem );
130 // an apex node to make common to all merged pyramids
131 SMDS_MeshNode* CommonNode = const_cast<SMDS_MeshNode*>(PrmI->GetNode(4));
132 if ( CommonNode == Nrem ) return; // already merged
133 int nbI = CommonNode->NbInverseElements( SMDSAbs_Volume );
134 SMESH_MeshEditor::TNodeXYZ Pi( CommonNode );
135 gp_XYZ Pnew = ( nbI*Pi + nbJ*Pj ) / (nbI+nbJ);
136 CommonNode->setXYZ( Pnew.X(), Pnew.Y(), Pnew.Z() );
138 nodesToMove.insert( CommonNode );
139 nodesToMove.erase ( Nrem );
141 // find and remove coincided faces of merged pyramids
142 SMDS_ElemIteratorPtr triItI = CommonNode->GetInverseElementIterator(SMDSAbs_Face);
143 while ( triItI->more() )
145 const SMDS_MeshElement* FI = triItI->next();
146 const SMDS_MeshElement* FJEqual = 0;
147 SMDS_ElemIteratorPtr triItJ = Nrem->GetInverseElementIterator(SMDSAbs_Face);
148 while ( !FJEqual && triItJ->more() )
150 const SMDS_MeshElement* FJ = triItJ->next();
151 if ( EqualTriangles( FJ, FI ))
156 ((Q2TAdaptor_Triangle*) FI)->MarkAsRemoved();
157 ((Q2TAdaptor_Triangle*) FJEqual)->MarkAsRemoved();
161 // set the common apex node to pyramids and triangles merged with J
162 SMDS_ElemIteratorPtr itJ = Nrem->GetInverseElementIterator();
163 while ( itJ->more() )
165 const SMDS_MeshElement* elem = itJ->next();
166 if ( elem->GetType() == SMDSAbs_Volume ) // pyramid
168 vector< const SMDS_MeshNode* > nodes( elem->begin_nodes(), elem->end_nodes() );
169 nodes[4] = CommonNode;
170 meshDS->ChangeElementNodes( elem, &nodes[0], nodes.size());
172 else if ( elem->GetEntityType() == Q2TAbs_TmpTriangle ) // tmp triangle
174 ((Q2TAdaptor_Triangle*) elem )->ChangeApex( CommonNode );
177 ASSERT( Nrem->NbInverseElements() == 0 );
178 meshDS->RemoveFreeNode( Nrem,
179 meshDS->MeshElements( Nrem->GetPosition()->GetShapeId()),
180 /*fromGroups=*/false);
183 //================================================================================
185 * \brief Return true if two adjacent pyramids are too close one to another
186 * so that a tetrahedron to built between them would have too poor quality
188 //================================================================================
190 bool TooCloseAdjacent( const SMDS_MeshElement* PrmI,
191 const SMDS_MeshElement* PrmJ,
194 const SMDS_MeshNode* nApexI = PrmI->GetNode(4);
195 const SMDS_MeshNode* nApexJ = PrmJ->GetNode(4);
196 if ( nApexI == nApexJ ||
197 nApexI->GetPosition()->GetShapeId() != nApexJ->GetPosition()->GetShapeId() )
200 // Find two common base nodes and their indices within PrmI and PrmJ
201 const SMDS_MeshNode* baseNodes[2] = { 0,0 };
202 int baseNodesIndI[2], baseNodesIndJ[2];
203 for ( int i = 0; i < 4 ; ++i )
205 int j = PrmJ->GetNodeIndex( PrmI->GetNode(i));
208 int ind = baseNodes[0] ? 1:0;
209 if ( baseNodes[ ind ])
210 return false; // pyramids with a common base face
211 baseNodes [ ind ] = PrmI->GetNode(i);
212 baseNodesIndI[ ind ] = i;
213 baseNodesIndJ[ ind ] = j;
216 if ( !baseNodes[1] ) return false; // not adjacent
218 // Get normals of triangles sharing baseNodes
219 gp_XYZ apexI = SMESH_MeshEditor::TNodeXYZ( nApexI );
220 gp_XYZ apexJ = SMESH_MeshEditor::TNodeXYZ( nApexJ );
221 gp_XYZ base1 = SMESH_MeshEditor::TNodeXYZ( baseNodes[0]);
222 gp_XYZ base2 = SMESH_MeshEditor::TNodeXYZ( baseNodes[1]);
223 gp_Vec baseVec( base1, base2 );
224 gp_Vec baI( base1, apexI );
225 gp_Vec baJ( base1, apexJ );
226 gp_Vec nI = baseVec.Crossed( baI );
227 gp_Vec nJ = baseVec.Crossed( baJ );
229 // Check angle between normals
230 double angle = nI.Angle( nJ );
231 bool tooClose = ( angle < 15 * PI180 );
233 // Check if pyramids collide
235 if ( !tooClose && baI * baJ > 0 )
237 // find out if nI points outside of PrmI or inside
238 int dInd = baseNodesIndI[1] - baseNodesIndI[0];
239 isOutI = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0;
241 // find out sign of projection of nJ to baI
242 double proj = baI * nJ;
244 tooClose = isOutI ? proj > 0 : proj < 0;
247 // Check if PrmI and PrmJ are in same domain
248 if ( tooClose && !hasShape )
250 // check order of baseNodes within pyramids, it must be opposite
251 int dInd = baseNodesIndJ[1] - baseNodesIndJ[0];
252 isOutJ = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0;
253 if ( isOutJ == isOutI )
254 return false; // other domain
256 // check absence of a face separating domains between pyramids
257 TIDSortedElemSet emptySet, avoidSet;
259 while ( const SMDS_MeshElement* f =
260 SMESH_MeshEditor::FindFaceInSet( baseNodes[0], baseNodes[1],
261 emptySet, avoidSet, &i1, &i2 ))
263 avoidSet.insert( f );
265 // face node other than baseNodes
266 int otherNodeInd = 0;
267 while ( otherNodeInd == i1 || otherNodeInd == i2 ) otherNodeInd++;
268 const SMDS_MeshNode* otherFaceNode = f->GetNode( otherNodeInd );
270 // check if f is a base face of either of pyramids
271 if ( f->NbCornerNodes() == 4 &&
272 ( PrmI->GetNodeIndex( otherFaceNode ) >= 0 ||
273 PrmJ->GetNodeIndex( otherFaceNode ) >= 0 ))
274 continue; // f is a base quadrangle
276 // check projections of face direction (baOFN) to triange normals (nI and nJ)
277 gp_Vec baOFN( base1, SMESH_MeshEditor::TNodeXYZ( otherFaceNode ));
278 ( isOutI ? nJ : nI ).Reverse();
279 if ( nI * baOFN > 0 && nJ * baOFN > 0 )
281 tooClose = false; // f is between pyramids
290 //================================================================================
292 * \brief Merges adjacent pyramids
294 //================================================================================
296 void MergeAdjacent(const SMDS_MeshElement* PrmI,
298 set<const SMDS_MeshNode*>& nodesToMove)
300 TIDSortedElemSet adjacentPyrams, mergedPyrams;
301 for(int k=0; k<4; k++) // loop on 4 base nodes of PrmI
303 const SMDS_MeshNode* n = PrmI->GetNode(k);
304 SMDS_ElemIteratorPtr vIt = n->GetInverseElementIterator( SMDSAbs_Volume );
305 while ( vIt->more() )
307 const SMDS_MeshElement* PrmJ = vIt->next();
308 if ( PrmJ->NbCornerNodes() != 5 || !adjacentPyrams.insert( PrmJ ).second )
310 if ( PrmI != PrmJ && TooCloseAdjacent( PrmI, PrmJ, mesh.HasShapeToMesh() ))
312 MergePiramids( PrmI, PrmJ, mesh.GetMeshDS(), nodesToMove );
313 mergedPyrams.insert( PrmJ );
317 if ( !mergedPyrams.empty() )
319 TIDSortedElemSet::iterator prm;
320 // for (prm = mergedPyrams.begin(); prm != mergedPyrams.end(); ++prm)
321 // MergeAdjacent( *prm, mesh, nodesToMove );
323 for (prm = adjacentPyrams.begin(); prm != adjacentPyrams.end(); ++prm)
324 MergeAdjacent( *prm, mesh, nodesToMove );
329 //================================================================================
333 //================================================================================
335 StdMeshers_QuadToTriaAdaptor::StdMeshers_QuadToTriaAdaptor():
336 myElemSearcher(0), myNbTriangles(0)
340 //================================================================================
344 //================================================================================
346 StdMeshers_QuadToTriaAdaptor::~StdMeshers_QuadToTriaAdaptor()
348 // delete temporary faces
349 TQuad2Trias::iterator f_f = myResMap.begin(), ffEnd = myResMap.end();
350 for ( ; f_f != ffEnd; ++f_f )
352 TTriaList& fList = f_f->second;
353 TTriaList::iterator f = fList.begin(), fEnd = fList.end();
354 for ( ; f != fEnd; ++f )
359 if ( myElemSearcher ) delete myElemSearcher;
364 //=======================================================================
365 //function : FindBestPoint
366 //purpose : Return a point P laying on the line (PC,V) so that triangle
367 // (P, P1, P2) to be equilateral as much as possible
368 // V - normal to (P1,P2,PC)
369 //=======================================================================
370 static gp_Pnt FindBestPoint(const gp_Pnt& P1, const gp_Pnt& P2,
371 const gp_Pnt& PC, const gp_Vec& V)
373 double a = P1.Distance(P2);
374 double b = P1.Distance(PC);
375 double c = P2.Distance(PC);
379 // find shift along V in order a to became equal to (b+c)/2
380 double shift = sqrt( a*a + (b*b-c*c)*(b*b-c*c)/16/a/a - (b*b+c*c)/2 );
382 gp_Pnt Pbest = PC.XYZ() + aDir.XYZ() * shift;
388 //=======================================================================
389 //function : HasIntersection3
390 //purpose : Auxilare for HasIntersection()
391 // find intersection point between triangle (P1,P2,P3)
392 // and segment [PC,P]
393 //=======================================================================
394 static bool HasIntersection3(const gp_Pnt& P, const gp_Pnt& PC, gp_Pnt& Pint,
395 const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3)
397 //cout<<"HasIntersection3"<<endl;
398 //cout<<" PC("<<PC.X()<<","<<PC.Y()<<","<<PC.Z()<<")"<<endl;
399 //cout<<" P("<<P.X()<<","<<P.Y()<<","<<P.Z()<<")"<<endl;
400 //cout<<" P1("<<P1.X()<<","<<P1.Y()<<","<<P1.Z()<<")"<<endl;
401 //cout<<" P2("<<P2.X()<<","<<P2.Y()<<","<<P2.Z()<<")"<<endl;
402 //cout<<" P3("<<P3.X()<<","<<P3.Y()<<","<<P3.Z()<<")"<<endl;
405 IntAna_Quadric IAQ(gp_Pln(P1,VP1.Crossed(VP2)));
406 IntAna_IntConicQuad IAICQ(gp_Lin(PC,gp_Dir(gp_Vec(PC,P))),IAQ);
408 if( IAICQ.IsInQuadric() )
410 if( IAICQ.NbPoints() == 1 ) {
411 gp_Pnt PIn = IAICQ.Point(1);
412 const double preci = 1.e-10 * P.Distance(PC);
413 // check if this point is internal for segment [PC,P]
415 ( (PC.X()-PIn.X())*(P.X()-PIn.X()) > preci ) ||
416 ( (PC.Y()-PIn.Y())*(P.Y()-PIn.Y()) > preci ) ||
417 ( (PC.Z()-PIn.Z())*(P.Z()-PIn.Z()) > preci );
421 // check if this point is internal for triangle (P1,P2,P3)
425 if( V1.Magnitude()<preci ||
426 V2.Magnitude()<preci ||
427 V3.Magnitude()<preci ) {
431 const double angularTol = 1e-6;
432 gp_Vec VC1 = V1.Crossed(V2);
433 gp_Vec VC2 = V2.Crossed(V3);
434 gp_Vec VC3 = V3.Crossed(V1);
435 if(VC1.Magnitude()<gp::Resolution()) {
436 if(VC2.IsOpposite(VC3,angularTol)) {
440 else if(VC2.Magnitude()<gp::Resolution()) {
441 if(VC1.IsOpposite(VC3,angularTol)) {
445 else if(VC3.Magnitude()<gp::Resolution()) {
446 if(VC1.IsOpposite(VC2,angularTol)) {
451 if( VC1.IsOpposite(VC2,angularTol) || VC1.IsOpposite(VC3,angularTol) ||
452 VC2.IsOpposite(VC3,angularTol) ) {
465 //=======================================================================
466 //function : HasIntersection
467 //purpose : Auxilare for CheckIntersection()
468 //=======================================================================
470 static bool HasIntersection(const gp_Pnt& P, const gp_Pnt& PC, gp_Pnt& Pint,
471 Handle(TColgp_HSequenceOfPnt)& aContour)
473 if(aContour->Length()==3) {
474 return HasIntersection3( P, PC, Pint, aContour->Value(1),
475 aContour->Value(2), aContour->Value(3) );
479 if( (aContour->Value(1).Distance(aContour->Value(2)) > 1.e-6) &&
480 (aContour->Value(1).Distance(aContour->Value(3)) > 1.e-6) &&
481 (aContour->Value(2).Distance(aContour->Value(3)) > 1.e-6) ) {
482 check = HasIntersection3( P, PC, Pint, aContour->Value(1),
483 aContour->Value(2), aContour->Value(3) );
485 if(check) return true;
486 if( (aContour->Value(1).Distance(aContour->Value(4)) > 1.e-6) &&
487 (aContour->Value(1).Distance(aContour->Value(3)) > 1.e-6) &&
488 (aContour->Value(4).Distance(aContour->Value(3)) > 1.e-6) ) {
489 check = HasIntersection3( P, PC, Pint, aContour->Value(1),
490 aContour->Value(3), aContour->Value(4) );
492 if(check) return true;
498 //================================================================================
500 * \brief Checks if a line segment (P,PC) intersects any mesh face.
501 * \param P - first segment end
502 * \param PC - second segment end (it is a gravity center of quadrangle)
503 * \param Pint - (out) intersection point
504 * \param aMesh - mesh
505 * \param aShape - shape to check faces on
506 * \param NotCheckedFace - mesh face not to check
507 * \retval bool - true if there is an intersection
509 //================================================================================
511 bool StdMeshers_QuadToTriaAdaptor::CheckIntersection (const gp_Pnt& P,
515 const TopoDS_Shape& aShape,
516 const SMDS_MeshElement* NotCheckedFace)
518 if ( !myElemSearcher )
519 myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher();
520 SMESH_ElementSearcher* searcher = const_cast<SMESH_ElementSearcher*>(myElemSearcher);
522 //SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
523 //cout<<" CheckIntersection: meshDS->NbFaces() = "<<meshDS->NbFaces()<<endl;
525 double dist = RealLast(); // find intersection closest to the segment
528 gp_Ax1 line( P, gp_Vec(P,PC));
529 vector< const SMDS_MeshElement* > suspectElems;
530 searcher->GetElementsNearLine( line, SMDSAbs_Face, suspectElems);
532 // for (TopExp_Explorer exp(aShape,TopAbs_FACE);exp.More();exp.Next()) {
533 // const TopoDS_Shape& aShapeFace = exp.Current();
534 // if(aShapeFace==NotCheckedFace)
536 // const SMESHDS_SubMesh * aSubMeshDSFace = meshDS->MeshElements(aShapeFace);
537 // if ( aSubMeshDSFace ) {
538 // SMDS_ElemIteratorPtr iteratorElem = aSubMeshDSFace->GetElements();
539 // while ( iteratorElem->more() ) { // loop on elements on a face
540 // const SMDS_MeshElement* face = iteratorElem->next();
541 for ( int i = 0; i < suspectElems.size(); ++i )
543 const SMDS_MeshElement* face = suspectElems[i];
544 if ( face == NotCheckedFace ) continue;
545 Handle(TColgp_HSequenceOfPnt) aContour = new TColgp_HSequenceOfPnt;
546 for ( int i = 0; i < face->NbCornerNodes(); ++i )
547 aContour->Append( SMESH_MeshEditor::TNodeXYZ( face->GetNode(i) ));
548 if( HasIntersection(P, PC, Pres, aContour) ) {
550 double tmp = PC.Distance(Pres);
560 //================================================================================
562 * \brief Prepare data for the given face
563 * \param PN - coordinates of face nodes
564 * \param VN - cross products of vectors (PC-PN(i)) ^ (PC-PN(i+1))
565 * \param FNodes - face nodes
566 * \param PC - gravity center of nodes
567 * \param VNorm - face normal (sum of VN)
568 * \param volumes - two volumes sharing the given face, the first is in VNorm direction
569 * \retval int - 0 if given face is not quad,
570 * 1 if given face is quad,
571 * 2 if given face is degenerate quad (two nodes are coincided)
573 //================================================================================
575 int StdMeshers_QuadToTriaAdaptor::Preparation(const SMDS_MeshElement* face,
576 Handle(TColgp_HArray1OfPnt)& PN,
577 Handle(TColgp_HArray1OfVec)& VN,
578 vector<const SMDS_MeshNode*>& FNodes,
581 const SMDS_MeshElement** volumes)
583 if( face->NbCornerNodes() != 4 )
585 myNbTriangles += int( face->NbCornerNodes() == 3 );
590 gp_XYZ xyzC(0., 0., 0.);
591 for ( i = 0; i < 4; ++i )
593 gp_XYZ p = SMESH_MeshEditor::TNodeXYZ( FNodes[i] = face->GetNode(i) );
594 PN->SetValue( i+1, p );
598 //cout<<" PC("<<PC.X()<<","<<PC.Y()<<","<<PC.Z()<<")"<<endl;
606 if( PN->Value(i).Distance(PN->Value(j)) < 1.e-6 )
611 //int deg_num = IsDegenarate(PN);
615 //cout<<"find degeneration"<<endl;
617 gp_Pnt Pdeg = PN->Value(i);
619 list< const SMDS_MeshNode* >::iterator itdg = myDegNodes.begin();
620 const SMDS_MeshNode* DegNode = 0;
621 for(; itdg!=myDegNodes.end(); itdg++) {
622 const SMDS_MeshNode* N = (*itdg);
623 gp_Pnt Ptmp(N->X(),N->Y(),N->Z());
624 if(Pdeg.Distance(Ptmp)<1.e-6) {
626 //DegNode = const_cast<SMDS_MeshNode*>(N);
631 DegNode = FNodes[i-1];
632 myDegNodes.push_back(DegNode);
635 FNodes[i-1] = DegNode;
638 PN->SetValue(i,PN->Value(i+1));
639 FNodes[i-1] = FNodes[i];
644 PN->SetValue(nbp+1,PN->Value(1));
645 FNodes[nbp] = FNodes[0];
646 // find normal direction
647 gp_Vec V1(PC,PN->Value(nbp));
648 gp_Vec V2(PC,PN->Value(1));
649 VNorm = V1.Crossed(V2);
650 VN->SetValue(nbp,VNorm);
651 for(i=1; i<nbp; i++) {
652 V1 = gp_Vec(PC,PN->Value(i));
653 V2 = gp_Vec(PC,PN->Value(i+1));
654 gp_Vec Vtmp = V1.Crossed(V2);
655 VN->SetValue(i,Vtmp);
659 // find volumes sharing the face
662 volumes[0] = volumes[1] = 0;
663 SMDS_ElemIteratorPtr vIt = FNodes[0]->GetInverseElementIterator( SMDSAbs_Volume );
664 while ( vIt->more() )
666 const SMDS_MeshElement* vol = vIt->next();
667 bool volSharesAllNodes = true;
668 for ( int i = 1; i < face->NbNodes() && volSharesAllNodes; ++i )
669 volSharesAllNodes = ( vol->GetNodeIndex( FNodes[i] ) >= 0 );
670 if ( volSharesAllNodes )
671 volumes[ volumes[0] ? 1 : 0 ] = vol;
672 // we could additionally check that vol has all FNodes in its one face using SMDS_VolumeTool
674 // define volume position relating to the face normal
678 SMDS_ElemIteratorPtr nodeIt = volumes[0]->nodesIterator();
680 volGC = accumulate( TXyzIterator(nodeIt), TXyzIterator(), volGC ) / volumes[0]->NbNodes();
682 if ( VNorm * gp_Vec( PC, volGC ) < 0 )
683 swap( volumes[0], volumes[1] );
687 //cout<<" VNorm("<<VNorm.X()<<","<<VNorm.Y()<<","<<VNorm.Z()<<")"<<endl;
688 return hasdeg ? DEGEN_QUAD : QUAD;
692 //=======================================================================
695 //=======================================================================
697 bool StdMeshers_QuadToTriaAdaptor::Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape)
704 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
705 SMESH_MesherHelper helper(aMesh);
706 helper.IsQuadraticSubMesh(aShape);
707 helper.SetElementsOnShape( true );
709 for (TopExp_Explorer exp(aShape,TopAbs_FACE);exp.More();exp.Next())
711 const TopoDS_Shape& aShapeFace = exp.Current();
712 const SMESHDS_SubMesh * aSubMeshDSFace = meshDS->MeshElements( aShapeFace );
713 if ( aSubMeshDSFace )
715 bool isRev = SMESH_Algo::IsReversedSubMesh( TopoDS::Face(aShapeFace), meshDS );
717 SMDS_ElemIteratorPtr iteratorElem = aSubMeshDSFace->GetElements();
718 while ( iteratorElem->more() ) // loop on elements on a geometrical face
720 const SMDS_MeshElement* face = iteratorElem->next();
721 //cout<<endl<<"================= face->GetID() = "<<face->GetID()<<endl;
722 // preparation step using face info
723 Handle(TColgp_HArray1OfPnt) PN = new TColgp_HArray1OfPnt(1,5);
724 Handle(TColgp_HArray1OfVec) VN = new TColgp_HArray1OfVec(1,4);
725 vector<const SMDS_MeshNode*> FNodes(5);
728 int stat = Preparation(face, PN, VN, FNodes, PC, VNorm);
735 // add triangles to result map
736 SMDS_MeshFace* NewFace;
738 NewFace = new Q2TAdaptor_Triangle( FNodes[0], FNodes[1], FNodes[2] );
740 NewFace = new Q2TAdaptor_Triangle( FNodes[0], FNodes[2], FNodes[1] );
741 TTriaList aList( 1, NewFace );
742 myResMap.insert(make_pair(face,aList));
746 if(!isRev) VNorm.Reverse();
747 double xc = 0., yc = 0., zc = 0.;
752 Pbest = FindBestPoint(PN->Value(i), PN->Value(i+1), PC, VN->Value(i).Reversed());
754 Pbest = FindBestPoint(PN->Value(i), PN->Value(i+1), PC, VN->Value(i));
759 gp_Pnt PCbest(xc/4., yc/4., zc/4.);
762 double height = PCbest.Distance(PC);
764 // create new PCbest using a bit shift along VNorm
765 PCbest = PC.XYZ() + VNorm.XYZ() * 0.001;
768 // check possible intersection with other faces
770 bool check = CheckIntersection(PCbest, PC, Pint, aMesh, aShape, face);
772 //cout<<"--PC("<<PC.X()<<","<<PC.Y()<<","<<PC.Z()<<")"<<endl;
773 //cout<<" PCbest("<<PCbest.X()<<","<<PCbest.Y()<<","<<PCbest.Z()<<")"<<endl;
774 double dist = PC.Distance(Pint)/3.;
775 gp_Dir aDir(gp_Vec(PC,PCbest));
776 PCbest = PC.XYZ() + aDir.XYZ() * dist;
779 gp_Vec VB(PC,PCbest);
780 gp_Pnt PCbestTmp = PC.XYZ() + VB.XYZ() * 3.0;
781 check = CheckIntersection(PCbestTmp, PC, Pint, aMesh, aShape, face);
783 double dist = PC.Distance(Pint)/3.;
785 gp_Dir aDir(gp_Vec(PC,PCbest));
786 PCbest = PC.XYZ() + aDir.XYZ() * dist;
791 // create node for PCbest
792 SMDS_MeshNode* NewNode = helper.AddNode( PCbest.X(), PCbest.Y(), PCbest.Z() );
794 // add triangles to result map
795 TTriaList& triaList = myResMap.insert( make_pair( face, TTriaList() ))->second;
797 triaList.push_back( new Q2TAdaptor_Triangle( NewNode, FNodes[i], FNodes[i+1] ));
800 if ( isRev ) swap( FNodes[1], FNodes[3]);
801 SMDS_MeshVolume* aPyram =
802 helper.AddVolume( FNodes[0], FNodes[1], FNodes[2], FNodes[3], NewNode );
803 myPyramids.push_back(aPyram);
805 } // end loop on elements on a face submesh
807 } // end for(TopExp_Explorer exp(aShape,TopAbs_FACE);exp.More();exp.Next()) {
809 return Compute2ndPart(aMesh);
812 //================================================================================
814 * \brief Computes pyramids in mesh with no shape
816 //================================================================================
818 bool StdMeshers_QuadToTriaAdaptor::Compute(SMESH_Mesh& aMesh)
822 SMESH_MesherHelper helper(aMesh);
823 helper.IsQuadraticSubMesh(aMesh.GetShapeToMesh());
824 helper.SetElementsOnShape( true );
826 if ( !myElemSearcher )
827 myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher();
828 SMESH_ElementSearcher* searcher = const_cast<SMESH_ElementSearcher*>(myElemSearcher);
830 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
832 SMDS_FaceIteratorPtr fIt = meshDS->facesIterator(/*idInceasingOrder=*/true);
835 const SMDS_MeshElement* face = fIt->next();
836 if ( !face ) continue;
837 //cout<<endl<<"================= face->GetID() = "<<face->GetID()<<endl;
838 // retrieve needed information about a face
839 Handle(TColgp_HArray1OfPnt) PN = new TColgp_HArray1OfPnt(1,5);
840 Handle(TColgp_HArray1OfVec) VN = new TColgp_HArray1OfVec(1,4);
841 vector<const SMDS_MeshNode*> FNodes(5);
844 const SMDS_MeshElement* volumes[2];
845 int what = Preparation(face, PN, VN, FNodes, PC, VNorm, volumes);
846 if ( what == NOT_QUAD )
848 if ( volumes[0] && volumes[1] )
849 continue; // face is shared by two volumes - no space for a pyramid
851 if ( what == DEGEN_QUAD )
854 // add triangles to result map
856 SMDS_MeshFace* NewFace;
859 double tmp = PN->Value(1).Distance(PN->Value(2)) + PN->Value(2).Distance(PN->Value(3));
860 // far points in VNorm direction
861 gp_Pnt Ptmp1 = PC.XYZ() + VNorm.XYZ() * tmp * 1.e6;
862 gp_Pnt Ptmp2 = PC.XYZ() - VNorm.XYZ() * tmp * 1.e6;
863 // check intersection for Ptmp1 and Ptmp2
867 double dist1 = RealLast();
868 double dist2 = RealLast();
871 gp_Ax1 line( PC, VNorm );
872 vector< const SMDS_MeshElement* > suspectElems;
873 searcher->GetElementsNearLine( line, SMDSAbs_Face, suspectElems);
875 for ( int iF = 0; iF < suspectElems.size(); ++iF ) {
876 const SMDS_MeshElement* F = suspectElems[iF];
877 if(F==face) continue;
878 Handle(TColgp_HSequenceOfPnt) aContour = new TColgp_HSequenceOfPnt;
879 for ( int i = 0; i < 4; ++i )
880 aContour->Append( SMESH_MeshEditor::TNodeXYZ( F->GetNode(i) ));
882 if( !volumes[0] && HasIntersection(Ptmp1, PC, PPP, aContour) ) {
884 double tmp = PC.Distance(PPP);
890 if( !volumes[1] && HasIntersection(Ptmp2, PC, PPP, aContour) ) {
892 double tmp = PC.Distance(PPP);
900 if( IsOK1 && !IsOK2 ) {
901 // using existed direction
903 else if( !IsOK1 && IsOK2 ) {
904 // using opposite direction
907 else { // IsOK1 && IsOK2
908 double tmp1 = PC.Distance(Pres1);
909 double tmp2 = PC.Distance(Pres2);
911 // using existed direction
914 // using opposite direction
919 NewFace = new Q2TAdaptor_Triangle( FNodes[0], FNodes[1], FNodes[2] );
921 NewFace = new Q2TAdaptor_Triangle( FNodes[0], FNodes[2], FNodes[1] );
922 aList.push_back(NewFace);
923 myResMap.insert(make_pair(face,aList));
929 gp_XYZ PCbest(0., 0., 0.); // pyramid peak
932 gp_Pnt Pbest = FindBestPoint(PN->Value(i), PN->Value(i+1), PC, VN->Value(i));
933 PCbest += Pbest.XYZ();
937 double height = PC.Distance(PCbest); // pyramid height to precise
939 // create new PCbest using a bit shift along VNorm
940 PCbest = PC.XYZ() + VNorm.XYZ() * 0.001;
941 height = PC.Distance(PCbest);
943 //cout<<" PCbest("<<PCbest.X()<<","<<PCbest.Y()<<","<<PCbest.Z()<<")"<<endl;
945 // Restrict pyramid height by intersection with other faces
946 gp_Vec tmpDir(PC,PCbest); tmpDir.Normalize();
947 double tmp = PN->Value(1).Distance(PN->Value(3)) + PN->Value(2).Distance(PN->Value(4));
948 // far points: in (PC, PCbest) direction and vice-versa
949 gp_Pnt farPnt[2] = { PC.XYZ() + tmpDir.XYZ() * tmp * 1.e6,
950 PC.XYZ() - tmpDir.XYZ() * tmp * 1.e6 };
951 // check intersection for farPnt1 and farPnt2
952 bool intersected[2] = { false, false };
953 double dist [2] = { RealLast(), RealLast() };
956 gp_Ax1 line( PC, tmpDir );
957 vector< const SMDS_MeshElement* > suspectElems;
958 searcher->GetElementsNearLine( line, SMDSAbs_Face, suspectElems);
960 for ( int iF = 0; iF < suspectElems.size(); ++iF )
962 const SMDS_MeshElement* F = suspectElems[iF];
963 if(F==face) continue;
964 Handle(TColgp_HSequenceOfPnt) aContour = new TColgp_HSequenceOfPnt;
965 int nbN = F->NbNodes() / ( F->IsQuadratic() ? 2 : 1 );
966 for ( i = 0; i < nbN; ++i )
967 aContour->Append( SMESH_MeshEditor::TNodeXYZ( F->GetNode(i) ));
969 for ( int isRev = 0; isRev < 2; ++isRev )
971 if( !volumes[isRev] && HasIntersection(farPnt[isRev], PC, intP, aContour) ) {
972 intersected[isRev] = true;
973 double d = PC.Distance( intP );
974 if( d < dist[isRev] )
976 intPnt[isRev] = intP;
983 // Create one or two pyramids
985 for ( int isRev = 0; isRev < 2; ++isRev )
987 if( !intersected[isRev] ) continue;
988 double pyramidH = Min( height, PC.Distance(intPnt[isRev])/3.);
989 PCbest = PC.XYZ() + tmpDir.XYZ() * (isRev ? -pyramidH : pyramidH);
991 // create node for PCbest
992 SMDS_MeshNode* NewNode = helper.AddNode( PCbest.X(), PCbest.Y(), PCbest.Z() );
994 // add triangles to result map
995 TTriaList& aList = myResMap.insert( make_pair( face, TTriaList()))->second;
997 SMDS_MeshFace* NewFace;
999 NewFace = new Q2TAdaptor_Triangle( NewNode, FNodes[i], FNodes[i+1] );
1001 NewFace = new Q2TAdaptor_Triangle( NewNode, FNodes[i+1], FNodes[i] );
1002 aList.push_back(NewFace);
1005 SMDS_MeshVolume* aPyram;
1007 aPyram = helper.AddVolume( FNodes[0], FNodes[1], FNodes[2], FNodes[3], NewNode );
1009 aPyram = helper.AddVolume( FNodes[0], FNodes[3], FNodes[2], FNodes[1], NewNode );
1010 myPyramids.push_back(aPyram);
1012 } // end loop on all faces
1014 return Compute2ndPart(aMesh);
1017 //================================================================================
1019 * \brief Update created pyramids and faces to avoid their intersection
1021 //================================================================================
1023 bool StdMeshers_QuadToTriaAdaptor::Compute2ndPart(SMESH_Mesh& aMesh)
1025 if(myPyramids.empty())
1028 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
1029 int i, j, k, myShapeID = myPyramids[0]->GetNode(4)->GetPosition()->GetShapeId();
1031 if ( !myElemSearcher )
1032 myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher();
1033 SMESH_ElementSearcher* searcher = const_cast<SMESH_ElementSearcher*>(myElemSearcher);
1035 set<const SMDS_MeshNode*> nodesToMove;
1037 // check adjacent pyramids
1039 for ( i = 0; i < myPyramids.size(); ++i )
1041 const SMDS_MeshElement* PrmI = myPyramids[i];
1042 MergeAdjacent( PrmI, aMesh, nodesToMove );
1045 // iterate on all pyramids
1046 for ( i = 0; i < myPyramids.size(); ++i )
1048 const SMDS_MeshElement* PrmI = myPyramids[i];
1050 // compare PrmI with all the rest pyramids
1052 // collect adjacent pyramids and nodes coordinates of PrmI
1053 set<const SMDS_MeshElement*> checkedPyrams;
1054 vector<gp_Pnt> PsI(5);
1055 for(k=0; k<5; k++) // loop on 4 base nodes of PrmI
1057 const SMDS_MeshNode* n = PrmI->GetNode(k);
1058 PsI[k] = SMESH_MeshEditor::TNodeXYZ( n );
1059 SMDS_ElemIteratorPtr vIt = n->GetInverseElementIterator( SMDSAbs_Volume );
1060 while ( vIt->more() )
1061 checkedPyrams.insert( vIt->next() );
1064 // check intersection with distant pyramids
1065 for(k=0; k<4; k++) // loop on 4 base nodes of PrmI
1067 gp_Vec Vtmp(PsI[k],PsI[4]);
1068 gp_Pnt Pshift = PsI[k].XYZ() + Vtmp.XYZ() * 0.01; // base node moved a bit to apex
1070 gp_Ax1 line( PsI[k], Vtmp );
1071 vector< const SMDS_MeshElement* > suspectPyrams;
1072 searcher->GetElementsNearLine( line, SMDSAbs_Volume, suspectPyrams);
1074 for ( j = 0; j < suspectPyrams.size(); ++j )
1076 const SMDS_MeshElement* PrmJ = suspectPyrams[j];
1077 if ( PrmJ == PrmI || PrmJ->NbCornerNodes() != 5 )
1079 if ( myShapeID != PrmJ->GetNode(4)->GetPosition()->GetShapeId())
1080 continue; // pyramid from other SOLID
1081 if ( PrmI->GetNode(4) == PrmJ->GetNode(4) )
1082 continue; // pyramids PrmI and PrmJ already merged
1083 if ( !checkedPyrams.insert( PrmJ ).second )
1084 continue; // already checked
1086 TXyzIterator xyzIt( PrmJ->nodesIterator() );
1087 vector<gp_Pnt> PsJ( xyzIt, TXyzIterator() );
1091 ( HasIntersection3( Pshift, PsI[4], Pint, PsJ[0], PsJ[1], PsJ[4]) ||
1092 HasIntersection3( Pshift, PsI[4], Pint, PsJ[1], PsJ[2], PsJ[4]) ||
1093 HasIntersection3( Pshift, PsI[4], Pint, PsJ[2], PsJ[3], PsJ[4]) ||
1094 HasIntersection3( Pshift, PsI[4], Pint, PsJ[3], PsJ[0], PsJ[4]) );
1096 for(k=0; k<4 && !hasInt; k++) {
1097 gp_Vec Vtmp(PsJ[k],PsJ[4]);
1098 gp_Pnt Pshift = PsJ[k].XYZ() + Vtmp.XYZ() * 0.01;
1100 ( HasIntersection3( Pshift, PsJ[4], Pint, PsI[0], PsI[1], PsI[4]) ||
1101 HasIntersection3( Pshift, PsJ[4], Pint, PsI[1], PsI[2], PsI[4]) ||
1102 HasIntersection3( Pshift, PsJ[4], Pint, PsI[2], PsI[3], PsI[4]) ||
1103 HasIntersection3( Pshift, PsJ[4], Pint, PsI[3], PsI[0], PsI[4]) );
1108 // count common nodes of base faces of two pyramids
1111 nbc += int ( PrmI->GetNodeIndex( PrmJ->GetNode(k) ) >= 0 );
1114 continue; // pyrams have a common base face
1118 // Merge the two pyramids and others already merged with them
1119 MergePiramids( PrmI, PrmJ, meshDS, nodesToMove );
1123 // decrease height of pyramids
1124 gp_XYZ PCi(0,0,0), PCj(0,0,0);
1125 for(k=0; k<4; k++) {
1126 PCi += PsI[k].XYZ();
1127 PCj += PsJ[k].XYZ();
1130 gp_Vec VN1(PCi,PsI[4]);
1131 gp_Vec VN2(PCj,PsJ[4]);
1132 gp_Vec VI1(PCi,Pint);
1133 gp_Vec VI2(PCj,Pint);
1134 double ang1 = fabs(VN1.Angle(VI1));
1135 double ang2 = fabs(VN2.Angle(VI2));
1136 double coef1 = 0.5 - (( ang1<PI/3 ) ? cos(ang1)*0.25 : 0 );
1137 double coef2 = 0.5 - (( ang2<PI/3 ) ? cos(ang2)*0.25 : 0 ); // cos(ang2) ?
1138 // double coef2 = 0.5;
1140 // coef2 -= cos(ang1)*0.25;
1144 SMDS_MeshNode* aNode1 = const_cast<SMDS_MeshNode*>(PrmI->GetNode(4));
1145 aNode1->setXYZ( PCi.X()+VN1.X(), PCi.Y()+VN1.Y(), PCi.Z()+VN1.Z() );
1146 SMDS_MeshNode* aNode2 = const_cast<SMDS_MeshNode*>(PrmJ->GetNode(4));
1147 aNode2->setXYZ( PCj.X()+VN2.X(), PCj.Y()+VN2.Y(), PCj.Z()+VN2.Z() );
1148 nodesToMove.insert( aNode1 );
1149 nodesToMove.insert( aNode2 );
1151 // fix intersections that could appear after apex movement
1152 MergeAdjacent( PrmI, aMesh, nodesToMove );
1153 MergeAdjacent( PrmJ, aMesh, nodesToMove );
1156 } // loop on suspectPyrams
1157 } // loop on 4 base nodes of PrmI
1159 } // loop on all pyramids
1161 if( !nodesToMove.empty() && !meshDS->IsEmbeddedMode() )
1163 set<const SMDS_MeshNode*>::iterator n = nodesToMove.begin();
1164 for ( ; n != nodesToMove.end(); ++n )
1165 meshDS->MoveNode( *n, (*n)->X(), (*n)->Y(), (*n)->Z() );
1168 // rebind triangles of pyramids sharing the same base quadrangle to the first
1169 // entrance of the base quadrangle
1170 TQuad2Trias::iterator q2t = myResMap.begin(), q2tPrev = q2t;
1171 for ( ++q2t; q2t != myResMap.end(); ++q2t, ++q2tPrev )
1173 if ( q2t->first == q2tPrev->first )
1174 q2tPrev->second.splice( q2tPrev->second.end(), q2t->second );
1176 // delete removed triangles and count resulting nb of triangles
1177 for ( q2t = myResMap.begin(); q2t != myResMap.end(); ++q2t )
1179 TTriaList & trias = q2t->second;
1180 for ( TTriaList::iterator tri = trias.begin(); tri != trias.end(); )
1181 if ( ((const Q2TAdaptor_Triangle*) *tri)->IsRemoved() )
1182 delete *tri, trias.erase( tri++ );
1184 tri++, myNbTriangles++;
1187 myPyramids.clear(); // no more needed
1190 delete myElemSearcher;
1196 //================================================================================
1198 * \brief Return list of created triangles for given face
1200 //================================================================================
1202 const list<const SMDS_MeshFace* >* StdMeshers_QuadToTriaAdaptor::GetTriangles (const SMDS_MeshElement* aQuad)
1204 TQuad2Trias::iterator it = myResMap.find(aQuad);
1205 return ( it != myResMap.end() ? & it->second : 0 );