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