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