Version: 6.3.1

src/StdMeshers/StdMeshers_QuadToTriaAdaptor.cxx

Go to the documentation of this file.
00001 // Copyright (C) 2007-2011  CEA/DEN, EDF R&D, OPEN CASCADE
00002 //
00003 // This library is free software; you can redistribute it and/or
00004 // modify it under the terms of the GNU Lesser General Public
00005 // License as published by the Free Software Foundation; either
00006 // version 2.1 of the License.
00007 //
00008 // This library is distributed in the hope that it will be useful,
00009 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00010 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00011 // Lesser General Public License for more details.
00012 //
00013 // You should have received a copy of the GNU Lesser General Public
00014 // License along with this library; if not, write to the Free Software
00015 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
00016 //
00017 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
00018 //
00019 
00020 // File      : StdMeshers_QuadToTriaAdaptor.cxx
00021 // Module    : SMESH
00022 // Created   : Wen May 07 16:37:07 2008
00023 // Author    : Sergey KUUL (skl)
00024 //
00025 #include "StdMeshers_QuadToTriaAdaptor.hxx"
00026 
00027 #include "SMDS_SetIterator.hxx"
00028 
00029 #include "SMESH_Algo.hxx"
00030 #include "SMESH_MesherHelper.hxx"
00031 
00032 #include <IntAna_IntConicQuad.hxx>
00033 #include <IntAna_Quadric.hxx>
00034 #include <TColgp_HArray1OfPnt.hxx>
00035 #include <TColgp_HArray1OfVec.hxx>
00036 #include <TColgp_HSequenceOfPnt.hxx>
00037 #include <TopExp_Explorer.hxx>
00038 #include <TopoDS.hxx>
00039 #include <gp_Lin.hxx>
00040 #include <gp_Pln.hxx>
00041 
00042 #include <numeric>
00043 
00044 using namespace std;
00045 
00046 enum EQuadNature { NOT_QUAD, QUAD, DEGEN_QUAD, PYRAM_APEX = 4, TRIA_APEX = 0 };
00047 
00048 // std-like iterator used to get coordinates of nodes of mesh element
00049 typedef SMDS_StdIterator< SMESH_TNodeXYZ, SMDS_ElemIteratorPtr > TXyzIterator;
00050 
00051 namespace
00052 {
00053   //================================================================================
00057   //================================================================================
00058 
00059   bool EqualTriangles(const SMDS_MeshElement* F1,const SMDS_MeshElement* F2)
00060   {
00061     return
00062       ( F1->GetNode(1)==F2->GetNode(2) && F1->GetNode(2)==F2->GetNode(1) ) ||
00063       ( F1->GetNode(1)==F2->GetNode(1) && F1->GetNode(2)==F2->GetNode(2) );
00064   }
00065   //================================================================================
00070   //================================================================================
00071 
00072   bool TooCloseAdjacent( const SMDS_MeshElement* PrmI,
00073                          const SMDS_MeshElement* PrmJ,
00074                          const bool              hasShape)
00075   {
00076     const SMDS_MeshNode* nApexI = PrmI->GetNode(4);
00077     const SMDS_MeshNode* nApexJ = PrmJ->GetNode(4);
00078     if ( nApexI == nApexJ ||
00079          nApexI->getshapeId() != nApexJ->getshapeId() )
00080       return false;
00081 
00082     // Find two common base nodes and their indices within PrmI and PrmJ
00083     const SMDS_MeshNode* baseNodes[2] = { 0,0 };
00084     int baseNodesIndI[2], baseNodesIndJ[2];
00085     for ( int i = 0; i < 4 ; ++i )
00086     {
00087       int j = PrmJ->GetNodeIndex( PrmI->GetNode(i));
00088       if ( j >= 0 )
00089       {
00090         int ind = baseNodes[0] ? 1:0;
00091         if ( baseNodes[ ind ])
00092           return false; // pyramids with a common base face
00093         baseNodes    [ ind ] = PrmI->GetNode(i);
00094         baseNodesIndI[ ind ] = i;
00095         baseNodesIndJ[ ind ] = j;
00096       }
00097     }
00098     if ( !baseNodes[1] ) return false; // not adjacent
00099 
00100     // Get normals of triangles sharing baseNodes
00101     gp_XYZ apexI = SMESH_TNodeXYZ( nApexI );
00102     gp_XYZ apexJ = SMESH_TNodeXYZ( nApexJ );
00103     gp_XYZ base1 = SMESH_TNodeXYZ( baseNodes[0]);
00104     gp_XYZ base2 = SMESH_TNodeXYZ( baseNodes[1]);
00105     gp_Vec baseVec( base1, base2 );
00106     gp_Vec baI( base1, apexI );
00107     gp_Vec baJ( base1, apexJ );
00108     gp_Vec nI = baseVec.Crossed( baI );
00109     gp_Vec nJ = baseVec.Crossed( baJ );
00110 
00111     // Check angle between normals
00112     double angle = nI.Angle( nJ );
00113     bool tooClose = ( angle < 15 * PI180 );
00114 
00115     // Check if pyramids collide
00116     if ( !tooClose && baI * baJ > 0 )
00117     {
00118       // find out if nI points outside of PrmI or inside
00119       int dInd = baseNodesIndI[1] - baseNodesIndI[0];
00120       bool isOutI = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0;
00121 
00122       // find out sign of projection of nJ to baI
00123       double proj = baI * nJ;
00124 
00125       tooClose = isOutI ? proj > 0 : proj < 0;
00126     }
00127 
00128     // Check if PrmI and PrmJ are in same domain
00129     if ( tooClose && !hasShape )
00130     {
00131       // check order of baseNodes within pyramids, it must be opposite
00132       int dInd;
00133       dInd = baseNodesIndI[1] - baseNodesIndI[0];
00134       bool isOutI = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0;
00135       dInd = baseNodesIndJ[1] - baseNodesIndJ[0];
00136       bool isOutJ = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0;
00137       if ( isOutJ == isOutI )
00138         return false; // other domain
00139 
00140       // direct both normals outside pyramid
00141       ( isOutI ? nJ : nI ).Reverse();
00142 
00143       // check absence of a face separating domains between pyramids
00144       TIDSortedElemSet emptySet, avoidSet;
00145       int i1, i2;
00146       while ( const SMDS_MeshElement* f =
00147               SMESH_MeshEditor::FindFaceInSet( baseNodes[0], baseNodes[1],
00148                                                emptySet, avoidSet, &i1, &i2 ))
00149       {
00150         avoidSet.insert( f );
00151 
00152         // face node other than baseNodes
00153         int otherNodeInd = 0;
00154         while ( otherNodeInd == i1 || otherNodeInd == i2 ) otherNodeInd++;
00155         const SMDS_MeshNode* otherFaceNode = f->GetNode( otherNodeInd );
00156 
00157         if ( otherFaceNode == nApexI || otherFaceNode == nApexJ )
00158           continue; // f is a temporary triangle
00159 
00160         // check if f is a base face of either of pyramids
00161         if ( f->NbCornerNodes() == 4 &&
00162              ( PrmI->GetNodeIndex( otherFaceNode ) >= 0 ||
00163                PrmJ->GetNodeIndex( otherFaceNode ) >= 0 ))
00164           continue; // f is a base quadrangle
00165 
00166         // check projections of face direction (baOFN) to triange normals (nI and nJ)
00167         gp_Vec baOFN( base1, SMESH_TNodeXYZ( otherFaceNode ));
00168         if ( nI * baOFN > 0 && nJ * baOFN > 0 )
00169         {
00170           tooClose = false; // f is between pyramids
00171           break;
00172         }
00173       }
00174     }
00175 
00176     return tooClose;
00177   }
00178 
00179   //================================================================================
00183   //================================================================================
00184 
00185   void UpdateQuadraticPyramids(const set<const SMDS_MeshNode*>& commonApex,
00186                                SMESHDS_Mesh*                    meshDS)
00187   {
00188     typedef SMDS_StdIterator< const SMDS_MeshElement*, SMDS_ElemIteratorPtr > TStdElemIterator;
00189     TStdElemIterator itEnd;
00190 
00191     // shift of node index to get medium nodes between the 4 base nodes and the apex
00192     const int base2MediumShift = 9;
00193 
00194     set<const SMDS_MeshNode*>::const_iterator nIt = commonApex.begin();
00195     for ( ; nIt != commonApex.end(); ++nIt )
00196     {
00197       SMESH_TNodeXYZ apex( *nIt );
00198 
00199       vector< const SMDS_MeshElement* > pyrams // pyramids sharing the apex node
00200         ( TStdElemIterator( apex._node->GetInverseElementIterator( SMDSAbs_Volume )), itEnd );
00201 
00202       // Select medium nodes to keep and medium nodes to remove
00203 
00204       typedef map < const SMDS_MeshNode*, const SMDS_MeshNode*, TIDCompare > TN2NMap;
00205       TN2NMap base2medium; // to keep
00206       vector< const SMDS_MeshNode* > nodesToRemove;
00207 
00208       for ( unsigned i = 0; i < pyrams.size(); ++i )
00209         for ( int baseIndex = 0; baseIndex < PYRAM_APEX; ++baseIndex )
00210         {
00211           SMESH_TNodeXYZ         base = pyrams[i]->GetNode( baseIndex );
00212           const SMDS_MeshNode* medium = pyrams[i]->GetNode( baseIndex + base2MediumShift );
00213           TN2NMap::iterator b2m = base2medium.insert( make_pair( base._node, medium )).first;
00214           if ( b2m->second != medium )
00215           {
00216             nodesToRemove.push_back( medium );
00217           }
00218           else
00219           {
00220             // move the kept medium node
00221             gp_XYZ newXYZ = 0.5 * ( apex + base );
00222             meshDS->MoveNode( medium, newXYZ.X(), newXYZ.Y(), newXYZ.Z() );
00223           }
00224         }
00225 
00226       // Within pyramids, replace nodes to remove by nodes to keep  
00227 
00228       for ( unsigned i = 0; i < pyrams.size(); ++i )
00229       {
00230         vector< const SMDS_MeshNode* > nodes( pyrams[i]->begin_nodes(),
00231                                               pyrams[i]->end_nodes() );
00232         for ( int baseIndex = 0; baseIndex < PYRAM_APEX; ++baseIndex )
00233         {
00234           const SMDS_MeshNode* base = pyrams[i]->GetNode( baseIndex );
00235           nodes[ baseIndex + base2MediumShift ] = base2medium[ base ];
00236         }
00237         meshDS->ChangeElementNodes( pyrams[i], &nodes[0], nodes.size());
00238       }
00239 
00240       // Remove the replaced nodes
00241 
00242       if ( !nodesToRemove.empty() )
00243       {
00244         SMESHDS_SubMesh * sm = meshDS->MeshElements( nodesToRemove[0]->getshapeId() );
00245         for ( unsigned i = 0; i < nodesToRemove.size(); ++i )
00246           meshDS->RemoveFreeNode( nodesToRemove[i], sm, /*fromGroups=*/false);
00247       }
00248     }
00249   }
00250 
00251 }
00252 
00253 //================================================================================
00257 //================================================================================
00258 
00259 void StdMeshers_QuadToTriaAdaptor::MergePiramids( const SMDS_MeshElement*     PrmI,
00260                                                   const SMDS_MeshElement*     PrmJ,
00261                                                   set<const SMDS_MeshNode*> & nodesToMove)
00262 {
00263   const SMDS_MeshNode* Nrem = PrmJ->GetNode(4); // node to remove
00264   //int nbJ = Nrem->NbInverseElements( SMDSAbs_Volume );
00265   SMESH_TNodeXYZ Pj( Nrem );
00266 
00267   // an apex node to make common to all merged pyramids
00268   SMDS_MeshNode* CommonNode = const_cast<SMDS_MeshNode*>(PrmI->GetNode(4));
00269   if ( CommonNode == Nrem ) return; // already merged
00270   //int nbI = CommonNode->NbInverseElements( SMDSAbs_Volume );
00271   SMESH_TNodeXYZ Pi( CommonNode );
00272   gp_XYZ Pnew = /*( nbI*Pi + nbJ*Pj ) / (nbI+nbJ);*/ 0.5 * ( Pi + Pj );
00273   CommonNode->setXYZ( Pnew.X(), Pnew.Y(), Pnew.Z() );
00274 
00275   nodesToMove.insert( CommonNode );
00276   nodesToMove.erase ( Nrem );
00277 
00278   typedef SMDS_StdIterator< const SMDS_MeshElement*, SMDS_ElemIteratorPtr > TStdElemIterator;
00279   TStdElemIterator itEnd;
00280 
00281   // find and remove coincided faces of merged pyramids
00282   vector< const SMDS_MeshElement* > inverseElems
00283     // copy inverse elements to avoid iteration on changing container 
00284     ( TStdElemIterator( CommonNode->GetInverseElementIterator(SMDSAbs_Face)), itEnd);
00285   for ( unsigned i = 0; i < inverseElems.size(); ++i )
00286   {
00287     const SMDS_MeshElement* FI = inverseElems[i];
00288     const SMDS_MeshElement* FJEqual = 0;
00289     SMDS_ElemIteratorPtr triItJ = Nrem->GetInverseElementIterator(SMDSAbs_Face);
00290     while ( !FJEqual && triItJ->more() )
00291     {
00292       const SMDS_MeshElement* FJ = triItJ->next();
00293       if ( EqualTriangles( FJ, FI ))
00294         FJEqual = FJ;
00295     }
00296     if ( FJEqual )
00297     {
00298       removeTmpElement( FI );
00299       removeTmpElement( FJEqual );
00300       myRemovedTrias.insert( FI );
00301       myRemovedTrias.insert( FJEqual );
00302     }
00303   }
00304 
00305   // set the common apex node to pyramids and triangles merged with J
00306   inverseElems.assign( TStdElemIterator( Nrem->GetInverseElementIterator()), itEnd );
00307   for ( unsigned i = 0; i < inverseElems.size(); ++i )
00308   {
00309     const SMDS_MeshElement* elem = inverseElems[i];
00310     vector< const SMDS_MeshNode* > nodes( elem->begin_nodes(), elem->end_nodes() );
00311     nodes[ elem->GetType() == SMDSAbs_Volume ? PYRAM_APEX : TRIA_APEX ] = CommonNode;
00312     GetMeshDS()->ChangeElementNodes( elem, &nodes[0], nodes.size());
00313   }
00314   ASSERT( Nrem->NbInverseElements() == 0 );
00315   GetMeshDS()->RemoveFreeNode( Nrem,
00316                                GetMeshDS()->MeshElements( Nrem->getshapeId()),
00317                                /*fromGroups=*/false);
00318 }
00319 
00320 //================================================================================
00324 //================================================================================
00325 
00326 void StdMeshers_QuadToTriaAdaptor::MergeAdjacent(const SMDS_MeshElement*    PrmI,
00327                                                  set<const SMDS_MeshNode*>& nodesToMove)
00328 {
00329   TIDSortedElemSet adjacentPyrams;
00330   bool mergedPyrams = false;
00331   for(int k=0; k<4; k++) // loop on 4 base nodes of PrmI
00332   {
00333     const SMDS_MeshNode* n = PrmI->GetNode(k);
00334     SMDS_ElemIteratorPtr vIt = n->GetInverseElementIterator( SMDSAbs_Volume );
00335     while ( vIt->more() )
00336     {
00337       const SMDS_MeshElement* PrmJ = vIt->next();
00338       if ( PrmJ->NbCornerNodes() != 5 || !adjacentPyrams.insert( PrmJ ).second  )
00339         continue;
00340       if ( PrmI != PrmJ && TooCloseAdjacent( PrmI, PrmJ, GetMesh()->HasShapeToMesh() ))
00341       {
00342         MergePiramids( PrmI, PrmJ, nodesToMove );
00343         mergedPyrams = true;
00344         // container of inverse elements can change
00345         vIt = n->GetInverseElementIterator( SMDSAbs_Volume );
00346       }
00347     }
00348   }
00349   if ( mergedPyrams )
00350   {
00351     TIDSortedElemSet::iterator prm;
00352     for (prm = adjacentPyrams.begin(); prm != adjacentPyrams.end(); ++prm)
00353       MergeAdjacent( *prm, nodesToMove );
00354   }
00355 }
00356 
00357 //================================================================================
00361 //================================================================================
00362 
00363 StdMeshers_QuadToTriaAdaptor::StdMeshers_QuadToTriaAdaptor():
00364   myElemSearcher(0)
00365 {
00366 }
00367 
00368 //================================================================================
00372 //================================================================================
00373 
00374 StdMeshers_QuadToTriaAdaptor::~StdMeshers_QuadToTriaAdaptor()
00375 {
00376   // temporary faces are deleted by ~SMESH_ProxyMesh()
00377   if ( myElemSearcher ) delete myElemSearcher;
00378   myElemSearcher=0;
00379 }
00380 
00381 
00382 //=======================================================================
00383 //function : FindBestPoint
00384 //purpose  : Return a point P laying on the line (PC,V) so that triangle
00385 //           (P, P1, P2) to be equilateral as much as possible
00386 //           V - normal to (P1,P2,PC)
00387 //=======================================================================
00388 static gp_Pnt FindBestPoint(const gp_Pnt& P1, const gp_Pnt& P2,
00389                             const gp_Pnt& PC, const gp_Vec& V)
00390 {
00391   double a = P1.Distance(P2);
00392   double b = P1.Distance(PC);
00393   double c = P2.Distance(PC);
00394   if( a < (b+c)/2 )
00395     return PC;
00396   else {
00397     // find shift along V in order a to became equal to (b+c)/2
00398     double shift = sqrt( a*a + (b*b-c*c)*(b*b-c*c)/16/a/a - (b*b+c*c)/2 );
00399     gp_Dir aDir(V);
00400     gp_Pnt Pbest = PC.XYZ() + aDir.XYZ() * shift;
00401     return Pbest;
00402   }
00403 }
00404 
00405 
00406 //=======================================================================
00407 //function : HasIntersection3
00408 //purpose  : Auxilare for HasIntersection()
00409 //           find intersection point between triangle (P1,P2,P3)
00410 //           and segment [PC,P]
00411 //=======================================================================
00412 static bool HasIntersection3(const gp_Pnt& P, const gp_Pnt& PC, gp_Pnt& Pint,
00413                              const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3)
00414 {
00415   //cout<<"HasIntersection3"<<endl;
00416   //cout<<"  PC("<<PC.X()<<","<<PC.Y()<<","<<PC.Z()<<")"<<endl;
00417   //cout<<"  P("<<P.X()<<","<<P.Y()<<","<<P.Z()<<")"<<endl;
00418   //cout<<"  P1("<<P1.X()<<","<<P1.Y()<<","<<P1.Z()<<")"<<endl;
00419   //cout<<"  P2("<<P2.X()<<","<<P2.Y()<<","<<P2.Z()<<")"<<endl;
00420   //cout<<"  P3("<<P3.X()<<","<<P3.Y()<<","<<P3.Z()<<")"<<endl;
00421   gp_Vec VP1(P1,P2);
00422   gp_Vec VP2(P1,P3);
00423   IntAna_Quadric IAQ(gp_Pln(P1,VP1.Crossed(VP2)));
00424   IntAna_IntConicQuad IAICQ(gp_Lin(PC,gp_Dir(gp_Vec(PC,P))),IAQ);
00425   if(IAICQ.IsDone()) {
00426     if( IAICQ.IsInQuadric() )
00427       return false;
00428     if( IAICQ.NbPoints() == 1 ) {
00429       gp_Pnt PIn = IAICQ.Point(1);
00430       const double preci = 1.e-10 * P.Distance(PC);
00431       // check if this point is internal for segment [PC,P]
00432       bool IsExternal =
00433         ( (PC.X()-PIn.X())*(P.X()-PIn.X()) > preci ) ||
00434         ( (PC.Y()-PIn.Y())*(P.Y()-PIn.Y()) > preci ) ||
00435         ( (PC.Z()-PIn.Z())*(P.Z()-PIn.Z()) > preci );
00436       if(IsExternal) {
00437         return false;
00438       }
00439       // check if this point is internal for triangle (P1,P2,P3)
00440       gp_Vec V1(PIn,P1);
00441       gp_Vec V2(PIn,P2);
00442       gp_Vec V3(PIn,P3);
00443       if( V1.Magnitude()<preci ||
00444           V2.Magnitude()<preci ||
00445           V3.Magnitude()<preci ) {
00446         Pint = PIn;
00447         return true;
00448       }
00449       const double angularTol = 1e-6;
00450       gp_Vec VC1 = V1.Crossed(V2);
00451       gp_Vec VC2 = V2.Crossed(V3);
00452       gp_Vec VC3 = V3.Crossed(V1);
00453       if(VC1.Magnitude()<gp::Resolution()) {
00454         if(VC2.IsOpposite(VC3,angularTol)) {
00455           return false;
00456         }
00457       }
00458       else if(VC2.Magnitude()<gp::Resolution()) {
00459         if(VC1.IsOpposite(VC3,angularTol)) {
00460           return false;
00461         }
00462       }
00463       else if(VC3.Magnitude()<gp::Resolution()) {
00464         if(VC1.IsOpposite(VC2,angularTol)) {
00465           return false;
00466         }
00467       }
00468       else {
00469         if( VC1.IsOpposite(VC2,angularTol) || VC1.IsOpposite(VC3,angularTol) ||
00470             VC2.IsOpposite(VC3,angularTol) ) {
00471           return false;
00472         }
00473       }
00474       Pint = PIn;
00475       return true;
00476     }
00477   }
00478 
00479   return false;
00480 }
00481 
00482 
00483 //=======================================================================
00484 //function : HasIntersection
00485 //purpose  : Auxilare for CheckIntersection()
00486 //=======================================================================
00487 
00488 static bool HasIntersection(const gp_Pnt& P, const gp_Pnt& PC, gp_Pnt& Pint,
00489                             Handle(TColgp_HSequenceOfPnt)& aContour)
00490 {
00491   if(aContour->Length()==3) {
00492     return HasIntersection3( P, PC, Pint, aContour->Value(1),
00493                              aContour->Value(2), aContour->Value(3) );
00494   }
00495   else {
00496     bool check = false;
00497     if( (aContour->Value(1).Distance(aContour->Value(2)) > 1.e-6) &&
00498         (aContour->Value(1).Distance(aContour->Value(3)) > 1.e-6) &&
00499         (aContour->Value(2).Distance(aContour->Value(3)) > 1.e-6) ) {
00500       check = HasIntersection3( P, PC, Pint, aContour->Value(1),
00501                                 aContour->Value(2), aContour->Value(3) );
00502     }
00503     if(check) return true;
00504     if( (aContour->Value(1).Distance(aContour->Value(4)) > 1.e-6) &&
00505         (aContour->Value(1).Distance(aContour->Value(3)) > 1.e-6) &&
00506         (aContour->Value(4).Distance(aContour->Value(3)) > 1.e-6) ) {
00507       check = HasIntersection3( P, PC, Pint, aContour->Value(1),
00508                                 aContour->Value(3), aContour->Value(4) );
00509     }
00510     if(check) return true;
00511   }
00512 
00513   return false;
00514 }
00515 
00516 //================================================================================
00527 //================================================================================
00528 
00529 bool StdMeshers_QuadToTriaAdaptor::CheckIntersection (const gp_Pnt&       P,
00530                                                       const gp_Pnt&       PC,
00531                                                       gp_Pnt&             Pint,
00532                                                       SMESH_Mesh&         aMesh,
00533                                                       const TopoDS_Shape& aShape,
00534                                                       const SMDS_MeshElement* NotCheckedFace)
00535 {
00536   if ( !myElemSearcher )
00537     myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher();
00538   SMESH_ElementSearcher* searcher = const_cast<SMESH_ElementSearcher*>(myElemSearcher);
00539 
00540   //SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
00541   //cout<<"    CheckIntersection: meshDS->NbFaces() = "<<meshDS->NbFaces()<<endl;
00542   bool res = false;
00543   double dist = RealLast(); // find intersection closest to the segment
00544   gp_Pnt Pres;
00545 
00546   gp_Ax1 line( P, gp_Vec(P,PC));
00547   vector< const SMDS_MeshElement* > suspectElems;
00548   searcher->GetElementsNearLine( line, SMDSAbs_Face, suspectElems);
00549   
00550   for ( int i = 0; i < suspectElems.size(); ++i )
00551   {
00552     const SMDS_MeshElement* face = suspectElems[i];
00553     if ( face == NotCheckedFace ) continue;
00554     Handle(TColgp_HSequenceOfPnt) aContour = new TColgp_HSequenceOfPnt;
00555     for ( int i = 0; i < face->NbCornerNodes(); ++i ) 
00556       aContour->Append( SMESH_TNodeXYZ( face->GetNode(i) ));
00557     if( HasIntersection(P, PC, Pres, aContour) ) {
00558       res = true;
00559       double tmp = PC.Distance(Pres);
00560       if(tmp<dist) {
00561         Pint = Pres;
00562         dist = tmp;
00563       }
00564     }
00565   }
00566   return res;
00567 }
00568 
00569 //================================================================================
00582 //================================================================================
00583 
00584 int StdMeshers_QuadToTriaAdaptor::Preparation(const SMDS_MeshElement*       face,
00585                                               Handle(TColgp_HArray1OfPnt)&  PN,
00586                                               Handle(TColgp_HArray1OfVec)&  VN,
00587                                               vector<const SMDS_MeshNode*>& FNodes,
00588                                               gp_Pnt&                       PC,
00589                                               gp_Vec&                       VNorm,
00590                                               const SMDS_MeshElement**      volumes)
00591 {
00592   if( face->NbCornerNodes() != 4 )
00593   {
00594     return NOT_QUAD;
00595   }
00596 
00597   int i = 0;
00598   gp_XYZ xyzC(0., 0., 0.);
00599   for ( i = 0; i < 4; ++i )
00600   {
00601     gp_XYZ p = SMESH_TNodeXYZ( FNodes[i] = face->GetNode(i) );
00602     PN->SetValue( i+1, p );
00603     xyzC += p;
00604   }
00605   PC = xyzC/4;
00606 
00607   int nbp = 4;
00608 
00609   int j = 0;
00610   for(i=1; i<4; i++) {
00611     j = i+1;
00612     for(; j<=4; j++) {
00613       if( PN->Value(i).Distance(PN->Value(j)) < 1.e-6 )
00614         break;
00615     }
00616     if(j<=4) break;
00617   }
00618   //int deg_num = IsDegenarate(PN);
00619   //if(deg_num>0) {
00620   bool hasdeg = false;
00621   if(i<4) {
00622     //cout<<"find degeneration"<<endl;
00623     hasdeg = true;
00624     gp_Pnt Pdeg = PN->Value(i);
00625 
00626     list< const SMDS_MeshNode* >::iterator itdg = myDegNodes.begin();
00627     const SMDS_MeshNode* DegNode = 0;
00628     for(; itdg!=myDegNodes.end(); itdg++) {
00629       const SMDS_MeshNode* N = (*itdg);
00630       gp_Pnt Ptmp(N->X(),N->Y(),N->Z());
00631       if(Pdeg.Distance(Ptmp)<1.e-6) {
00632         DegNode = N;
00633         //DegNode = const_cast<SMDS_MeshNode*>(N);
00634         break;
00635       }
00636     }
00637     if(!DegNode) {
00638       DegNode = FNodes[i-1];
00639       myDegNodes.push_back(DegNode);
00640     }
00641     else {
00642       FNodes[i-1] = DegNode;
00643     }
00644     for(i=j; i<4; i++) {
00645       PN->SetValue(i,PN->Value(i+1));
00646       FNodes[i-1] = FNodes[i];
00647     }
00648     nbp = 3;
00649   }
00650 
00651   PN->SetValue(nbp+1,PN->Value(1));
00652   FNodes[nbp] = FNodes[0];
00653   // find normal direction
00654   gp_Vec V1(PC,PN->Value(nbp));
00655   gp_Vec V2(PC,PN->Value(1));
00656   VNorm = V1.Crossed(V2);
00657   VN->SetValue(nbp,VNorm);
00658   for(i=1; i<nbp; i++) {
00659     V1 = gp_Vec(PC,PN->Value(i));
00660     V2 = gp_Vec(PC,PN->Value(i+1));
00661     gp_Vec Vtmp = V1.Crossed(V2);
00662     VN->SetValue(i,Vtmp);
00663     VNorm += Vtmp;
00664   }
00665 
00666   // find volumes sharing the face
00667   if ( volumes )
00668   {
00669     volumes[0] = volumes[1] = 0;
00670     SMDS_ElemIteratorPtr vIt = FNodes[0]->GetInverseElementIterator( SMDSAbs_Volume );
00671     while ( vIt->more() )
00672     {
00673       const SMDS_MeshElement* vol = vIt->next();
00674       bool volSharesAllNodes = true;
00675       for ( int i = 1; i < face->NbNodes() && volSharesAllNodes; ++i )
00676         volSharesAllNodes = ( vol->GetNodeIndex( FNodes[i] ) >= 0 );
00677       if ( volSharesAllNodes )
00678         volumes[ volumes[0] ? 1 : 0 ] = vol;
00679       // we could additionally check that vol has all FNodes in its one face using SMDS_VolumeTool
00680     }
00681     // define volume position relating to the face normal
00682     if ( volumes[0] )
00683     {
00684       // get volume gc
00685       SMDS_ElemIteratorPtr nodeIt = volumes[0]->nodesIterator();
00686       gp_XYZ volGC(0,0,0);
00687       volGC = accumulate( TXyzIterator(nodeIt), TXyzIterator(), volGC ) / volumes[0]->NbNodes();
00688 
00689       if ( VNorm * gp_Vec( PC, volGC ) < 0 )
00690         swap( volumes[0], volumes[1] );
00691     }
00692   }
00693 
00694   //cout<<"  VNorm("<<VNorm.X()<<","<<VNorm.Y()<<","<<VNorm.Z()<<")"<<endl;
00695   return hasdeg ? DEGEN_QUAD : QUAD;
00696 }
00697 
00698 
00699 //=======================================================================
00700 //function : Compute
00701 //purpose  : 
00702 //=======================================================================
00703 
00704 bool StdMeshers_QuadToTriaAdaptor::Compute(SMESH_Mesh&           aMesh,
00705                                            const TopoDS_Shape&   aShape,
00706                                            SMESH_ProxyMesh* aProxyMesh)
00707 {
00708   SMESH_ProxyMesh::setMesh( aMesh );
00709 
00710   if ( aShape.ShapeType() != TopAbs_SOLID &&
00711        aShape.ShapeType() != TopAbs_SHELL )
00712     return false;
00713 
00714   vector<const SMDS_MeshElement*> myPyramids;
00715 
00716   SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
00717   SMESH_MesherHelper helper(aMesh);
00718   helper.IsQuadraticSubMesh(aShape);
00719   helper.SetElementsOnShape( true );
00720 
00721   if ( myElemSearcher ) delete myElemSearcher;
00722   if ( aProxyMesh )
00723     myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher( aProxyMesh->GetFaces(aShape));
00724   else
00725     myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher();
00726 
00727   const SMESHDS_SubMesh * aSubMeshDSFace;
00728   Handle(TColgp_HArray1OfPnt) PN = new TColgp_HArray1OfPnt(1,5);
00729   Handle(TColgp_HArray1OfVec) VN = new TColgp_HArray1OfVec(1,4);
00730   vector<const SMDS_MeshNode*> FNodes(5);
00731   gp_Pnt PC;
00732   gp_Vec VNorm;
00733 
00734   for (TopExp_Explorer exp(aShape,TopAbs_FACE);exp.More();exp.Next())
00735   {
00736     const TopoDS_Shape& aShapeFace = exp.Current();
00737     if ( aProxyMesh )
00738       aSubMeshDSFace = aProxyMesh->GetSubMesh( aShapeFace );
00739     else
00740       aSubMeshDSFace = meshDS->MeshElements( aShapeFace );
00741 
00742     vector<const SMDS_MeshElement*> trias, quads;
00743     bool hasNewTrias = false;
00744 
00745     if ( aSubMeshDSFace )
00746     {
00747       bool isRev = false;
00748       if ( helper.NbAncestors( aShapeFace, aMesh, aShape.ShapeType() ) > 1 )
00749         isRev = SMESH_Algo::IsReversedSubMesh( TopoDS::Face(aShapeFace), meshDS );
00750 
00751       SMDS_ElemIteratorPtr iteratorElem = aSubMeshDSFace->GetElements();
00752       while ( iteratorElem->more() ) // loop on elements on a geometrical face
00753       {
00754         const SMDS_MeshElement* face = iteratorElem->next();
00755         // preparation step to get face info
00756         int stat = Preparation(face, PN, VN, FNodes, PC, VNorm);
00757         switch ( stat )
00758         {
00759         case NOT_QUAD:
00760 
00761           trias.push_back( face );
00762           break;
00763 
00764         case DEGEN_QUAD:
00765           {
00766             // degenerate face
00767             // add triangles to result map
00768             SMDS_MeshFace* NewFace;
00769             if(!isRev)
00770               NewFace = meshDS->AddFace( FNodes[0], FNodes[1], FNodes[2] );
00771             else
00772               NewFace = meshDS->AddFace( FNodes[0], FNodes[2], FNodes[1] );
00773             storeTmpElement( NewFace );
00774             trias.push_back ( NewFace );
00775             quads.push_back( face );
00776             hasNewTrias = true;
00777             break;
00778           }
00779 
00780         case QUAD:
00781           {
00782             if(!isRev) VNorm.Reverse();
00783             double xc = 0., yc = 0., zc = 0.;
00784             int i = 1;
00785             for(; i<=4; i++) {
00786               gp_Pnt Pbest;
00787               if(!isRev)
00788                 Pbest = FindBestPoint(PN->Value(i), PN->Value(i+1), PC, VN->Value(i).Reversed());
00789               else
00790                 Pbest = FindBestPoint(PN->Value(i), PN->Value(i+1), PC, VN->Value(i));
00791               xc += Pbest.X();
00792               yc += Pbest.Y();
00793               zc += Pbest.Z();
00794             }
00795             gp_Pnt PCbest(xc/4., yc/4., zc/4.);
00796 
00797             // check PCbest
00798             double height = PCbest.Distance(PC);
00799             if(height<1.e-6) {
00800               // create new PCbest using a bit shift along VNorm
00801               PCbest = PC.XYZ() + VNorm.XYZ() * 0.001;
00802             }
00803             else {
00804               // check possible intersection with other faces
00805               gp_Pnt Pint;
00806               bool check = CheckIntersection(PCbest, PC, Pint, aMesh, aShape, face);
00807               if(check) {
00808                 //cout<<"--PC("<<PC.X()<<","<<PC.Y()<<","<<PC.Z()<<")"<<endl;
00809                 //cout<<"  PCbest("<<PCbest.X()<<","<<PCbest.Y()<<","<<PCbest.Z()<<")"<<endl;
00810                 double dist = PC.Distance(Pint)/3.;
00811                 gp_Dir aDir(gp_Vec(PC,PCbest));
00812                 PCbest = PC.XYZ() + aDir.XYZ() * dist;
00813               }
00814               else {
00815                 gp_Vec VB(PC,PCbest);
00816                 gp_Pnt PCbestTmp = PC.XYZ() + VB.XYZ() * 3.0;
00817                 check = CheckIntersection(PCbestTmp, PC, Pint, aMesh, aShape, face);
00818                 if(check) {
00819                   double dist = PC.Distance(Pint)/3.;
00820                   if(dist<height) {
00821                     gp_Dir aDir(gp_Vec(PC,PCbest));
00822                     PCbest = PC.XYZ() + aDir.XYZ() * dist;
00823                   }
00824                 }
00825               }
00826             }
00827             // create node for PCbest
00828             SMDS_MeshNode* NewNode = helper.AddNode( PCbest.X(), PCbest.Y(), PCbest.Z() );
00829 
00830             // add triangles to result map
00831             for(i=0; i<4; i++)
00832             {
00833               trias.push_back ( meshDS->AddFace( NewNode, FNodes[i], FNodes[i+1] ));
00834               storeTmpElement( trias.back() );
00835             }
00836             // create a pyramid
00837             if ( isRev ) swap( FNodes[1], FNodes[3]);
00838             SMDS_MeshVolume* aPyram =
00839               helper.AddVolume( FNodes[0], FNodes[1], FNodes[2], FNodes[3], NewNode );
00840             myPyramids.push_back(aPyram);
00841 
00842             quads.push_back( face );
00843             hasNewTrias = true;
00844             break;
00845 
00846           } // case QUAD:
00847 
00848         } // switch ( stat )
00849       } // end loop on elements on a face submesh
00850 
00851       bool sourceSubMeshIsProxy = false;
00852       if ( aProxyMesh )
00853       {
00854         // move proxy sub-mesh from other proxy mesh to this
00855         sourceSubMeshIsProxy = takeProxySubMesh( aShapeFace, aProxyMesh );
00856         // move also tmp elements added in mesh
00857         takeTmpElemsInMesh( aProxyMesh );
00858       }
00859       if ( hasNewTrias )
00860       {
00861         SMESH_ProxyMesh::SubMesh* prxSubMesh = getProxySubMesh( aShapeFace );
00862         prxSubMesh->ChangeElements( trias.begin(), trias.end() );
00863 
00864         // delete tmp quadrangles removed from aProxyMesh
00865         if ( sourceSubMeshIsProxy )
00866         {
00867           for ( unsigned i = 0; i < quads.size(); ++i )
00868             removeTmpElement( quads[i] );
00869 
00870           delete myElemSearcher;
00871           myElemSearcher =
00872             SMESH_MeshEditor(&aMesh).GetElementSearcher( aProxyMesh->GetFaces(aShape));
00873         }
00874       }
00875     }
00876   } // end for(TopExp_Explorer exp(aShape,TopAbs_FACE);exp.More();exp.Next()) {
00877 
00878   return Compute2ndPart(aMesh, myPyramids);
00879 }
00880 
00881 //================================================================================
00885 //================================================================================
00886 
00887 bool StdMeshers_QuadToTriaAdaptor::Compute(SMESH_Mesh& aMesh)
00888 {
00889   SMESH_ProxyMesh::setMesh( aMesh );
00890   SMESH_ProxyMesh::_allowedTypes.push_back( SMDSEntity_Triangle );
00891   SMESH_ProxyMesh::_allowedTypes.push_back( SMDSEntity_Quad_Triangle );
00892   if ( aMesh.NbQuadrangles() < 1 )
00893     return false;
00894 
00895   vector<const SMDS_MeshElement*> myPyramids;
00896   SMESH_MesherHelper helper(aMesh);
00897   helper.IsQuadraticSubMesh(aMesh.GetShapeToMesh());
00898   helper.SetElementsOnShape( true );
00899 
00900   if ( !myElemSearcher )
00901     myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher();
00902   SMESH_ElementSearcher* searcher = const_cast<SMESH_ElementSearcher*>(myElemSearcher);
00903 
00904   SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
00905   SMESH_ProxyMesh::SubMesh* prxSubMesh = getProxySubMesh();
00906 
00907   SMDS_FaceIteratorPtr fIt = meshDS->facesIterator(/*idInceasingOrder=*/true);
00908   while( fIt->more()) 
00909   {
00910     const SMDS_MeshElement* face = fIt->next();
00911     if ( !face ) continue;
00912     // retrieve needed information about a face
00913     Handle(TColgp_HArray1OfPnt) PN = new TColgp_HArray1OfPnt(1,5);
00914     Handle(TColgp_HArray1OfVec) VN = new TColgp_HArray1OfVec(1,4);
00915     vector<const SMDS_MeshNode*> FNodes(5);
00916     gp_Pnt PC;
00917     gp_Vec VNorm;
00918     const SMDS_MeshElement* volumes[2];
00919     int what = Preparation(face, PN, VN, FNodes, PC, VNorm, volumes);
00920     if ( what == NOT_QUAD )
00921       continue;
00922     if ( volumes[0] && volumes[1] )
00923       continue; // face is shared by two volumes - no space for a pyramid
00924 
00925     if ( what == DEGEN_QUAD )
00926     {
00927       // degenerate face
00928       // add a triangle to the proxy mesh
00929       SMDS_MeshFace* NewFace;
00930 
00931       // check orientation
00932       double tmp = PN->Value(1).Distance(PN->Value(2)) + PN->Value(2).Distance(PN->Value(3));
00933       // far points in VNorm direction
00934       gp_Pnt Ptmp1 = PC.XYZ() + VNorm.XYZ() * tmp * 1.e6;
00935       gp_Pnt Ptmp2 = PC.XYZ() - VNorm.XYZ() * tmp * 1.e6;
00936       // check intersection for Ptmp1 and Ptmp2
00937       bool IsRev = false;
00938       bool IsOK1 = false;
00939       bool IsOK2 = false;
00940       double dist1 = RealLast();
00941       double dist2 = RealLast();
00942       gp_Pnt Pres1,Pres2;
00943 
00944       gp_Ax1 line( PC, VNorm );
00945       vector< const SMDS_MeshElement* > suspectElems;
00946       searcher->GetElementsNearLine( line, SMDSAbs_Face, suspectElems);
00947 
00948       for ( int iF = 0; iF < suspectElems.size(); ++iF ) {
00949         const SMDS_MeshElement* F = suspectElems[iF];
00950         if(F==face) continue;
00951         Handle(TColgp_HSequenceOfPnt) aContour = new TColgp_HSequenceOfPnt;
00952         for ( int i = 0; i < 4; ++i )
00953           aContour->Append( SMESH_TNodeXYZ( F->GetNode(i) ));
00954         gp_Pnt PPP;
00955         if( !volumes[0] && HasIntersection(Ptmp1, PC, PPP, aContour) ) {
00956           IsOK1 = true;
00957           double tmp = PC.Distance(PPP);
00958           if(tmp<dist1) {
00959             Pres1 = PPP;
00960             dist1 = tmp;
00961           }
00962         }
00963         if( !volumes[1] && HasIntersection(Ptmp2, PC, PPP, aContour) ) {
00964           IsOK2 = true;
00965           double tmp = PC.Distance(PPP);
00966           if(tmp<dist2) {
00967             Pres2 = PPP;
00968             dist2 = tmp;
00969           }
00970         }
00971       }
00972 
00973       if( IsOK1 && !IsOK2 ) {
00974         // using existed direction
00975       }
00976       else if( !IsOK1 && IsOK2 ) {
00977         // using opposite direction
00978         IsRev = true;
00979       }
00980       else { // IsOK1 && IsOK2
00981         double tmp1 = PC.Distance(Pres1);
00982         double tmp2 = PC.Distance(Pres2);
00983         if(tmp1<tmp2) {
00984           // using existed direction
00985         }
00986         else {
00987           // using opposite direction
00988           IsRev = true;
00989         }
00990       }
00991       if(!IsRev)
00992         NewFace = meshDS->AddFace( FNodes[0], FNodes[1], FNodes[2] );
00993       else
00994         NewFace = meshDS->AddFace( FNodes[0], FNodes[2], FNodes[1] );
00995       storeTmpElement( NewFace );
00996       prxSubMesh->AddElement( NewFace );
00997       continue;
00998     }
00999 
01000     // Case of non-degenerated quadrangle 
01001 
01002     // Find pyramid peak
01003 
01004     gp_XYZ PCbest(0., 0., 0.); // pyramid peak
01005     int i = 1;
01006     for(; i<=4; i++) {
01007       gp_Pnt Pbest = FindBestPoint(PN->Value(i), PN->Value(i+1), PC, VN->Value(i));
01008       PCbest += Pbest.XYZ();
01009     }
01010     PCbest /= 4;
01011 
01012     double height = PC.Distance(PCbest); // pyramid height to precise
01013     if(height<1.e-6) {
01014       // create new PCbest using a bit shift along VNorm
01015       PCbest = PC.XYZ() + VNorm.XYZ() * 0.001;
01016       height = PC.Distance(PCbest);
01017     }
01018 
01019     // Restrict pyramid height by intersection with other faces
01020     gp_Vec tmpDir(PC,PCbest); tmpDir.Normalize();
01021     double tmp = PN->Value(1).Distance(PN->Value(3)) + PN->Value(2).Distance(PN->Value(4));
01022     // far points: in (PC, PCbest) direction and vice-versa
01023     gp_Pnt farPnt[2] = { PC.XYZ() + tmpDir.XYZ() * tmp * 1.e6,
01024                          PC.XYZ() - tmpDir.XYZ() * tmp * 1.e6 };
01025     // check intersection for farPnt1 and farPnt2
01026     bool   intersected[2] = { false, false };
01027     double dist       [2] = { RealLast(), RealLast() };
01028     gp_Pnt intPnt[2];
01029 
01030     gp_Ax1 line( PC, tmpDir );
01031     vector< const SMDS_MeshElement* > suspectElems;
01032     searcher->GetElementsNearLine( line, SMDSAbs_Face, suspectElems);
01033 
01034     for ( int iF = 0; iF < suspectElems.size(); ++iF )
01035     {
01036       const SMDS_MeshElement* F = suspectElems[iF];
01037       if(F==face) continue;
01038       Handle(TColgp_HSequenceOfPnt) aContour = new TColgp_HSequenceOfPnt;
01039       int nbN = F->NbNodes() / ( F->IsQuadratic() ? 2 : 1 );
01040       for ( i = 0; i < nbN; ++i )
01041         aContour->Append( SMESH_TNodeXYZ( F->GetNode(i) ));
01042       gp_Pnt intP;
01043       for ( int isRev = 0; isRev < 2; ++isRev )
01044       {
01045         if( !volumes[isRev] && HasIntersection(farPnt[isRev], PC, intP, aContour) ) {
01046           intersected[isRev] = true;
01047           double d = PC.Distance( intP );
01048           if( d < dist[isRev] )
01049           {
01050             intPnt[isRev] = intP;
01051             dist  [isRev] = d;
01052           }
01053         }
01054       }
01055     }
01056 
01057     // Create one or two pyramids
01058 
01059     for ( int isRev = 0; isRev < 2; ++isRev )
01060     {
01061       if( !intersected[isRev] ) continue;
01062       double pyramidH = Min( height, PC.Distance(intPnt[isRev])/3.);
01063       PCbest = PC.XYZ() + tmpDir.XYZ() * (isRev ? -pyramidH : pyramidH);
01064 
01065       // create node for PCbest
01066       SMDS_MeshNode* NewNode = helper.AddNode( PCbest.X(), PCbest.Y(), PCbest.Z() );
01067 
01068       // add triangles to result map
01069       for(i=0; i<4; i++) {
01070         SMDS_MeshFace* NewFace;
01071         if(isRev)
01072           NewFace = meshDS->AddFace( NewNode, FNodes[i], FNodes[i+1] );
01073         else
01074           NewFace = meshDS->AddFace( NewNode, FNodes[i+1], FNodes[i] );
01075         storeTmpElement( NewFace );
01076         prxSubMesh->AddElement( NewFace );
01077       }
01078       // create a pyramid
01079       SMDS_MeshVolume* aPyram;
01080       if(isRev)
01081         aPyram = helper.AddVolume( FNodes[0], FNodes[1], FNodes[2], FNodes[3], NewNode );
01082       else
01083         aPyram = helper.AddVolume( FNodes[0], FNodes[3], FNodes[2], FNodes[1], NewNode );
01084       myPyramids.push_back(aPyram);
01085     }
01086   } // end loop on all faces
01087 
01088   return Compute2ndPart(aMesh, myPyramids);
01089 }
01090 
01091 //================================================================================
01095 //================================================================================
01096 
01097 bool StdMeshers_QuadToTriaAdaptor::Compute2ndPart(SMESH_Mesh&                            aMesh,
01098                                                   const vector<const SMDS_MeshElement*>& myPyramids)
01099 {
01100   if(myPyramids.empty())
01101     return true;
01102 
01103   SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
01104   int i, j, k, myShapeID = myPyramids[0]->GetNode(4)->getshapeId();
01105 
01106   if ( myElemSearcher ) delete myElemSearcher;
01107   myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher();
01108   SMESH_ElementSearcher* searcher = const_cast<SMESH_ElementSearcher*>(myElemSearcher);
01109 
01110   set<const SMDS_MeshNode*> nodesToMove;
01111 
01112   // check adjacent pyramids
01113 
01114   for ( i = 0; i <  myPyramids.size(); ++i )
01115   {
01116     const SMDS_MeshElement* PrmI = myPyramids[i];
01117     MergeAdjacent( PrmI, nodesToMove );
01118   }
01119 
01120   // iterate on all pyramids
01121   for ( i = 0; i <  myPyramids.size(); ++i )
01122   {
01123     const SMDS_MeshElement* PrmI = myPyramids[i];
01124 
01125     // compare PrmI with all the rest pyramids
01126 
01127     // collect adjacent pyramids and nodes coordinates of PrmI
01128     set<const SMDS_MeshElement*> checkedPyrams;
01129     vector<gp_Pnt> PsI(5);
01130     for(k=0; k<5; k++) // loop on 4 base nodes of PrmI
01131     {
01132       const SMDS_MeshNode* n = PrmI->GetNode(k);
01133       PsI[k] = SMESH_TNodeXYZ( n );
01134       SMDS_ElemIteratorPtr vIt = n->GetInverseElementIterator( SMDSAbs_Volume );
01135       while ( vIt->more() )
01136       {
01137         const SMDS_MeshElement* PrmJ = vIt->next();
01138         if ( SMESH_Algo::GetCommonNodes( PrmI, PrmJ ).size() > 1 )
01139           checkedPyrams.insert( PrmJ );
01140       }
01141     }
01142 
01143     // check intersection with distant pyramids
01144     for(k=0; k<4; k++) // loop on 4 base nodes of PrmI
01145     {
01146       gp_Vec Vtmp(PsI[k],PsI[4]);
01147       gp_Ax1 line( PsI[k], Vtmp );
01148       vector< const SMDS_MeshElement* > suspectPyrams;
01149       searcher->GetElementsNearLine( line, SMDSAbs_Volume, suspectPyrams);
01150 
01151       for ( j = 0; j < suspectPyrams.size(); ++j )
01152       {
01153         const SMDS_MeshElement* PrmJ = suspectPyrams[j];
01154         if ( PrmJ == PrmI || PrmJ->NbCornerNodes() != 5 )
01155           continue;
01156         if ( myShapeID != PrmJ->GetNode(4)->getshapeId())
01157           continue; // pyramid from other SOLID
01158         if ( PrmI->GetNode(4) == PrmJ->GetNode(4) )
01159           continue; // pyramids PrmI and PrmJ already merged
01160         if ( !checkedPyrams.insert( PrmJ ).second )
01161           continue; // already checked
01162 
01163         TXyzIterator xyzIt( PrmJ->nodesIterator() );
01164         vector<gp_Pnt> PsJ( xyzIt, TXyzIterator() );
01165 
01166         gp_Pnt Pint;
01167         bool hasInt=false;
01168         for(k=0; k<4 && !hasInt; k++) {
01169           gp_Vec Vtmp(PsI[k],PsI[4]);
01170           gp_Pnt Pshift = PsI[k].XYZ() + Vtmp.XYZ() * 0.01; // base node moved a bit to apex
01171           hasInt = 
01172           ( HasIntersection3( Pshift, PsI[4], Pint, PsJ[0], PsJ[1], PsJ[4]) ||
01173             HasIntersection3( Pshift, PsI[4], Pint, PsJ[1], PsJ[2], PsJ[4]) ||
01174             HasIntersection3( Pshift, PsI[4], Pint, PsJ[2], PsJ[3], PsJ[4]) ||
01175             HasIntersection3( Pshift, PsI[4], Pint, PsJ[3], PsJ[0], PsJ[4]) );
01176         }
01177         for(k=0; k<4 && !hasInt; k++) {
01178           gp_Vec Vtmp(PsJ[k],PsJ[4]);
01179           gp_Pnt Pshift = PsJ[k].XYZ() + Vtmp.XYZ() * 0.01;
01180           hasInt = 
01181             ( HasIntersection3( Pshift, PsJ[4], Pint, PsI[0], PsI[1], PsI[4]) ||
01182               HasIntersection3( Pshift, PsJ[4], Pint, PsI[1], PsI[2], PsI[4]) ||
01183               HasIntersection3( Pshift, PsJ[4], Pint, PsI[2], PsI[3], PsI[4]) ||
01184               HasIntersection3( Pshift, PsJ[4], Pint, PsI[3], PsI[0], PsI[4]) );
01185         }
01186 
01187         if ( hasInt )
01188         {
01189           // count common nodes of base faces of two pyramids
01190           int nbc = 0;
01191           for (k=0; k<4; k++)
01192             nbc += int ( PrmI->GetNodeIndex( PrmJ->GetNode(k) ) >= 0 );
01193 
01194           if ( nbc == 4 )
01195             continue; // pyrams have a common base face
01196 
01197           if(nbc>0)
01198           {
01199             // Merge the two pyramids and others already merged with them
01200             MergePiramids( PrmI, PrmJ, nodesToMove );
01201           }
01202           else { // nbc==0
01203 
01204             // decrease height of pyramids
01205             gp_XYZ PCi(0,0,0), PCj(0,0,0);
01206             for(k=0; k<4; k++) {
01207               PCi += PsI[k].XYZ();
01208               PCj += PsJ[k].XYZ();
01209             }
01210             PCi /= 4; PCj /= 4; 
01211             gp_Vec VN1(PCi,PsI[4]);
01212             gp_Vec VN2(PCj,PsJ[4]);
01213             gp_Vec VI1(PCi,Pint);
01214             gp_Vec VI2(PCj,Pint);
01215             double ang1 = fabs(VN1.Angle(VI1));
01216             double ang2 = fabs(VN2.Angle(VI2));
01217             double coef1 = 0.5 - (( ang1<PI/3 ) ? cos(ang1)*0.25 : 0 );
01218             double coef2 = 0.5 - (( ang2<PI/3 ) ? cos(ang2)*0.25 : 0 ); // cos(ang2) ?
01219 //             double coef2 = 0.5;
01220 //             if(ang2<PI/3)
01221 //               coef2 -= cos(ang1)*0.25;
01222 
01223             VN1.Scale(coef1);
01224             VN2.Scale(coef2);
01225             SMDS_MeshNode* aNode1 = const_cast<SMDS_MeshNode*>(PrmI->GetNode(4));
01226             aNode1->setXYZ( PCi.X()+VN1.X(), PCi.Y()+VN1.Y(), PCi.Z()+VN1.Z() );
01227             SMDS_MeshNode* aNode2 = const_cast<SMDS_MeshNode*>(PrmJ->GetNode(4));
01228             aNode2->setXYZ( PCj.X()+VN2.X(), PCj.Y()+VN2.Y(), PCj.Z()+VN2.Z() );
01229             nodesToMove.insert( aNode1 );
01230             nodesToMove.insert( aNode2 );
01231           }
01232           // fix intersections that could appear after apex movement
01233           MergeAdjacent( PrmI, nodesToMove );
01234           MergeAdjacent( PrmJ, nodesToMove );
01235 
01236         } // end if(hasInt)
01237       } // loop on suspectPyrams
01238     }  // loop on 4 base nodes of PrmI
01239 
01240   } // loop on all pyramids
01241 
01242   if( !nodesToMove.empty() && !meshDS->IsEmbeddedMode() )
01243   {
01244     set<const SMDS_MeshNode*>::iterator n = nodesToMove.begin();
01245     for ( ; n != nodesToMove.end(); ++n )
01246       meshDS->MoveNode( *n, (*n)->X(), (*n)->Y(), (*n)->Z() );
01247   }
01248 
01249   // move medium nodes of merged quadratic pyramids
01250   if ( myPyramids[0]->IsQuadratic() )
01251     UpdateQuadraticPyramids( nodesToMove, GetMeshDS() );
01252 
01253   // erase removed triangles from the proxy mesh
01254   if ( !myRemovedTrias.empty() )
01255   {
01256     for ( int i = 0; i <= meshDS->MaxShapeIndex(); ++i )
01257       if ( SMESH_ProxyMesh::SubMesh* sm = findProxySubMesh(i))
01258       {
01259         vector<const SMDS_MeshElement *> faces;
01260         faces.reserve( sm->NbElements() );
01261         SMDS_ElemIteratorPtr fIt = sm->GetElements();
01262         while ( fIt->more() )
01263         {
01264           const SMDS_MeshElement* tria = fIt->next();
01265           set<const SMDS_MeshElement*>::iterator rmTria = myRemovedTrias.find( tria );
01266           if ( rmTria != myRemovedTrias.end() )
01267             myRemovedTrias.erase( rmTria );
01268           else
01269             faces.push_back( tria );
01270         }
01271         sm->ChangeElements( faces.begin(), faces.end() );
01272       }
01273   }
01274 
01275   myDegNodes.clear();
01276 
01277   delete myElemSearcher;
01278   myElemSearcher=0;
01279 
01280   return true;
01281 }
01282 
Copyright © 2007-2011 CEA/DEN, EDF R&D, OPEN CASCADE
Copyright © 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN, CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS