Salome HOME
Update of CheckDone
[modules/smesh.git] / src / StdMeshers / StdMeshers_QuadToTriaAdaptor.cxx
1 // Copyright (C) 2007-2022  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
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, or (at your option) any later version.
7 //
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.
12 //
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
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 // File      : StdMeshers_QuadToTriaAdaptor.cxx
21 // Module    : SMESH
22 // Created   : Wen May 07 16:37:07 2008
23 // Author    : Sergey KUUL (skl)
24
25 #include "StdMeshers_QuadToTriaAdaptor.hxx"
26
27 #include "SMDS_IteratorOnIterators.hxx"
28 #include "SMDS_SetIterator.hxx"
29 #include "SMESHDS_GroupBase.hxx"
30 #include "SMESHDS_Mesh.hxx"
31 #include "SMESH_Algo.hxx"
32 #include "SMESH_Group.hxx"
33 #include "SMESH_Mesh.hxx"
34 #include "SMESH_MeshAlgos.hxx"
35 #include "SMESH_MesherHelper.hxx"
36 #include "SMESH_subMesh.hxx"
37
38 #include <IntAna_IntConicQuad.hxx>
39 #include <IntAna_Quadric.hxx>
40 #include <TColgp_Array1OfPnt.hxx>
41 #include <TColgp_Array1OfVec.hxx>
42 #include <TColgp_SequenceOfPnt.hxx>
43 #include <TopExp_Explorer.hxx>
44 #include <TopoDS.hxx>
45 #include <TopoDS_Iterator.hxx>
46 #include <gp_Lin.hxx>
47 #include <gp_Pln.hxx>
48
49 #include "utilities.h"
50
51 #include <string>
52 #include <numeric>
53 #include <limits>
54
55 using namespace std;
56
57 enum EQuadNature { NOT_QUAD, QUAD, DEGEN_QUAD, PYRAM_APEX = 4, TRIA_APEX = 0 };
58
59 // std-like iterator used to get coordinates of nodes of mesh element
60 typedef SMDS_StdIterator< SMESH_TNodeXYZ, SMDS_ElemIteratorPtr > TXyzIterator;
61
62 //================================================================================
63 /*!
64  * \brief Return ID of pyramid base face, for debug
65  */
66 //================================================================================
67
68 int getFaceID(const SMDS_MeshElement* pyram)
69 {
70   if ( pyram )
71     if ( const SMDS_MeshElement* f = SMDS_Mesh::FindFace( pyram->GetNode(0),
72                                                           pyram->GetNode(1),
73                                                           pyram->GetNode(2),
74                                                           pyram->GetNode(3)))
75       return f->GetID();
76   return -1;
77 }
78
79 namespace
80 {
81   //================================================================================
82   /*!
83    * \brief Return true if two nodes of triangles are equal
84    */
85   //================================================================================
86
87   bool EqualTriangles(const SMDS_MeshElement* F1,const SMDS_MeshElement* F2)
88   {
89     return
90       ( F1->GetNode(1)==F2->GetNode(2) && F1->GetNode(2)==F2->GetNode(1) ) ||
91       ( F1->GetNode(1)==F2->GetNode(1) && F1->GetNode(2)==F2->GetNode(2) );
92   }
93   //================================================================================
94   /*!
95    * \brief Return true if two adjacent pyramids are too close one to another
96    * so that a tetrahedron to built between them would have too poor quality
97    */
98   //================================================================================
99
100   bool TooCloseAdjacent( const SMDS_MeshElement* PrmI,
101                          const SMDS_MeshElement* PrmJ,
102                          const bool              hasShape)
103   {
104     const SMDS_MeshNode* nApexI = PrmI->GetNode(4);
105     const SMDS_MeshNode* nApexJ = PrmJ->GetNode(4);
106     if ( nApexI == nApexJ ||
107          nApexI->getshapeId() != nApexJ->getshapeId() )
108       return false;
109
110     // Find two common base nodes and their indices within PrmI and PrmJ
111     const SMDS_MeshNode* baseNodes[2] = { 0,0 };
112     int baseNodesIndI[2], baseNodesIndJ[2];
113     for ( int i = 0; i < 4 ; ++i )
114     {
115       int j = PrmJ->GetNodeIndex( PrmI->GetNode(i));
116       if ( j >= 0 )
117       {
118         int ind = baseNodes[0] ? 1:0;
119         if ( baseNodes[ ind ])
120           return false; // pyramids with a common base face
121         baseNodes    [ ind ] = PrmI->GetNode(i);
122         baseNodesIndI[ ind ] = i;
123         baseNodesIndJ[ ind ] = j;
124       }
125     }
126     if ( !baseNodes[1] ) return false; // not adjacent
127
128     // Get normals of triangles sharing baseNodes
129     gp_XYZ apexI = SMESH_TNodeXYZ( nApexI );
130     gp_XYZ apexJ = SMESH_TNodeXYZ( nApexJ );
131     gp_XYZ base1 = SMESH_TNodeXYZ( baseNodes[0]);
132     gp_XYZ base2 = SMESH_TNodeXYZ( baseNodes[1]);
133     gp_Vec baseVec( base1, base2 );
134     gp_Vec baI( base1, apexI );
135     gp_Vec baJ( base1, apexJ );
136     gp_Vec nI = baseVec.Crossed( baI );
137     gp_Vec nJ = baseVec.Crossed( baJ );
138
139     // Check angle between normals
140     double  angle = nI.Angle( nJ );
141     bool tooClose = ( angle < 15. * M_PI / 180. );
142
143     // Check if pyramids collide
144     if ( !tooClose && ( baI * baJ > 0 ) && ( nI * nJ > 0 ))
145     {
146       // find out if nI points outside of PrmI or inside
147       int    dInd = baseNodesIndI[1] - baseNodesIndI[0];
148       bool isOutI = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0;
149
150       // find out sign of projection of baI to nJ
151       double proj = baI * nJ;
152
153       tooClose = ( isOutI ? proj > 0 : proj < 0 );
154     }
155
156     // Check if PrmI and PrmJ are in same domain
157     if ( tooClose && !hasShape )
158     {
159       // check order of baseNodes within pyramids, it must be opposite
160       int dInd;
161       dInd = baseNodesIndI[1] - baseNodesIndI[0];
162       bool isOutI = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0;
163       dInd = baseNodesIndJ[1] - baseNodesIndJ[0];
164       bool isOutJ = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0;
165       if ( isOutJ == isOutI )
166         return false; // other domain
167
168       // direct both normals outside pyramid
169       ( isOutI ? nJ : nI ).Reverse();
170
171       // check absence of a face separating domains between pyramids
172       TIDSortedElemSet emptySet, avoidSet;
173       int i1, i2;
174       while ( const SMDS_MeshElement* f =
175               SMESH_MeshAlgos::FindFaceInSet( baseNodes[0], baseNodes[1],
176                                               emptySet, avoidSet, &i1, &i2 ))
177       {
178         avoidSet.insert( f );
179
180         // face node other than baseNodes
181         int otherNodeInd = 0;
182         while ( otherNodeInd == i1 || otherNodeInd == i2 ) otherNodeInd++;
183         const SMDS_MeshNode* otherFaceNode = f->GetNode( otherNodeInd );
184
185         if ( otherFaceNode == nApexI || otherFaceNode == nApexJ )
186           continue; // f is a temporary triangle
187
188         // check if f is a base face of either of pyramids
189         if ( f->NbCornerNodes() == 4 &&
190              ( PrmI->GetNodeIndex( otherFaceNode ) >= 0 ||
191                PrmJ->GetNodeIndex( otherFaceNode ) >= 0 ))
192           continue; // f is a base quadrangle
193
194         // check projections of face direction (baOFN) to triangle normals (nI and nJ)
195         gp_Vec baOFN( base2, SMESH_TNodeXYZ( otherFaceNode ));
196         if ( nI * baOFN > 0 && nJ * baOFN > 0 &&
197              baI* baOFN > 0 && baJ* baOFN > 0 ) // issue 0023212
198         {
199           tooClose = false; // f is between pyramids
200           break;
201         }
202       }
203     }
204
205     return tooClose;
206   }
207
208   //================================================================================
209   /*!
210    * \brief Move medium nodes of merged quadratic pyramids
211    */
212   //================================================================================
213
214   void UpdateQuadraticPyramids(const set<const SMDS_MeshNode*>& commonApex,
215                                SMESHDS_Mesh*                    meshDS)
216   {
217     typedef SMDS_StdIterator< const SMDS_MeshElement*, SMDS_ElemIteratorPtr > TStdElemIterator;
218     TStdElemIterator itEnd;
219
220     // shift of node index to get medium nodes between the 4 base nodes and the apex
221     const int base2MediumShift = 9;
222
223     set<const SMDS_MeshNode*>::const_iterator nIt = commonApex.begin();
224     for ( ; nIt != commonApex.end(); ++nIt )
225     {
226       SMESH_TNodeXYZ apex( *nIt );
227
228       vector< const SMDS_MeshElement* > pyrams // pyramids sharing the apex node
229         ( TStdElemIterator( apex._node->GetInverseElementIterator( SMDSAbs_Volume )), itEnd );
230
231       // Select medium nodes to keep and medium nodes to remove
232
233       typedef map < const SMDS_MeshNode*, const SMDS_MeshNode*, TIDCompare > TN2NMap;
234       TN2NMap base2medium; // to keep
235       vector< const SMDS_MeshNode* > nodesToRemove;
236
237       for ( unsigned i = 0; i < pyrams.size(); ++i )
238         for ( int baseIndex = 0; baseIndex < PYRAM_APEX; ++baseIndex )
239         {
240           SMESH_TNodeXYZ         base = pyrams[i]->GetNode( baseIndex );
241           const SMDS_MeshNode* medium = pyrams[i]->GetNode( baseIndex + base2MediumShift );
242           TN2NMap::iterator b2m = base2medium.insert( make_pair( base._node, medium )).first;
243           if ( b2m->second != medium )
244           {
245             nodesToRemove.push_back( medium );
246           }
247           else
248           {
249             // move the kept medium node
250             gp_XYZ newXYZ = 0.5 * ( apex + base );
251             meshDS->MoveNode( medium, newXYZ.X(), newXYZ.Y(), newXYZ.Z() );
252           }
253         }
254
255       // Within pyramids, replace nodes to remove by nodes to keep
256
257       for ( unsigned i = 0; i < pyrams.size(); ++i )
258       {
259         vector< const SMDS_MeshNode* > nodes( pyrams[i]->begin_nodes(),
260                                               pyrams[i]->end_nodes() );
261         for ( int baseIndex = 0; baseIndex < PYRAM_APEX; ++baseIndex )
262         {
263           const SMDS_MeshNode* base = pyrams[i]->GetNode( baseIndex );
264           nodes[ baseIndex + base2MediumShift ] = base2medium[ base ];
265         }
266         meshDS->ChangeElementNodes( pyrams[i], &nodes[0], nodes.size());
267       }
268
269       // Remove the replaced nodes
270
271       if ( !nodesToRemove.empty() )
272       {
273         SMESHDS_SubMesh * sm = meshDS->MeshElements( nodesToRemove[0]->getshapeId() );
274         for ( unsigned i = 0; i < nodesToRemove.size(); ++i )
275           meshDS->RemoveFreeNode( nodesToRemove[i], sm, /*fromGroups=*/false);
276       }
277     }
278     return;
279   }
280
281   //================================================================================
282   /*!
283    * \brief Store an error about overlapping faces
284    */
285   //================================================================================
286
287   bool overlapError( SMESH_Mesh&             mesh,
288                      const SMDS_MeshElement* face1,
289                      const SMDS_MeshElement* face2,
290                      const TopoDS_Shape&     shape = TopoDS_Shape())
291   {
292     if ( !face1 || !face2 ) return false;
293
294     SMESH_Comment msg;
295     msg << "face " << face1->GetID() << " overlaps face " << face2->GetID();
296
297     SMESH_subMesh * sm = 0;
298     if ( shape.IsNull() )
299     {
300       sm = mesh.GetSubMesh( mesh.GetShapeToMesh() );
301     }
302     else if ( shape.ShapeType() >= TopAbs_SOLID )
303     {
304       sm = mesh.GetSubMesh( shape );
305     }
306     else
307     {
308       TopoDS_Iterator it ( shape );
309       if ( it.More() )
310         sm = mesh.GetSubMesh( it.Value() );
311     }
312     if ( sm )
313     {
314       SMESH_ComputeErrorPtr& err = sm->GetComputeError();
315       if ( !err || err->IsOK() )
316       {
317         SMESH_BadInputElements* badElems =
318           new SMESH_BadInputElements( mesh.GetMeshDS(),COMPERR_BAD_INPUT_MESH, msg, sm->GetAlgo() );
319         badElems->add( face1 );
320         badElems->add( face2 );
321         err.reset( badElems );
322       }
323     }
324
325     return false; // == "algo fails"
326   }
327
328   //================================================================================
329   /*!
330    * \brief Check if a face is in a SOLID
331    */
332   //================================================================================
333
334   bool isInSolid( vector<const SMDS_MeshNode*> & faceNodes,
335                   const int                      nbNodes,
336                   const int                      solidID )
337   {
338     if ( !faceNodes[0] )
339       return true; // NOT_QUAD
340     for ( int i = 0; i < nbNodes; ++i )
341     {
342       int shapeID = faceNodes[i]->GetShapeID();
343       if ( shapeID == solidID )
344         return true;
345     }
346     faceNodes.resize( nbNodes );
347     std::vector<const SMDS_MeshElement*> vols;
348     SMDS_Mesh::GetElementsByNodes( faceNodes, vols, SMDSAbs_Volume );
349     bool inSolid = false;
350     for ( size_t i = 0; i < vols.size() && !inSolid; ++i )
351     {
352       int shapeID = vols[i]->GetShapeID();
353       inSolid = ( shapeID == solidID );
354     }
355     faceNodes.push_back( faceNodes[0] );
356     return inSolid;
357   }
358 }
359
360 //================================================================================
361 /*!
362  * \brief Merge the two pyramids (i.e. fuse their apex) and others already merged with them
363  */
364 //================================================================================
365
366 void StdMeshers_QuadToTriaAdaptor::MergePiramids( const SMDS_MeshElement*     PrmI,
367                                                   const SMDS_MeshElement*     PrmJ,
368                                                   set<const SMDS_MeshNode*> & nodesToMove)
369 {
370   const SMDS_MeshNode* Nrem = PrmJ->GetNode(4); // node to remove
371   //int nbJ = Nrem->NbInverseElements( SMDSAbs_Volume );
372   SMESH_TNodeXYZ Pj( Nrem );
373
374   // an apex node to make common to all merged pyramids
375   SMDS_MeshNode* CommonNode = const_cast<SMDS_MeshNode*>(PrmI->GetNode(4));
376   if ( CommonNode == Nrem ) return; // already merged
377   //int nbI = CommonNode->NbInverseElements( SMDSAbs_Volume );
378   SMESH_TNodeXYZ Pi( CommonNode );
379   gp_XYZ Pnew = /*( nbI*Pi + nbJ*Pj ) / (nbI+nbJ);*/ 0.5 * ( Pi + Pj );
380   CommonNode->setXYZ( Pnew.X(), Pnew.Y(), Pnew.Z() );
381
382   nodesToMove.insert( CommonNode );
383   nodesToMove.erase ( Nrem );
384
385   //cout << "MergePiramids F" << getFaceID( PrmI ) << " - F" << getFaceID( PrmJ ) << endl;
386
387   typedef SMDS_StdIterator< const SMDS_MeshElement*, SMDS_ElemIteratorPtr > TStdElemIterator;
388   TStdElemIterator itEnd;
389
390   typedef std::map< const SMDS_MeshNode*, const SMDS_MeshNode* > TNNMap;
391   TNNMap mediumReplaceMap;
392
393   // find and remove coincided faces of merged pyramids
394   vector< const SMDS_MeshElement* > inverseElems
395     // copy inverse elements to avoid iteration on changing container
396     ( TStdElemIterator( CommonNode->GetInverseElementIterator(SMDSAbs_Face)), itEnd);
397   for ( size_t i = 0; i < inverseElems.size(); ++i )
398   {
399     const SMDS_MeshElement* FI = inverseElems[i];
400     const SMDS_MeshElement* FJEqual = 0;
401     SMDS_ElemIteratorPtr triItJ = Nrem->GetInverseElementIterator(SMDSAbs_Face);
402     while ( !FJEqual && triItJ->more() )
403     {
404       const SMDS_MeshElement* FJ = triItJ->next();
405       if ( EqualTriangles( FJ, FI ))
406         FJEqual = FJ;
407     }
408     if ( FJEqual )
409     {
410       if ( FJEqual->NbNodes() == 6 ) // find medium nodes to replace
411       {
412         mediumReplaceMap.insert( std::make_pair( FJEqual->GetNode(3), FI->GetNode(5) ));
413         mediumReplaceMap.insert( std::make_pair( FJEqual->GetNode(5), FI->GetNode(3) ));
414       }
415       removeTmpElement( FI );
416       removeTmpElement( FJEqual );
417       myRemovedTrias.insert( FI );
418       myRemovedTrias.insert( FJEqual );
419     }
420   }
421
422   // set the common apex node to pyramids and triangles merged with J
423   vector< const SMDS_MeshNode* > nodes;
424   inverseElems.assign( TStdElemIterator( Nrem->GetInverseElementIterator()), itEnd );
425   for ( size_t i = 0; i < inverseElems.size(); ++i )
426   {
427     const SMDS_MeshElement* elem = inverseElems[i];
428     nodes.assign( elem->begin_nodes(), elem->end_nodes() );
429     nodes[ elem->GetType() == SMDSAbs_Volume ? PYRAM_APEX : TRIA_APEX ] = CommonNode;
430     if ( !mediumReplaceMap.empty() )
431       for ( size_t iN = elem->NbCornerNodes(); iN < nodes.size(); ++iN )
432       {
433         TNNMap::iterator n2n = mediumReplaceMap.find( nodes[iN] );
434         if ( n2n != mediumReplaceMap.end() )
435           nodes[iN] = n2n->second;
436       }
437     GetMeshDS()->ChangeElementNodes( elem, &nodes[0], nodes.size());
438   }
439   ASSERT( Nrem->NbInverseElements() == 0 );
440   GetMeshDS()->RemoveFreeNode( Nrem,
441                                GetMeshDS()->MeshElements( Nrem->getshapeId()),
442                                /*fromGroups=*/false);
443   if ( !mediumReplaceMap.empty() )
444     for ( TNNMap::iterator n2n = mediumReplaceMap.begin(); n2n != mediumReplaceMap.end(); ++n2n )
445     {
446       const SMDS_MeshNode* remNode = n2n->first;
447       if ( !remNode->IsNull() && remNode->NbInverseElements() == 0 )
448         GetMeshDS()->RemoveFreeNode( remNode, 0, /*fromGroups=*/false);
449     }
450   return;
451 }
452
453 //================================================================================
454 /*!
455  * \brief Merges adjacent pyramids
456  */
457 //================================================================================
458
459 void StdMeshers_QuadToTriaAdaptor::MergeAdjacent(const SMDS_MeshElement*    PrmI,
460                                                  set<const SMDS_MeshNode*>& nodesToMove,
461                                                  const bool                 isRecursion)
462 {
463   TIDSortedElemSet adjacentPyrams;
464   bool mergedPyrams = false;
465   for ( int k = 0; k < 4; k++ ) // loop on 4 base nodes of PrmI
466   {
467     const SMDS_MeshNode*   n = PrmI->GetNode(k);
468     SMDS_ElemIteratorPtr vIt = n->GetInverseElementIterator( SMDSAbs_Volume );
469     while ( vIt->more() )
470     {
471       const SMDS_MeshElement* PrmJ = vIt->next();
472       if ( PrmJ == PrmI || PrmJ->NbCornerNodes() != 5 || !adjacentPyrams.insert( PrmJ ).second  )
473         continue;
474       if ( TooCloseAdjacent( PrmI, PrmJ, GetMesh()->HasShapeToMesh() ))
475       {
476         MergePiramids( PrmI, PrmJ, nodesToMove );
477         mergedPyrams = true;
478         // container of inverse elements can change
479         // vIt = n->GetInverseElementIterator( SMDSAbs_Volume ); -- iterator re-implemented
480       }
481     }
482   }
483   if ( mergedPyrams && !isRecursion )
484   {
485     TIDSortedElemSet::iterator prm;
486     for (prm = adjacentPyrams.begin(); prm != adjacentPyrams.end(); ++prm)
487       MergeAdjacent( *prm, nodesToMove, true );
488   }
489   return;
490 }
491
492 //================================================================================
493 /*!
494  * \brief Decrease height of a given or adjacent pyramids if height difference
495  *        is too large
496  *  \param [in] pyram - a pyramid to treat
497  *  \param [inout] h2 - pyramid's square height
498  *  \return bool - true if the height changes
499  */
500 //================================================================================
501
502 bool StdMeshers_QuadToTriaAdaptor::DecreaseHeightDifference( const SMDS_MeshElement* thePyram,
503                                                              const double            theH2 )
504 {
505   const double allowedFactor2 = 2. * 2.;
506
507   bool modif = false;
508   myNodes[0] = thePyram->GetNode( 3 );
509   for ( int i = 0; i < 4; ++i )
510   {
511     myNodes[1] = thePyram->GetNode( i );
512     SMDS_Mesh::GetElementsByNodes( myNodes, myAdjPyrams, SMDSAbs_Volume );
513     myNodes[0] = myNodes[1];
514
515     for ( const SMDS_MeshElement* pyramAdj : myAdjPyrams )
516     {
517       if ( pyramAdj == thePyram )
518         continue;
519       if ( !myPyramHeight2.IsBound( pyramAdj ))
520         continue;
521       double h2Adj = Abs( myPyramHeight2( pyramAdj ));
522       double h2    = Abs( theH2 );
523       if ( h2Adj > h2 )
524       {
525         if ( h2 * allowedFactor2 < h2Adj )
526         {
527           // bind negative value to allow finding pyramids whose height must change
528           myPyramHeight2.Bind( pyramAdj, - h2 * allowedFactor2 );
529           modif = true;
530         }
531       }
532       else
533       {
534         if ( h2Adj * allowedFactor2 < h2 )
535         {
536           // bind negative value to allow finding pyramids whose height must change
537           myPyramHeight2.Bind( thePyram, - h2Adj * allowedFactor2 );
538           modif = true;
539         }
540       }
541     }
542   }
543   return modif;
544 }
545
546
547 //================================================================================
548 /*!
549  * \brief Constructor
550  */
551 //================================================================================
552
553 StdMeshers_QuadToTriaAdaptor::StdMeshers_QuadToTriaAdaptor():
554   myElemSearcher(0)
555 {
556 }
557
558 //================================================================================
559 /*!
560  * \brief Destructor
561  */
562 //================================================================================
563
564 StdMeshers_QuadToTriaAdaptor::~StdMeshers_QuadToTriaAdaptor()
565 {
566   // temporary faces are deleted by ~SMESH_ProxyMesh()
567   if ( myElemSearcher ) delete myElemSearcher;
568   myElemSearcher=0;
569 }
570
571 //=======================================================================
572 //function : FindBestPoint
573 //purpose  : Return a point P laying on the line (PC,V) so that triangle
574 //           (P, P1, P2) to be equilateral as much as possible
575 //           V - normal to (P1,P2,PC)
576 //=======================================================================
577
578 static gp_Pnt FindBestPoint(const gp_Pnt& P1, const gp_Pnt& P2,
579                             const gp_Pnt& PC, const gp_Vec& V,
580                             double & shift)
581 {
582   gp_Pnt Pbest = PC;
583   shift = 0;
584   const double a2 = P1.SquareDistance(P2);
585   const double b2 = P1.SquareDistance(PC);
586   const double c2 = P2.SquareDistance(PC);
587   if ( a2 < ( b2 + Sqrt( 4 * b2 * c2 ) + c2 ) / 4 ) // ( a < (b+c)/2 )
588     return Pbest;
589   else {
590     // find shift along V in order a to became equal to (b+c)/2
591     const double Vsize = V.Magnitude();
592     if ( fabs( Vsize ) > std::numeric_limits<double>::min() )
593     {
594       shift = sqrt( a2 + (b2-c2)*(b2-c2)/16/a2 - (b2+c2)/2 );
595       Pbest.ChangeCoord() += shift * V.XYZ() / Vsize;
596     }
597   }
598   return Pbest;
599 }
600
601 //=======================================================================
602 //function : HasIntersection3
603 //purpose  : Find intersection point between a triangle (P1,P2,P3)
604 //           and a segment [PC,P]
605 //=======================================================================
606
607 static bool HasIntersection3(const gp_Pnt& P,  const gp_Pnt& PC, gp_Pnt&       Pint,
608                              const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3)
609 {
610   const double EPSILON = 1e-6;
611   double segLen = P.Distance( PC );
612
613   gp_XYZ  orig = PC.XYZ();
614   gp_XYZ   dir = ( P.XYZ() - PC.XYZ() ) / segLen;
615   gp_XYZ vert0 = P1.XYZ();
616   gp_XYZ vert1 = P2.XYZ();
617   gp_XYZ vert2 = P3.XYZ();
618
619   gp_XYZ edge1 = vert1 - vert0;
620   gp_XYZ edge2 = vert2 - vert0;
621
622   /* begin calculating determinant - also used to calculate U parameter */
623   gp_XYZ pvec = dir ^ edge2;
624
625   /* if determinant is near zero, ray lies in plane of triangle */
626   double det = edge1 * pvec;
627
628   const double ANGL_EPSILON = 1e-12;
629   if ( det > -ANGL_EPSILON && det < ANGL_EPSILON )
630     return false;
631
632   /* calculate distance from vert0 to ray origin */
633   gp_XYZ  tvec = orig - vert0;
634
635   /* calculate U parameter and test bounds */
636   double u = ( tvec * pvec ) / det;
637   //if (u < 0.0 || u > 1.0)
638   if (u < -EPSILON || u > 1.0 + EPSILON)
639     return false;
640
641   /* prepare to test V parameter */
642   gp_XYZ qvec = tvec ^ edge1;
643
644   /* calculate V parameter and test bounds */
645   double v = (dir * qvec) / det;
646   //if ( v < 0.0 || u + v > 1.0 )
647   if ( v < -EPSILON || u + v > 1.0 + EPSILON)
648     return false;
649
650   /* calculate t, ray intersects triangle */
651   double t = (edge2 * qvec) / det;
652
653   Pint = orig + dir * t;
654
655   bool hasInt = ( t > 0.  &&  t < segLen );
656
657   if ( hasInt && det < EPSILON ) // t is inaccurate, additionally check
658   {
659     gp_XYZ triNorm = edge1 ^ edge2;
660     gp_XYZ int0vec = Pint.XYZ() - vert0;
661     gp_XYZ in      = triNorm ^ edge1; // dir inside triangle from edge1
662     double dot     = int0vec * in;
663     if ( dot < 0 && dot / triNorm.Modulus() < -EPSILON )
664       return false;
665     in  = edge2 ^ triNorm;
666     dot = int0vec * in;
667     if ( dot < 0 && dot / triNorm.Modulus() < -EPSILON )
668       return false;
669     gp_XYZ int1vec = Pint.XYZ() - vert1;
670     in  = triNorm ^ ( vert2 - vert1 );
671     dot = int1vec * in;
672     if ( dot < 0 && dot / triNorm.Modulus() < -EPSILON )
673       return false;
674   }
675   return hasInt;
676 }
677
678 //=======================================================================
679 //function : HasIntersection
680 //purpose  : Auxiliary for CheckIntersection()
681 //=======================================================================
682
683 static bool HasIntersection(const gp_Pnt& P, const gp_Pnt& PC, gp_Pnt& Pint,
684                             TColgp_SequenceOfPnt& aContour)
685 {
686   if ( aContour.Length() == 3 ) {
687     return HasIntersection3( P, PC, Pint, aContour(1), aContour(2), aContour(3) );
688   }
689   else {
690     bool check = false;
691     if( (aContour(1).SquareDistance(aContour(2)) > 1.e-12) &&
692         (aContour(1).SquareDistance(aContour(3)) > 1.e-12) &&
693         (aContour(2).SquareDistance(aContour(3)) > 1.e-12) ) {
694       check = HasIntersection3( P, PC, Pint, aContour(1), aContour(2), aContour(3) );
695     }
696     if(check) return true;
697     if( (aContour(1).SquareDistance(aContour(4)) > 1.e-12) &&
698         (aContour(1).SquareDistance(aContour(3)) > 1.e-12) &&
699         (aContour(4).SquareDistance(aContour(3)) > 1.e-12) ) {
700       check = HasIntersection3( P, PC, Pint, aContour(1), aContour(3), aContour(4) );
701     }
702     if(check) return true;
703   }
704
705   return false;
706 }
707
708 //================================================================================
709 /*!
710  * \brief Return allowed height of a pyramid
711  *  \param Papex - optimal pyramid apex
712  *  \param PC - gravity center of a quadrangle
713  *  \param PN - four nodes of the quadrangle
714  *  \param aMesh - mesh
715  *  \param NotCheckedFace - the quadrangle face
716  *  \param Shape - the shape being meshed
717  *  \retval false if mesh invalidity detected
718  */
719 //================================================================================
720
721 bool StdMeshers_QuadToTriaAdaptor::LimitHeight (gp_Pnt&                             Papex,
722                                                 const gp_Pnt&                       PC,
723                                                 const TColgp_Array1OfPnt&           PN,
724                                                 const vector<const SMDS_MeshNode*>& FNodes,
725                                                 SMESH_Mesh&                         aMesh,
726                                                 const SMDS_MeshElement*             NotCheckedFace,
727                                                 const bool                          UseApexRay,
728                                                 const TopoDS_Shape&                 Shape)
729 {
730   if ( !myElemSearcher )
731     myElemSearcher = SMESH_MeshAlgos::GetElementSearcher( *aMesh.GetMeshDS() );
732   SMESH_ElementSearcher* searcher = const_cast<SMESH_ElementSearcher*>(myElemSearcher);
733
734   // Find intersection of faces with (P,PC) segment elongated 3 times
735
736   double height = Papex.Distance( PC );
737   gp_Ax1 line( PC, gp_Vec( PC, Papex ));
738   gp_Pnt Pint, Ptest;
739   vector< const SMDS_MeshElement* > suspectFaces;
740   TColgp_SequenceOfPnt aContour;
741
742   if ( UseApexRay )
743   {
744     double idealHeight = height;
745     const SMDS_MeshElement* intFace = 0;
746
747     // find intersection closest to PC
748     Ptest = PC.XYZ() + line.Direction().XYZ() * height * 3;
749
750     searcher->GetElementsNearLine( line, SMDSAbs_Face, suspectFaces );
751     for ( size_t iF = 0; iF < suspectFaces.size(); ++iF )
752     {
753       const SMDS_MeshElement* face = suspectFaces[iF];
754       if ( face == NotCheckedFace ) continue;
755
756       aContour.Clear();
757       for ( int i = 0, nb = face->NbCornerNodes(); i < nb; ++i )
758         aContour.Append( SMESH_TNodeXYZ( face->GetNode(i) ));
759
760       if ( HasIntersection( Ptest, PC, Pint, aContour ))
761       {
762         double dInt = PC.Distance( Pint ) / 3.;
763         if ( dInt < height )
764         {
765           height = dInt;
766           intFace = face;
767         }
768       }
769     }
770     if ( height < 1e-2 * idealHeight && intFace )
771       return overlapError( aMesh, NotCheckedFace, intFace, Shape );
772   }
773
774   // Find faces intersecting triangular facets of the pyramid (issue 23212)
775
776   gp_XYZ center   = PC.XYZ() + line.Direction().XYZ() * height * 0.5;
777   double diameter = Max( PN(1).Distance(PN(3)), PN(2).Distance(PN(4)));
778   suspectFaces.clear();
779   searcher->GetElementsInSphere( center, diameter * 0.6, SMDSAbs_Face, suspectFaces);
780
781   const double upShift = 1.5;
782   Ptest = PC.XYZ() + line.Direction().XYZ() * height * upShift; // tmp apex
783
784   for ( size_t iF = 0; iF < suspectFaces.size(); ++iF )
785   {
786     const SMDS_MeshElement* face = suspectFaces[iF];
787     if ( face == NotCheckedFace ) continue;
788     if ( face->GetNodeIndex( FNodes[0] ) >= 0 ||
789          face->GetNodeIndex( FNodes[1] ) >= 0 ||
790          face->GetNodeIndex( FNodes[2] ) >= 0 ||
791          face->GetNodeIndex( FNodes[3] ) >= 0 )
792       continue; // neighbor face of the quadrangle
793
794     // limit height using points of intersection of face links with pyramid facets
795     int   nbN = face->NbCornerNodes();
796     gp_Pnt P1 = SMESH_TNodeXYZ( face->GetNode( nbN-1 )); // 1st link end
797     for ( int i = 0; i < nbN; ++i )
798     {
799       gp_Pnt P2 = SMESH_TNodeXYZ( face->GetNode(i) );    // 2nd link end
800
801       for ( int iN = 1; iN <= 4; ++iN ) // loop on pyramid facets
802       {
803         if ( HasIntersection3( P1, P2, Pint, PN(iN), PN(iN+1), Ptest ))
804         {
805           height = Min( height, gp_Vec( PC, Pint ) * line.Direction() );
806           //Ptest = PC.XYZ() + line.Direction().XYZ() * height * upShift; // new tmp apex
807         }
808       }
809       P1 = P2;
810     }
811   }
812
813   Papex  = PC.XYZ() + line.Direction().XYZ() * height;
814
815   return true;
816 }
817
818 //================================================================================
819 /*!
820  * \brief Retrieve data of the given face
821  *  \param PN - coordinates of face nodes
822  *  \param VN - cross products of vectors (PC-PN(i)) ^ (PC-PN(i+1))
823  *  \param FNodes - face nodes
824  *  \param PC - gravity center of nodes
825  *  \param VNorm - face normal (sum of VN)
826  *  \param volumes - two volumes sharing the given face, the first is in VNorm direction
827  *  \retval int - 0 if given face is not quad,
828  *                1 if given face is quad,
829  *                2 if given face is degenerate quad (two nodes are coincided)
830  */
831 //================================================================================
832
833 int StdMeshers_QuadToTriaAdaptor::Preparation(const SMDS_MeshElement*       face,
834                                               TColgp_Array1OfPnt&           PN,
835                                               TColgp_Array1OfVec&           VN,
836                                               vector<const SMDS_MeshNode*>& FNodes,
837                                               gp_Pnt&                       PC,
838                                               gp_Vec&                       VNorm,
839                                               const SMDS_MeshElement**      volumes)
840 {
841   if( face->NbCornerNodes() != 4 )
842   {
843     return NOT_QUAD;
844   }
845
846   int i = 0;
847   gp_XYZ xyzC(0., 0., 0.);
848   for ( i = 0; i < 4; ++i )
849   {
850     gp_XYZ p = SMESH_TNodeXYZ( FNodes[i] = face->GetNode(i) );
851     PN.SetValue( i+1, p );
852     xyzC += p;
853   }
854   PC = xyzC/4;
855
856   int nbp = 4;
857
858   int j = 0;
859   for ( i = 1; i < 4; i++ )
860   {
861     j = i+1;
862     for(; j<=4; j++) {
863       if( PN(i).Distance(PN(j)) < 1.e-6 )
864         break;
865     }
866     if(j<=4) break;
867   }
868
869   bool hasdeg = false;
870   if ( i < 4 )
871   {
872     hasdeg = true;
873     gp_Pnt Pdeg = PN(i);
874
875     list< const SMDS_MeshNode* >::iterator itdg = myDegNodes.begin();
876     const SMDS_MeshNode* DegNode = 0;
877     for(; itdg!=myDegNodes.end(); itdg++) {
878       const SMDS_MeshNode* N = (*itdg);
879       gp_Pnt Ptmp(N->X(),N->Y(),N->Z());
880       if(Pdeg.Distance(Ptmp)<1.e-6) {
881         DegNode = N;
882         break;
883       }
884     }
885     if(!DegNode) {
886       DegNode = FNodes[i-1];
887       myDegNodes.push_back(DegNode);
888     }
889     else {
890       FNodes[i-1] = DegNode;
891     }
892     for(i=j; i<4; i++) {
893       PN.SetValue(i,PN.Value(i+1));
894       FNodes[i-1] = FNodes[i];
895     }
896     nbp = 3;
897   }
898
899   PN.SetValue(nbp+1,PN(1));
900   FNodes[nbp] = FNodes[0];
901
902   // find normal direction
903   gp_Vec V1(PC,PN(nbp));
904   gp_Vec V2(PC,PN(1));
905   VNorm = V1.Crossed(V2);
906   VN.SetValue(nbp,VNorm);
907   for(i=1; i<nbp; i++) {
908     V1 = gp_Vec(PC,PN(i));
909     V2 = gp_Vec(PC,PN(i+1));
910     gp_Vec Vtmp = V1.Crossed(V2);
911     VN.SetValue(i,Vtmp);
912     VNorm += Vtmp;
913   }
914
915   // find volumes sharing the face
916   if ( volumes )
917   {
918     volumes[0] = volumes[1] = 0;
919     SMDS_ElemIteratorPtr vIt = FNodes[0]->GetInverseElementIterator( SMDSAbs_Volume );
920     while ( vIt->more() )
921     {
922       const SMDS_MeshElement* vol = vIt->next();
923       bool volSharesAllNodes = true;
924       for ( int i = 1; i < face->NbNodes() && volSharesAllNodes; ++i )
925         volSharesAllNodes = ( vol->GetNodeIndex( FNodes[i] ) >= 0 );
926       if ( volSharesAllNodes )
927         volumes[ volumes[0] ? 1 : 0 ] = vol;
928       // we could additionally check that vol has all FNodes in its one face using SMDS_VolumeTool
929     }
930     // define volume position relating to the face normal
931     if ( volumes[0] )
932     {
933       // get volume gc
934       SMDS_ElemIteratorPtr nodeIt = volumes[0]->nodesIterator();
935       gp_XYZ volGC(0,0,0);
936       volGC = accumulate( TXyzIterator(nodeIt), TXyzIterator(), volGC ) / volumes[0]->NbNodes();
937
938       if ( VNorm * gp_Vec( PC, volGC ) < 0 )
939         swap( volumes[0], volumes[1] );
940     }
941   }
942
943   return hasdeg ? DEGEN_QUAD : QUAD;
944 }
945
946
947 //=======================================================================
948 //function : Compute
949 //purpose  :
950 //=======================================================================
951
952 bool StdMeshers_QuadToTriaAdaptor::Compute(SMESH_Mesh&         aMesh,
953                                            const TopoDS_Shape& aShape,
954                                            SMESH_ProxyMesh*    aProxyMesh)
955 {
956   SMESH_ProxyMesh::setMesh( aMesh );
957
958   if ( aShape.ShapeType() != TopAbs_SOLID )
959     return false;
960
961   myShape = aShape;
962
963   vector<const SMDS_MeshElement*> myPyramids;
964
965   const SMESHDS_SubMesh * aSubMeshDSFace;
966   SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
967   SMESH_MesherHelper helper1(aMesh);
968   helper1.IsQuadraticSubMesh(aShape);
969
970   if ( myElemSearcher ) delete myElemSearcher;
971   vector< SMDS_ElemIteratorPtr > itVec;
972   if ( aProxyMesh )
973   {
974     itVec.push_back( aProxyMesh->GetFaces( aShape ));
975   }
976   else
977   {
978     for ( TopExp_Explorer exp(aShape,TopAbs_FACE); exp.More(); exp.Next() )
979       if (( aSubMeshDSFace = meshDS->MeshElements( exp.Current() )))
980         itVec.push_back( aSubMeshDSFace->GetElements() );
981   }
982   typedef
983     SMDS_IteratorOnIterators< const SMDS_MeshElement*, vector< SMDS_ElemIteratorPtr > > TIter;
984   SMDS_ElemIteratorPtr faceIt( new TIter( itVec ));
985   myElemSearcher = SMESH_MeshAlgos::GetElementSearcher( *meshDS, faceIt );
986
987   TColgp_Array1OfPnt PN(1,5);
988   TColgp_Array1OfVec VN(1,4);
989   vector<const SMDS_MeshNode*> FNodes(5);
990   gp_Pnt PC;
991   gp_Vec VNorm;
992   const int solidID = meshDS->ShapeToIndex( aShape );
993
994   for ( TopExp_Explorer exp(aShape,TopAbs_FACE); exp.More(); exp.Next() )
995   {
996     const TopoDS_Shape& aShapeFace = exp.Current();
997     if ( aProxyMesh )
998       aSubMeshDSFace = aProxyMesh->GetSubMesh( aShapeFace );
999     else
1000       aSubMeshDSFace = meshDS->MeshElements( aShapeFace );
1001     if ( !aSubMeshDSFace )
1002       continue;
1003
1004     vector<const SMDS_MeshElement*> trias, quads;
1005     bool hasNewTrias = false;
1006
1007     const bool toCheckFaceInSolid =
1008       aProxyMesh ? aProxyMesh->HasPrismsOnTwoSides( meshDS->MeshElements( aShapeFace )) : false;
1009     if ( toCheckFaceInSolid && !dynamic_cast< const SMESH_ProxyMesh::SubMesh* >( aSubMeshDSFace ))
1010       continue; // no room for pyramids as prisms are on both sides 
1011
1012     {
1013       bool isRevGlob = false;
1014       SMESH_MesherHelper helper2( aMesh );
1015       PShapeIteratorPtr sIt = helper2.GetAncestors( aShapeFace, aMesh, aShape.ShapeType() );
1016       while ( const TopoDS_Shape * solid = sIt->next() )
1017         if ( !solid->IsSame( aShape ))
1018         {
1019           isRevGlob = helper2.IsReversedSubMesh( TopoDS::Face( aShapeFace ));
1020           if ( toCheckFaceInSolid )
1021             helper2.IsQuadraticSubMesh( *solid );
1022           break;
1023         }
1024
1025       SMDS_ElemIteratorPtr iteratorElem = aSubMeshDSFace->GetElements();
1026       while ( iteratorElem->more() ) // loop on elements on a geometrical face
1027       {
1028         const SMDS_MeshElement* face = iteratorElem->next();
1029
1030         // preparation step to get face info
1031         int stat = Preparation( face, PN, VN, FNodes, PC, VNorm );
1032
1033         bool isRev = isRevGlob;
1034         SMESH_MesherHelper* helper = &helper1;
1035         if ( toCheckFaceInSolid && !isInSolid( FNodes, face->NbCornerNodes(), solidID ))
1036         {
1037           isRev  = !isRevGlob;
1038           helper = &helper2;
1039         }
1040
1041         switch ( stat )
1042         {
1043         case NOT_QUAD:
1044
1045           trias.push_back( face );
1046           break;
1047
1048         case DEGEN_QUAD:
1049           {
1050             // degenerate face
1051             // add triangles to result map
1052             SMDS_MeshFace* NewFace;
1053             helper->SetElementsOnShape( false );
1054             if(!isRev)
1055               NewFace = helper->AddFace( FNodes[0], FNodes[1], FNodes[2] );
1056             else
1057               NewFace = helper->AddFace( FNodes[0], FNodes[2], FNodes[1] );
1058             storeTmpElement( NewFace );
1059             trias.push_back ( NewFace );
1060             quads.push_back( face );
1061             hasNewTrias = true;
1062             break;
1063           }
1064
1065         case QUAD:
1066           {
1067             if(!isRev) VNorm.Reverse();
1068             //double xc = 0., yc = 0., zc = 0.;
1069             double h, hMin = Precision::Infinite();
1070             gp_Pnt PCbest = PC;
1071             int i = 1;
1072             for(; i<=4; i++) {
1073               gp_Pnt Pbest;
1074               if(!isRev)
1075                 Pbest = FindBestPoint(PN(i), PN(i+1), PC, VN(i).Reversed(), h);
1076               else
1077                 Pbest = FindBestPoint(PN(i), PN(i+1), PC, VN(i), h);
1078               if ( 0 < h && h < hMin )
1079               {
1080                 PCbest = Pbest;
1081                 hMin = h;
1082               }
1083             }
1084             //gp_Pnt PCbest(xc/4., yc/4., zc/4.);
1085
1086             // check PCbest
1087             double height = PCbest.Distance(PC);
1088             if ( height < 1.e-6 ) {
1089               // create new PCbest using a bit shift along VNorm
1090               PCbest = PC.XYZ() + VNorm.XYZ() * 0.001;
1091             }
1092             else {
1093               // check possible intersection with other faces
1094               if ( !LimitHeight( PCbest, PC, PN, FNodes, aMesh, face, /*UseApexRay=*/true, aShape ))
1095                 return false;
1096             }
1097             // create node at PCbest
1098             helper->SetElementsOnShape( true );
1099             SMDS_MeshNode* NewNode = helper->AddNode( PCbest.X(), PCbest.Y(), PCbest.Z() );
1100
1101             // create a pyramid
1102             SMDS_MeshVolume* aPyram;
1103             if ( isRev )
1104               aPyram = helper->AddVolume( FNodes[0], FNodes[3], FNodes[2], FNodes[1], NewNode );
1105             else
1106               aPyram = helper->AddVolume( FNodes[0], FNodes[1], FNodes[2], FNodes[3], NewNode );
1107             myPyramids.push_back(aPyram);
1108             //cout << "F" << face->GetID() << " - V" << aPyram->GetID() << endl;
1109
1110             myPyramHeight2.Bind( aPyram, PCbest.SquareDistance( PC ));
1111
1112             // add triangles to result map
1113             helper->SetElementsOnShape( false );
1114             for ( i = 0; i < 4; i++ )
1115             {
1116               trias.push_back ( helper->AddFace( NewNode, FNodes[i], FNodes[i+1] ));
1117               storeTmpElement( trias.back() );
1118             }
1119
1120             quads.push_back( face );
1121             hasNewTrias = true;
1122             break;
1123
1124           } // case QUAD:
1125
1126         } // switch ( stat )
1127       } // end loop on elements on a face submesh
1128
1129       bool sourceSubMeshIsProxy = false;
1130       if ( aProxyMesh )
1131       {
1132         // move proxy sub-mesh from other proxy mesh to this
1133         sourceSubMeshIsProxy = takeProxySubMesh( aShapeFace, aProxyMesh );
1134         // move also tmp elements added in mesh
1135         takeTmpElemsInMesh( aProxyMesh );
1136       }
1137       if ( hasNewTrias )
1138       {
1139         SMESH_ProxyMesh::SubMesh* prxSubMesh = getProxySubMesh( aShapeFace );
1140         prxSubMesh->ChangeElements( trias.begin(), trias.end() );
1141
1142         // delete tmp quadrangles removed from aProxyMesh
1143         if ( sourceSubMeshIsProxy )
1144         {
1145           for ( unsigned i = 0; i < quads.size(); ++i )
1146             removeTmpElement( quads[i] );
1147
1148           delete myElemSearcher;
1149           myElemSearcher =
1150             SMESH_MeshAlgos::GetElementSearcher( *meshDS, aProxyMesh->GetFaces(aShape));
1151         }
1152       }
1153     }
1154   } // end for(TopExp_Explorer exp(aShape,TopAbs_FACE);exp.More();exp.Next()) {
1155
1156   return Compute2ndPart(aMesh, myPyramids);
1157 }
1158
1159 //================================================================================
1160 /*!
1161  * \brief Computes pyramids in mesh with no shape
1162  */
1163 //================================================================================
1164
1165 bool StdMeshers_QuadToTriaAdaptor::Compute(SMESH_Mesh& aMesh)
1166 {
1167   SMESH_ProxyMesh::setMesh( aMesh );
1168   SMESH_ProxyMesh::_allowedTypes.push_back( SMDSEntity_Triangle );
1169   SMESH_ProxyMesh::_allowedTypes.push_back( SMDSEntity_Quad_Triangle );
1170   if ( aMesh.NbQuadrangles() < 1 )
1171     return false;
1172
1173   // find if there is a group of faces identified as skin faces, with normal going outside the volume
1174   std::string groupName = "skinFaces";
1175   SMESHDS_GroupBase* groupDS = 0;
1176   SMESH_Mesh::GroupIteratorPtr groupIt = aMesh.GetGroups();
1177   while ( groupIt->more() )
1178   {
1179     groupDS = 0;
1180     SMESH_Group * group = groupIt->next();
1181     if ( !group ) continue;
1182     groupDS = group->GetGroupDS();
1183     if ( !groupDS || groupDS->IsEmpty() )
1184     {
1185       groupDS = 0;
1186       continue;
1187     }
1188     if (groupDS->GetType() != SMDSAbs_Face)
1189     {
1190       groupDS = 0;
1191       continue;
1192     }
1193     std::string grpName = group->GetName();
1194     if (grpName == groupName)
1195     {
1196       break;
1197     }
1198     else
1199       groupDS = 0;
1200   }
1201
1202   const bool toFindVolumes = aMesh.NbVolumes() > 0;
1203
1204   vector<const SMDS_MeshElement*> myPyramids;
1205   SMESH_MesherHelper helper(aMesh);
1206   helper.IsQuadraticSubMesh(aMesh.GetShapeToMesh());
1207
1208   SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
1209   SMESH_ProxyMesh::SubMesh* prxSubMesh = getProxySubMesh();
1210
1211   if ( !myElemSearcher )
1212     myElemSearcher = SMESH_MeshAlgos::GetElementSearcher( *meshDS );
1213   SMESH_ElementSearcher* searcher = const_cast<SMESH_ElementSearcher*>( myElemSearcher );
1214   SMESHUtils::Deleter<SMESH_ElementSearcher>
1215     volSearcher( SMESH_MeshAlgos::GetElementSearcher( *meshDS ));
1216   vector< const SMDS_MeshElement* > suspectFaces, foundVolumes;
1217
1218   TColgp_Array1OfPnt PN(1,5);
1219   TColgp_Array1OfVec VN(1,4);
1220   vector<const SMDS_MeshNode*> FNodes(5);
1221   TColgp_SequenceOfPnt aContour;
1222
1223   SMDS_FaceIteratorPtr fIt = meshDS->facesIterator();
1224   while( fIt->more())
1225   {
1226     const SMDS_MeshElement* face = fIt->next();
1227     if ( !face ) continue;
1228     // retrieve needed information about a face
1229     gp_Pnt PC;
1230     gp_Vec VNorm;
1231     const SMDS_MeshElement* volumes[2];
1232     int what = Preparation( face, PN, VN, FNodes, PC, VNorm, volumes );
1233     if ( what == NOT_QUAD )
1234       continue;
1235     if ( volumes[0] && volumes[1] )
1236       continue; // face is shared by two volumes - no room for a pyramid
1237
1238     if ( what == DEGEN_QUAD )
1239     {
1240       // degenerate face
1241       // add a triangle to the proxy mesh
1242       SMDS_MeshFace* NewFace;
1243
1244       // check orientation
1245       double tmp = PN(1).Distance(PN(2)) + PN(2).Distance(PN(3));
1246       // far points in VNorm direction
1247       gp_Pnt Ptmp1 = PC.XYZ() + VNorm.XYZ() * tmp * 1.e6;
1248       gp_Pnt Ptmp2 = PC.XYZ() - VNorm.XYZ() * tmp * 1.e6;
1249       // check intersection for Ptmp1 and Ptmp2
1250       bool IsRev = false;
1251       bool IsOK1 = false;
1252       bool IsOK2 = false;
1253       double dist1 = RealLast();
1254       double dist2 = RealLast();
1255       gp_Pnt Pres1,Pres2;
1256
1257       gp_Ax1 line( PC, VNorm );
1258       vector< const SMDS_MeshElement* > suspectFaces;
1259       searcher->GetElementsNearLine( line, SMDSAbs_Face, suspectFaces);
1260
1261       for ( size_t iF = 0; iF < suspectFaces.size(); ++iF ) {
1262         const SMDS_MeshElement* F = suspectFaces[iF];
1263         if ( F == face ) continue;
1264         aContour.Clear();
1265         for ( int i = 0; i < 4; ++i )
1266           aContour.Append( SMESH_TNodeXYZ( F->GetNode(i) ));
1267         gp_Pnt PPP;
1268         if ( !volumes[0] && HasIntersection( Ptmp1, PC, PPP, aContour )) {
1269           IsOK1 = true;
1270           double tmp = PC.Distance(PPP);
1271           if ( tmp < dist1 ) {
1272             Pres1 = PPP;
1273             dist1 = tmp;
1274           }
1275         }
1276         if ( !volumes[1] && HasIntersection( Ptmp2, PC, PPP, aContour )) {
1277           IsOK2 = true;
1278           double tmp = PC.Distance(PPP);
1279           if ( tmp < dist2 ) {
1280             Pres2 = PPP;
1281             dist2 = tmp;
1282           }
1283         }
1284       }
1285
1286       if( IsOK1 && !IsOK2 ) {
1287         // using existed direction
1288       }
1289       else if( !IsOK1 && IsOK2 ) {
1290         // using opposite direction
1291         IsRev = true;
1292       }
1293       else { // IsOK1 && IsOK2
1294         double tmp1 = PC.Distance(Pres1);
1295         double tmp2 = PC.Distance(Pres2);
1296         if(tmp1<tmp2) {
1297           // using existed direction
1298         }
1299         else {
1300           // using opposite direction
1301           IsRev = true;
1302         }
1303       }
1304       helper.SetElementsOnShape( false );
1305       if(!IsRev)
1306         NewFace = helper.AddFace( FNodes[0], FNodes[1], FNodes[2] );
1307       else
1308         NewFace = helper.AddFace( FNodes[0], FNodes[2], FNodes[1] );
1309       storeTmpElement( NewFace );
1310       prxSubMesh->AddElement( NewFace );
1311       continue;
1312     }
1313
1314     // -----------------------------------
1315     // Case of non-degenerated quadrangle
1316     // -----------------------------------
1317
1318     // Find pyramid peak
1319
1320     gp_XYZ PCbest = PC.XYZ();//(0., 0., 0.); // pyramid peak
1321     double h, hMin = Precision::Infinite();
1322     int i = 1;
1323     for ( ; i <= 4; i++ ) {
1324       gp_Pnt Pbest = FindBestPoint(PN(i), PN(i+1), PC, VN(i), h);
1325       if ( 0 < h && h < hMin )
1326       {
1327         PCbest = Pbest.XYZ();
1328         h = hMin;
1329       }
1330       //PCbest += Pbest.XYZ();
1331     }
1332     //PCbest /= 4;
1333
1334     double height = PC.Distance(PCbest); // pyramid height to precise
1335     if ( height < 1.e-6 ) {
1336       // create new PCbest using a bit shift along VNorm
1337       PCbest = PC.XYZ() + VNorm.XYZ() * 0.001;
1338       height = PC.Distance(PCbest);
1339       if ( height < std::numeric_limits<double>::min() )
1340         return false; // batterfly element
1341     }
1342
1343     // Restrict pyramid height by intersection with other faces
1344
1345     gp_Vec tmpDir(PC,PCbest); tmpDir.Normalize();
1346     double tmp = PN(1).Distance(PN(3)) + PN(2).Distance(PN(4));
1347     // far points: in (PC, PCbest) direction and vice-versa
1348     gp_Pnt farPnt[2] = { PC.XYZ() + tmpDir.XYZ() * tmp * 1.e6,
1349                          PC.XYZ() - tmpDir.XYZ() * tmp * 1.e6 };
1350     // check intersection for farPnt1 and farPnt2
1351     bool   intersected[2] = { false, false };
1352     double dist2int   [2] = { RealLast(), RealLast() };
1353     gp_Pnt intPnt     [2];
1354     int    intFaceInd [2] = { 0, 0 };
1355
1356     if ( toFindVolumes && 0 ) // non-conformal mesh is not suitable for any mesher so far
1357     {
1358       // there are volumes in the mesh, in a non-conformal mesh a neighbor
1359       // volume can be not found yet
1360       for ( int isRev = 0; isRev < 2; ++isRev )
1361       {
1362         if ( volumes[isRev] ) continue;
1363         gp_Pnt testPnt = PC.XYZ() + tmpDir.XYZ() * height * ( isRev ? -0.1: 0.1 );
1364         foundVolumes.clear();
1365         if ( volSearcher->FindElementsByPoint( testPnt, SMDSAbs_Volume, foundVolumes ))
1366           volumes[isRev] = foundVolumes[0];
1367       }
1368       if ( volumes[0] && volumes[1] )
1369         continue; // no room for a pyramid
1370     }
1371
1372     gp_Ax1 line( PC, tmpDir );
1373     suspectFaces.clear();
1374     searcher->GetElementsNearLine( line, SMDSAbs_Face, suspectFaces);
1375
1376     for ( size_t iF = 0; iF < suspectFaces.size(); ++iF )
1377     {
1378       const SMDS_MeshElement* F = suspectFaces[iF];
1379       if ( F == face ) continue;
1380       aContour.Clear();
1381       int nbN = F->NbCornerNodes();
1382       for ( i = 0; i < nbN; ++i )
1383         aContour.Append( SMESH_TNodeXYZ( F->GetNode(i) ));
1384       gp_Pnt intP;
1385       for ( int isRev = 0; isRev < 2; ++isRev )
1386       {
1387         if( !volumes[isRev] && HasIntersection( farPnt[isRev], PC, intP, aContour ))
1388         {
1389           double d = PC.Distance( intP );
1390           if ( d < dist2int[isRev] )
1391           {
1392             intersected[isRev] = true;
1393             intPnt     [isRev] = intP;
1394             dist2int   [isRev] = d;
1395             intFaceInd [isRev] = iF;
1396           }
1397         }
1398       }
1399     }
1400
1401     // if the face belong to the group of skinFaces, do not build a pyramid outside
1402     if ( groupDS && groupDS->Contains(face) )
1403     {
1404       intersected[0] = false;
1405     }
1406     else if ( intersected[0] && intersected[1] ) // check if one of pyramids is in a hole
1407     {
1408       gp_Pnt P ( PC.XYZ() + tmpDir.XYZ() * 0.5 * dist2int[0] );
1409       if ( searcher->GetPointState( P ) == TopAbs_OUT )
1410         intersected[0] = false;
1411       else
1412       {
1413         P = ( PC.XYZ() - tmpDir.XYZ() * 0.5 * dist2int[1] );
1414         if ( searcher->GetPointState( P ) == TopAbs_OUT )
1415           intersected[1] = false;
1416       }
1417     }
1418
1419     // Create one or two pyramids
1420
1421     for ( int isRev = 0; isRev < 2; ++isRev )
1422     {
1423       if ( !intersected[isRev] ) continue;
1424       double pyramidH = Min( height, dist2int[isRev]/3. );
1425       gp_Pnt    Papex = PC.XYZ() + tmpDir.XYZ() * (isRev ? -pyramidH : pyramidH);
1426       if ( pyramidH < 1e-2 * height )
1427         return overlapError( aMesh, face, suspectFaces[ intFaceInd[isRev] ] );
1428
1429       if ( !LimitHeight( Papex, PC, PN, FNodes, aMesh, face, /*UseApexRay=*/false ))
1430         return false;
1431
1432       // create node for Papex
1433       helper.SetElementsOnShape( true );
1434       SMDS_MeshNode* NewNode = helper.AddNode( Papex.X(), Papex.Y(), Papex.Z() );
1435
1436       // create a pyramid
1437       SMDS_MeshVolume* aPyram;
1438       if(isRev)
1439         aPyram = helper.AddVolume( FNodes[0], FNodes[1], FNodes[2], FNodes[3], NewNode );
1440       else
1441         aPyram = helper.AddVolume( FNodes[0], FNodes[3], FNodes[2], FNodes[1], NewNode );
1442       myPyramids.push_back(aPyram);
1443
1444       //myPyramHeight2.Bind( aPyram, Papex.SquareDistance( PC ));
1445
1446       // add triangles to result map
1447       helper.SetElementsOnShape( false );
1448       for ( i = 0; i < 4; i++) {
1449         SMDS_MeshFace* NewFace;
1450         if(isRev)
1451           NewFace = helper.AddFace( NewNode, FNodes[i], FNodes[i+1] );
1452         else
1453           NewFace = helper.AddFace( NewNode, FNodes[i+1], FNodes[i] );
1454         storeTmpElement( NewFace );
1455         prxSubMesh->AddElement( NewFace );
1456       }
1457     }
1458   } // end loop on all faces
1459
1460   return Compute2ndPart(aMesh, myPyramids);
1461 }
1462
1463 //================================================================================
1464 /*!
1465  * \brief Update created pyramids and faces to avoid their intersection
1466  */
1467 //================================================================================
1468
1469 bool StdMeshers_QuadToTriaAdaptor::Compute2ndPart(SMESH_Mesh&                            aMesh,
1470                                                   const vector<const SMDS_MeshElement*>& myPyramids)
1471 {
1472   if ( myPyramids.empty() )
1473     return true;
1474
1475   SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
1476   size_t i, j, k;
1477   //int myShapeID = myPyramids[0]->GetNode(4)->getshapeId();
1478   {
1479     SMDS_ElemIteratorPtr
1480       pyramIt( new SMDS_ElementVectorIterator( myPyramids.begin(), myPyramids.end() ));
1481     if ( myElemSearcher ) delete myElemSearcher;
1482     myElemSearcher = SMESH_MeshAlgos::GetElementSearcher( *meshDS, pyramIt );
1483   }
1484   SMESH_ElementSearcher* searcher = const_cast<SMESH_ElementSearcher*>( myElemSearcher );
1485
1486   set<const SMDS_MeshNode*> nodesToMove;
1487
1488   // check adjacent pyramids
1489
1490   // for ( i = 0; i <  myPyramids.size(); ++i )
1491   // {
1492   //   const SMDS_MeshElement* PrmI = myPyramids[i];
1493   //   MergeAdjacent( PrmI, nodesToMove );
1494   // }
1495
1496   // Fix adjacent pyramids whose heights differ too much
1497
1498   myNodes.resize(2);
1499   bool modifHeight = true;
1500   typedef NCollection_DataMap< const SMDS_MeshElement*, double >::Iterator TPyramToH2Iter;
1501   while ( modifHeight )
1502   {
1503     modifHeight = false;
1504     for ( TPyramToH2Iter pyramToH2( myPyramHeight2 ); pyramToH2.More(); pyramToH2.Next() )
1505       modifHeight |= DecreaseHeightDifference( pyramToH2.Key(), pyramToH2.Value() );
1506   }
1507   for ( TPyramToH2Iter pyramToH2( myPyramHeight2 ); pyramToH2.More(); pyramToH2.Next() )
1508   {
1509     if ( pyramToH2.Value() > 0 )
1510       continue; // not changed
1511     const double                h = Sqrt( - pyramToH2.Value() );
1512     const SMDS_MeshElement* pyram = pyramToH2.Key();
1513     SMESH_NodeXYZ Papex = pyram->GetNode( PYRAM_APEX );
1514     gp_XYZ PC( 0,0,0 );
1515     for ( int i = 0; i < 4; ++i )
1516       PC += SMESH_NodeXYZ( pyram->GetNode( i ));
1517     PC /= 4;
1518     gp_Vec V( PC, Papex );
1519     gp_Pnt newApex = gp_Pnt( PC ).Translated( h * V.Normalized() );
1520     meshDS->MoveNode( Papex.Node(), newApex.X(), newApex.Y(), newApex.Z() );
1521   }
1522
1523   // iterate on all new pyramids
1524   vector< const SMDS_MeshElement* > suspectPyrams;
1525   for ( i = 0; i <  myPyramids.size(); ++i )
1526   {
1527     const SMDS_MeshElement*  PrmI = myPyramids[i];
1528     const SMDS_MeshNode*    apexI = PrmI->GetNode( PYRAM_APEX );
1529
1530     // compare PrmI with all the rest pyramids
1531
1532     // collect adjacent pyramids and nodes coordinates of PrmI
1533     set<const SMDS_MeshElement*> checkedPyrams;
1534     gp_Pnt PsI[5];
1535     for ( k = 0; k < 5; k++ )
1536     {
1537       const SMDS_MeshNode* n = PrmI->GetNode(k);
1538       PsI[k] = SMESH_TNodeXYZ( n );
1539       SMDS_ElemIteratorPtr vIt = n->GetInverseElementIterator( SMDSAbs_Volume );
1540       while ( vIt->more() )
1541       {
1542         const SMDS_MeshElement* PrmJ = vIt->next();
1543         if ( SMESH_MeshAlgos::NbCommonNodes( PrmI, PrmJ ) > 1 )
1544           checkedPyrams.insert( PrmJ );
1545       }
1546     }
1547
1548     // get pyramids to check
1549     gp_XYZ       PC = ( PsI[0].XYZ() + PsI[1].XYZ() + PsI[2].XYZ() + PsI[3].XYZ() ) / 4.;
1550     gp_XYZ      ray = PsI[4].XYZ() - PC;
1551     gp_XYZ   center = PC + 0.5 * ray;
1552     double diameter = Max( PsI[0].Distance(PsI[2]), PsI[1].Distance(PsI[3]));
1553     suspectPyrams.clear();
1554     searcher->GetElementsInSphere( center, diameter * 0.6, SMDSAbs_Volume, suspectPyrams);
1555
1556     // check intersection with distant pyramids
1557     for ( j = 0; j < suspectPyrams.size(); ++j )
1558     {
1559       const SMDS_MeshElement* PrmJ = suspectPyrams[j];
1560       if ( PrmJ == PrmI )
1561         continue;
1562       if ( apexI == PrmJ->GetNode( PYRAM_APEX ))
1563         continue; // pyramids PrmI and PrmJ already merged
1564       if ( !checkedPyrams.insert( PrmJ ).second )
1565         continue; // already checked
1566
1567       gp_Pnt PsJ[5];
1568       for ( k = 0; k < 5; k++ )
1569         PsJ[k] = SMESH_TNodeXYZ( PrmJ->GetNode(k) );
1570
1571       if ( ray * ( PsJ[4].XYZ() - PC ) < 0. )
1572         continue; // PrmJ is below PrmI
1573
1574       for ( k = 0; k < 4; k++ ) // loop on 4 base nodes of PrmI
1575       {
1576         gp_Pnt Pint;
1577         bool hasInt=false;
1578         for ( k = 0; k < 4  &&  !hasInt; k++ )
1579         {
1580           gp_Vec Vtmp( PsI[k], PsI[ PYRAM_APEX ]);
1581           gp_Pnt Pshift = PsI[k].XYZ() + Vtmp.XYZ() * 0.01; // base node moved a bit to apex
1582           hasInt =
1583             ( HasIntersection3( Pshift, PsI[4], Pint, PsJ[0], PsJ[1], PsJ[PYRAM_APEX]) ||
1584               HasIntersection3( Pshift, PsI[4], Pint, PsJ[1], PsJ[2], PsJ[PYRAM_APEX]) ||
1585               HasIntersection3( Pshift, PsI[4], Pint, PsJ[2], PsJ[3], PsJ[PYRAM_APEX]) ||
1586               HasIntersection3( Pshift, PsI[4], Pint, PsJ[3], PsJ[0], PsJ[PYRAM_APEX]) );
1587         }
1588         for ( k = 0; k < 4  &&  !hasInt; k++ )
1589         {
1590           gp_Vec Vtmp( PsJ[k], PsJ[ PYRAM_APEX ]);
1591           gp_Pnt Pshift = PsJ[k].XYZ() + Vtmp.XYZ() * 0.01;
1592           hasInt =
1593             ( HasIntersection3( Pshift, PsJ[4], Pint, PsI[0], PsI[1], PsI[PYRAM_APEX]) ||
1594               HasIntersection3( Pshift, PsJ[4], Pint, PsI[1], PsI[2], PsI[PYRAM_APEX]) ||
1595               HasIntersection3( Pshift, PsJ[4], Pint, PsI[2], PsI[3], PsI[PYRAM_APEX]) ||
1596               HasIntersection3( Pshift, PsJ[4], Pint, PsI[3], PsI[0], PsI[PYRAM_APEX]) );
1597         }
1598
1599         if ( hasInt )
1600         {
1601           // count common nodes of base faces of two pyramids
1602           int nbc = 0;
1603           for ( k = 0; k < 4; k++ )
1604             nbc += int ( PrmI->GetNodeIndex( PrmJ->GetNode(k) ) >= 0 );
1605
1606           if ( nbc == 4 )
1607             continue; // pyrams have a common base face
1608
1609           if ( nbc > 1 )
1610           {
1611             // Merge the two pyramids and others already merged with them
1612             MergePiramids( PrmI, PrmJ, nodesToMove );
1613           }
1614           else  // nbc==0
1615           {
1616             // decrease height of pyramids
1617             gp_XYZ PCi(0,0,0), PCj(0,0,0);
1618             for ( k = 0; k < 4; k++ ) {
1619               PCi += PsI[k].XYZ();
1620               PCj += PsJ[k].XYZ();
1621             }
1622             PCi /= 4; PCj /= 4;
1623             gp_Vec VN1(PCi,PsI[4]);
1624             gp_Vec VN2(PCj,PsJ[4]);
1625             gp_Vec VI1(PCi,Pint);
1626             gp_Vec VI2(PCj,Pint);
1627             double ang1 = fabs(VN1.Angle(VI1));
1628             double ang2 = fabs(VN2.Angle(VI2));
1629             double coef1 = 0.5 - (( ang1 < M_PI/3. ) ? cos(ang1)*0.25 : 0 );
1630             double coef2 = 0.5 - (( ang2 < M_PI/3. ) ? cos(ang2)*0.25 : 0 ); // cos(ang2) ?
1631             //             double coef2 = 0.5;
1632             //             if(ang2<PI/3)
1633             //               coef2 -= cos(ang1)*0.25;
1634
1635             VN1.Scale(coef1);
1636             VN2.Scale(coef2);
1637             SMDS_MeshNode* aNode1 = const_cast<SMDS_MeshNode*>( apexI );
1638             aNode1->setXYZ( PCi.X()+VN1.X(), PCi.Y()+VN1.Y(), PCi.Z()+VN1.Z() );
1639             SMDS_MeshNode* aNode2 = const_cast<SMDS_MeshNode*>(PrmJ->GetNode( PYRAM_APEX ));
1640             aNode2->setXYZ( PCj.X()+VN2.X(), PCj.Y()+VN2.Y(), PCj.Z()+VN2.Z() );
1641             nodesToMove.insert( aNode1 );
1642             nodesToMove.insert( aNode2 );
1643             //cout << "Limit H  F" << getFaceID( PrmI ) << " - F" << getFaceID( PrmJ ) << endl;
1644           }
1645           // fix intersections that can appear after apex movement
1646           //MergeAdjacent( PrmI, nodesToMove );
1647           //MergeAdjacent( PrmJ, nodesToMove );
1648
1649           apexI = PrmI->GetNode( PYRAM_APEX ); // apexI can be removed by merge
1650
1651         } // end if(hasInt)
1652       } // loop on 4 base nodes of PrmI
1653     }  // loop on suspectPyrams
1654
1655   } // loop on all pyramids
1656
1657   //smIdType nbNodes = aMesh.NbNodes();
1658   for ( i = 0; i <  myPyramids.size(); ++i )
1659   {
1660     const SMDS_MeshElement* PrmI = myPyramids[i];
1661     MergeAdjacent( PrmI, nodesToMove );
1662   }
1663
1664   if( !nodesToMove.empty() && !meshDS->IsEmbeddedMode() )
1665   {
1666     set<const SMDS_MeshNode*>::iterator n = nodesToMove.begin();
1667     for ( ; n != nodesToMove.end(); ++n )
1668       meshDS->MoveNode( *n, (*n)->X(), (*n)->Y(), (*n)->Z() );
1669   }
1670
1671   // move medium nodes of merged quadratic pyramids
1672   if ( myPyramids[0]->IsQuadratic() )
1673     UpdateQuadraticPyramids( nodesToMove, GetMeshDS() );
1674
1675   // erase removed triangles from the proxy mesh
1676   if ( !myRemovedTrias.empty() )
1677   {
1678     for ( int i = 0; i <= meshDS->MaxShapeIndex(); ++i )
1679       if ( SMESH_ProxyMesh::SubMesh* sm = findProxySubMesh(i))
1680       {
1681         vector<const SMDS_MeshElement *> faces;
1682         faces.reserve( sm->NbElements() );
1683         SMDS_ElemIteratorPtr fIt = sm->GetElements();
1684         while ( fIt->more() )
1685         {
1686           const SMDS_MeshElement* tria = fIt->next();
1687           set<const SMDS_MeshElement*>::iterator rmTria = myRemovedTrias.find( tria );
1688           if ( rmTria != myRemovedTrias.end() )
1689             myRemovedTrias.erase( rmTria );
1690           else
1691             faces.push_back( tria );
1692         }
1693         sm->ChangeElements( faces.begin(), faces.end() );
1694       }
1695   }
1696
1697   myDegNodes.clear();
1698
1699   delete myElemSearcher;
1700   myElemSearcher=0;
1701
1702   return true;
1703 }