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