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