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