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