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