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