Salome HOME
Update copyright information
[modules/smesh.git] / src / StdMeshers / StdMeshers_QuadToTriaAdaptor.cxx
index c2a54ebb0cb061cbefd53da304f0dafaa58c1314..008efb414cd46f81473160e97211e945b54466bd 100644 (file)
@@ -1,28 +1,27 @@
-//  Copyright (C) 2007-2010  CEA/DEN, EDF R&D, OPEN CASCADE
+// Copyright (C) 2007-2012  CEA/DEN, EDF R&D, OPEN CASCADE
 //
-//  This library is free software; you can redistribute it and/or
-//  modify it under the terms of the GNU Lesser General Public
-//  License as published by the Free Software Foundation; either
-//  version 2.1 of the License.
+// This library is free software; you can redistribute it and/or
+// modify it under the terms of the GNU Lesser General Public
+// License as published by the Free Software Foundation; either
+// version 2.1 of the License.
 //
-//  This library is distributed in the hope that it will be useful,
-//  but WITHOUT ANY WARRANTY; without even the implied warranty of
-//  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
-//  Lesser General Public License for more details.
+// This library is distributed in the hope that it will be useful,
+// but WITHOUT ANY WARRANTY; without even the implied warranty of
+// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+// Lesser General Public License for more details.
 //
-//  You should have received a copy of the GNU Lesser General Public
-//  License along with this library; if not, write to the Free Software
-//  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
+// You should have received a copy of the GNU Lesser General Public
+// License along with this library; if not, write to the Free Software
+// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
 //
-//  See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
+// See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
 //
 
-//  SMESH SMESH : implementaion of SMESH idl descriptions
 // File      : StdMeshers_QuadToTriaAdaptor.cxx
 // Module    : SMESH
 // Created   : Wen May 07 16:37:07 2008
 // Author    : Sergey KUUL (skl)
-//
+
 #include "StdMeshers_QuadToTriaAdaptor.hxx"
 
 #include "SMDS_SetIterator.hxx"
@@ -91,7 +90,7 @@ namespace
         int ind = baseNodes[0] ? 1:0;
         if ( baseNodes[ ind ])
           return false; // pyramids with a common base face
-        baseNodes   [ ind ] = PrmI->GetNode(i);
+        baseNodes    [ ind ] = PrmI->GetNode(i);
         baseNodesIndI[ ind ] = i;
         baseNodesIndJ[ ind ] = j;
       }
@@ -111,15 +110,14 @@ namespace
 
     // Check angle between normals
     double angle = nI.Angle( nJ );
-    bool tooClose = ( angle < 15 * PI180 );
+    bool tooClose = ( angle < 15. * M_PI / 180. );
 
     // Check if pyramids collide
-    bool isOutI, isOutJ;
     if ( !tooClose && baI * baJ > 0 )
     {
       // find out if nI points outside of PrmI or inside
       int dInd = baseNodesIndI[1] - baseNodesIndI[0];
-      isOutI = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0;
+      bool isOutI = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0;
 
       // find out sign of projection of nJ to baI
       double proj = baI * nJ;
@@ -133,12 +131,15 @@ namespace
       // check order of baseNodes within pyramids, it must be opposite
       int dInd;
       dInd = baseNodesIndI[1] - baseNodesIndI[0];
-      isOutI = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0;
+      bool isOutI = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0;
       dInd = baseNodesIndJ[1] - baseNodesIndJ[0];
-      isOutJ = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0;
+      bool isOutJ = ( abs(dInd)==1 ) ? dInd < 0 : dInd > 0;
       if ( isOutJ == isOutI )
         return false; // other domain
 
+      // direct both normals outside pyramid
+      ( isOutI ? nJ : nI ).Reverse();
+
       // check absence of a face separating domains between pyramids
       TIDSortedElemSet emptySet, avoidSet;
       int i1, i2;
@@ -153,6 +154,9 @@ namespace
         while ( otherNodeInd == i1 || otherNodeInd == i2 ) otherNodeInd++;
         const SMDS_MeshNode* otherFaceNode = f->GetNode( otherNodeInd );
 
+        if ( otherFaceNode == nApexI || otherFaceNode == nApexJ )
+          continue; // f is a temporary triangle
+
         // check if f is a base face of either of pyramids
         if ( f->NbCornerNodes() == 4 &&
              ( PrmI->GetNodeIndex( otherFaceNode ) >= 0 ||
@@ -161,7 +165,6 @@ namespace
 
         // check projections of face direction (baOFN) to triange normals (nI and nJ)
         gp_Vec baOFN( base1, SMESH_TNodeXYZ( otherFaceNode ));
-        ( isOutI ? nJ : nI ).Reverse();
         if ( nI * baOFN > 0 && nJ * baOFN > 0 )
         {
           tooClose = false; // f is between pyramids
@@ -172,6 +175,79 @@ namespace
 
     return tooClose;
   }
+
+  //================================================================================
+  /*!
+   * \brief Move medium nodes of merged quadratic pyramids
+   */
+  //================================================================================
+
+  void UpdateQuadraticPyramids(const set<const SMDS_MeshNode*>& commonApex,
+                               SMESHDS_Mesh*                    meshDS)
+  {
+    typedef SMDS_StdIterator< const SMDS_MeshElement*, SMDS_ElemIteratorPtr > TStdElemIterator;
+    TStdElemIterator itEnd;
+
+    // shift of node index to get medium nodes between the 4 base nodes and the apex
+    const int base2MediumShift = 9;
+
+    set<const SMDS_MeshNode*>::const_iterator nIt = commonApex.begin();
+    for ( ; nIt != commonApex.end(); ++nIt )
+    {
+      SMESH_TNodeXYZ apex( *nIt );
+
+      vector< const SMDS_MeshElement* > pyrams // pyramids sharing the apex node
+        ( TStdElemIterator( apex._node->GetInverseElementIterator( SMDSAbs_Volume )), itEnd );
+
+      // Select medium nodes to keep and medium nodes to remove
+
+      typedef map < const SMDS_MeshNode*, const SMDS_MeshNode*, TIDCompare > TN2NMap;
+      TN2NMap base2medium; // to keep
+      vector< const SMDS_MeshNode* > nodesToRemove;
+
+      for ( unsigned i = 0; i < pyrams.size(); ++i )
+        for ( int baseIndex = 0; baseIndex < PYRAM_APEX; ++baseIndex )
+        {
+          SMESH_TNodeXYZ         base = pyrams[i]->GetNode( baseIndex );
+          const SMDS_MeshNode* medium = pyrams[i]->GetNode( baseIndex + base2MediumShift );
+          TN2NMap::iterator b2m = base2medium.insert( make_pair( base._node, medium )).first;
+          if ( b2m->second != medium )
+          {
+            nodesToRemove.push_back( medium );
+          }
+          else
+          {
+            // move the kept medium node
+            gp_XYZ newXYZ = 0.5 * ( apex + base );
+            meshDS->MoveNode( medium, newXYZ.X(), newXYZ.Y(), newXYZ.Z() );
+          }
+        }
+
+      // Within pyramids, replace nodes to remove by nodes to keep  
+
+      for ( unsigned i = 0; i < pyrams.size(); ++i )
+      {
+        vector< const SMDS_MeshNode* > nodes( pyrams[i]->begin_nodes(),
+                                              pyrams[i]->end_nodes() );
+        for ( int baseIndex = 0; baseIndex < PYRAM_APEX; ++baseIndex )
+        {
+          const SMDS_MeshNode* base = pyrams[i]->GetNode( baseIndex );
+          nodes[ baseIndex + base2MediumShift ] = base2medium[ base ];
+        }
+        meshDS->ChangeElementNodes( pyrams[i], &nodes[0], nodes.size());
+      }
+
+      // Remove the replaced nodes
+
+      if ( !nodesToRemove.empty() )
+      {
+        SMESHDS_SubMesh * sm = meshDS->MeshElements( nodesToRemove[0]->getshapeId() );
+        for ( unsigned i = 0; i < nodesToRemove.size(); ++i )
+          meshDS->RemoveFreeNode( nodesToRemove[i], sm, /*fromGroups=*/false);
+      }
+    }
+  }
+
 }
 
 //================================================================================
@@ -185,15 +261,15 @@ void StdMeshers_QuadToTriaAdaptor::MergePiramids( const SMDS_MeshElement*     Pr
                                                   set<const SMDS_MeshNode*> & nodesToMove)
 {
   const SMDS_MeshNode* Nrem = PrmJ->GetNode(4); // node to remove
-  int nbJ = Nrem->NbInverseElements( SMDSAbs_Volume );
+  //int nbJ = Nrem->NbInverseElements( SMDSAbs_Volume );
   SMESH_TNodeXYZ Pj( Nrem );
 
   // an apex node to make common to all merged pyramids
   SMDS_MeshNode* CommonNode = const_cast<SMDS_MeshNode*>(PrmI->GetNode(4));
   if ( CommonNode == Nrem ) return; // already merged
-  int nbI = CommonNode->NbInverseElements( SMDSAbs_Volume );
+  //int nbI = CommonNode->NbInverseElements( SMDSAbs_Volume );
   SMESH_TNodeXYZ Pi( CommonNode );
-  gp_XYZ Pnew = ( nbI*Pi + nbJ*Pj ) / (nbI+nbJ);
+  gp_XYZ Pnew = /*( nbI*Pi + nbJ*Pj ) / (nbI+nbJ);*/ 0.5 * ( Pi + Pj );
   CommonNode->setXYZ( Pnew.X(), Pnew.Y(), Pnew.Z() );
 
   nodesToMove.insert( CommonNode );
@@ -204,7 +280,7 @@ void StdMeshers_QuadToTriaAdaptor::MergePiramids( const SMDS_MeshElement*     Pr
 
   // find and remove coincided faces of merged pyramids
   vector< const SMDS_MeshElement* > inverseElems
-    // copy inverse elements to avoid iteration on changing conainer 
+    // copy inverse elements to avoid iteration on changing container 
     ( TStdElemIterator( CommonNode->GetInverseElementIterator(SMDSAbs_Face)), itEnd);
   for ( unsigned i = 0; i < inverseElems.size(); ++i )
   {
@@ -1027,8 +1103,8 @@ bool StdMeshers_QuadToTriaAdaptor::Compute2ndPart(SMESH_Mesh&
   SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
   int i, j, k, myShapeID = myPyramids[0]->GetNode(4)->getshapeId();
 
-  if ( !myElemSearcher )
-    myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher();
+  if ( myElemSearcher ) delete myElemSearcher;
+  myElemSearcher = SMESH_MeshEditor(&aMesh).GetElementSearcher();
   SMESH_ElementSearcher* searcher = const_cast<SMESH_ElementSearcher*>(myElemSearcher);
 
   set<const SMDS_MeshNode*> nodesToMove;
@@ -1057,15 +1133,17 @@ bool StdMeshers_QuadToTriaAdaptor::Compute2ndPart(SMESH_Mesh&
       PsI[k] = SMESH_TNodeXYZ( n );
       SMDS_ElemIteratorPtr vIt = n->GetInverseElementIterator( SMDSAbs_Volume );
       while ( vIt->more() )
-        checkedPyrams.insert( vIt->next() );
+      {
+        const SMDS_MeshElement* PrmJ = vIt->next();
+        if ( SMESH_Algo::GetCommonNodes( PrmI, PrmJ ).size() > 1 )
+          checkedPyrams.insert( PrmJ );
+      }
     }
 
     // check intersection with distant pyramids
     for(k=0; k<4; k++) // loop on 4 base nodes of PrmI
     {
       gp_Vec Vtmp(PsI[k],PsI[4]);
-      gp_Pnt Pshift = PsI[k].XYZ() + Vtmp.XYZ() * 0.01; // base node moved a bit to apex
-
       gp_Ax1 line( PsI[k], Vtmp );
       vector< const SMDS_MeshElement* > suspectPyrams;
       searcher->GetElementsNearLine( line, SMDSAbs_Volume, suspectPyrams);
@@ -1086,12 +1164,16 @@ bool StdMeshers_QuadToTriaAdaptor::Compute2ndPart(SMESH_Mesh&
         vector<gp_Pnt> PsJ( xyzIt, TXyzIterator() );
 
         gp_Pnt Pint;
-        bool hasInt = 
+        bool hasInt=false;
+        for(k=0; k<4 && !hasInt; k++) {
+          gp_Vec Vtmp(PsI[k],PsI[4]);
+          gp_Pnt Pshift = PsI[k].XYZ() + Vtmp.XYZ() * 0.01; // base node moved a bit to apex
+          hasInt = 
           ( HasIntersection3( Pshift, PsI[4], Pint, PsJ[0], PsJ[1], PsJ[4]) ||
             HasIntersection3( Pshift, PsI[4], Pint, PsJ[1], PsJ[2], PsJ[4]) ||
             HasIntersection3( Pshift, PsI[4], Pint, PsJ[2], PsJ[3], PsJ[4]) ||
             HasIntersection3( Pshift, PsI[4], Pint, PsJ[3], PsJ[0], PsJ[4]) );
-
+        }
         for(k=0; k<4 && !hasInt; k++) {
           gp_Vec Vtmp(PsJ[k],PsJ[4]);
           gp_Pnt Pshift = PsJ[k].XYZ() + Vtmp.XYZ() * 0.01;
@@ -1132,8 +1214,8 @@ bool StdMeshers_QuadToTriaAdaptor::Compute2ndPart(SMESH_Mesh&
             gp_Vec VI2(PCj,Pint);
             double ang1 = fabs(VN1.Angle(VI1));
             double ang2 = fabs(VN2.Angle(VI2));
-            double coef1 = 0.5 - (( ang1<PI/3 ) ? cos(ang1)*0.25 : 0 );
-            double coef2 = 0.5 - (( ang2<PI/3 ) ? cos(ang2)*0.25 : 0 ); // cos(ang2) ?
+            double coef1 = 0.5 - (( ang1 < M_PI/3. ) ? cos(ang1)*0.25 : 0 );
+            double coef2 = 0.5 - (( ang2 < M_PI/3. ) ? cos(ang2)*0.25 : 0 ); // cos(ang2) ?
 //             double coef2 = 0.5;
 //             if(ang2<PI/3)
 //               coef2 -= cos(ang1)*0.25;
@@ -1164,6 +1246,10 @@ bool StdMeshers_QuadToTriaAdaptor::Compute2ndPart(SMESH_Mesh&
       meshDS->MoveNode( *n, (*n)->X(), (*n)->Y(), (*n)->Z() );
   }
 
+  // move medium nodes of merged quadratic pyramids
+  if ( myPyramids[0]->IsQuadratic() )
+    UpdateQuadraticPyramids( nodesToMove, GetMeshDS() );
+
   // erase removed triangles from the proxy mesh
   if ( !myRemovedTrias.empty() )
   {
@@ -1193,15 +1279,3 @@ bool StdMeshers_QuadToTriaAdaptor::Compute2ndPart(SMESH_Mesh&
 
   return true;
 }
-
-//================================================================================
-/*!
- * \brief Return list of created triangles for given face
- */
-//================================================================================
-
-// const list<const SMDS_MeshFace* >* StdMeshers_QuadToTriaAdaptor::GetTriangles (const SMDS_MeshElement* aQuad)
-// {
-//   TQuad2Trias::iterator it = myResMap.find(aQuad);
-//   return ( it != myResMap.end() ?  & it->second : 0 );
-// }