]> SALOME platform Git repositories - modules/smesh.git/commitdiff
Salome HOME
Merge remote branch 'origin/V7_dev' into V8_0_0_BR
authorvsr <vsr@opencascade.com>
Wed, 3 Feb 2016 14:09:39 +0000 (17:09 +0300)
committervsr <vsr@opencascade.com>
Wed, 3 Feb 2016 14:09:39 +0000 (17:09 +0300)
src/Controls/SMESH_Controls.cxx
src/SMESH/SMESH_MeshEditor.cxx
src/SMESH_I/SMESH_MeshEditor_i.cxx
src/StdMeshers/StdMeshers_Prism_3D.cxx
src/StdMeshers/StdMeshers_Prism_3D.hxx

index 3f612a24fc823cc2b4750d44f01095709f07ba70..0623399e9dba8246b5958f7f32117d1937ec5324 100644 (file)
@@ -3212,14 +3212,13 @@ bool RangeOfIds::SetRangeStr( const TCollection_AsciiString& theStr )
   myIds.Clear();
 
   TCollection_AsciiString aStr = theStr;
-  //aStr.RemoveAll( ' ' );
-  //aStr.RemoveAll( '\t' );
   for ( int i = 1; i <= aStr.Length(); ++i )
-    if ( isspace( aStr.Value( i )))
-      aStr.SetValue( i, ',');
-
-  for ( int aPos = aStr.Search( ",," ); aPos != -1; aPos = aStr.Search( ",," ) )
-    aStr.Remove( aPos, 1 );
+  {
+    char c = aStr.Value( i );
+    if ( !isdigit( c ) && c != ',' && c != '-' )
+      aStr.SetValue( i, ' ');
+  }
+  aStr.RemoveAll( ' ' );
 
   TCollection_AsciiString tmpStr = aStr.Token( ",", 1 );
   int i = 1;
index 56371d7112930ca2ba49711fe92bed76519cbf84..0f2f1ad710c39fa4159264fa4c5290fa8a0d5230 100644 (file)
@@ -1601,13 +1601,13 @@ void SMESH_MeshEditor::QuadTo4Tri (TIDSortedElemSet & theElems)
 
     // create 4 triangles
 
-    GetMeshDS()->RemoveFreeElement( quad, subMeshDS, /*fromGroups=*/false );
-    
     helper.SetIsQuadratic  ( nodes.size() > 4 );
     helper.SetIsBiQuadratic( nodes.size() == 9 );
     if ( helper.GetIsQuadratic() )
       helper.AddTLinks( static_cast< const SMDS_MeshFace*>( quad ));
 
+    GetMeshDS()->RemoveFreeElement( quad, subMeshDS, /*fromGroups=*/false );
+
     for ( int i = 0; i < 4; ++i )
     {
       SMDS_MeshElement* tria = helper.AddFace( nodes[ i ],
@@ -2979,7 +2979,7 @@ bool SMESH_MeshEditor::QuadToTri (TIDSortedElemSet & theElems,
 
     // Quadratic quadrangle
 
-    else if ( elem->NbNodes() == 8 )
+    else if ( elem->NbNodes() >= 8 )
     {
       // get surface elem is on
       int aShapeId = FindShape( elem );
@@ -2996,52 +2996,34 @@ bool SMESH_MeshEditor::QuadToTri (TIDSortedElemSet & theElems,
         }
       }
 
-      const SMDS_MeshNode* aNodes [8];
-      const SMDS_MeshNode* inFaceNode = 0;
+      const SMDS_MeshNode* aNodes [9]; aNodes[8] = 0;
       SMDS_ElemIteratorPtr itN = elem->nodesIterator();
-      int i = 0;
-      if ( helper.GetNodeUVneedInFaceNode() )
-        while ( itN->more() && !inFaceNode ) {
-          aNodes[ i++ ] = static_cast<const SMDS_MeshNode*>( itN->next() );
-          if ( aNodes[ i-1 ]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_FACE )
-          {
-            inFaceNode = aNodes[ i-1 ];
-          }
-        }
+      for ( int i = 0; itN->more(); ++i )
+        aNodes[ i ] = static_cast<const SMDS_MeshNode*>( itN->next() );
 
-      // find middle point for (0,1,2,3)
-      // and create a node in this point;
-      gp_XYZ p( 0,0,0 );
-      if ( surface.IsNull() ) {
-        for(i=0; i<4; i++)
-          p += gp_XYZ(aNodes[i]->X(), aNodes[i]->Y(), aNodes[i]->Z() );
-        p /= 4;
-      }
-      else {
-        TopoDS_Face geomFace = TopoDS::Face( helper.GetSubShape() );
-        gp_XY uv( 0,0 );
-        for(i=0; i<4; i++)
-          uv += helper.GetNodeUV( geomFace, aNodes[i], inFaceNode );
-        uv /= 4.;
-        p = surface->Value( uv.X(), uv.Y() ).XYZ();
+      const SMDS_MeshNode* centrNode = aNodes[8];
+      if ( centrNode == 0 )
+      {
+        centrNode = helper.GetCentralNode( aNodes[0], aNodes[1], aNodes[2], aNodes[3],
+                                           aNodes[4], aNodes[5], aNodes[6], aNodes[7],
+                                           surface.IsNull() );
+        myLastCreatedNodes.Append(centrNode);
       }
-      const SMDS_MeshNode* newN = aMesh->AddNode( p.X(), p.Y(), p.Z() );
-      myLastCreatedNodes.Append(newN);
 
       // create a new element
       const SMDS_MeshElement* newElem1 = 0;
       const SMDS_MeshElement* newElem2 = 0;
       if ( the13Diag ) {
         newElem1 = aMesh->AddFace(aNodes[2], aNodes[3], aNodes[0],
-                                  aNodes[6], aNodes[7], newN );
+                                  aNodes[6], aNodes[7], centrNode );
         newElem2 = aMesh->AddFace(aNodes[2], aNodes[0], aNodes[1],
-                                  newN,      aNodes[4], aNodes[5] );
+                                  centrNode, aNodes[4], aNodes[5] );
       }
       else {
         newElem1 = aMesh->AddFace(aNodes[3], aNodes[0], aNodes[1],
-                                  aNodes[7], aNodes[4], newN );
+                                  aNodes[7], aNodes[4], centrNode );
         newElem2 = aMesh->AddFace(aNodes[3], aNodes[1], aNodes[2],
-                                  newN,      aNodes[5], aNodes[6] );
+                                  centrNode, aNodes[5], aNodes[6] );
       }
       myLastCreatedElems.Append(newElem1);
       myLastCreatedElems.Append(newElem2);
index f31e583ebe935941469bc286e6a5755fa6469fd7..bd08f9d67785800357512b06fb0b24fe07213f70 100644 (file)
@@ -5304,7 +5304,7 @@ bool SMESH_MeshEditor_i::idSourceToSet(SMESH::SMESH_IDSource_ptr  theIDSource,
 {
   if ( error ) *error = IDSource_OK;
 
-  if ( CORBA::is_nil( theIDSource ) )
+  if ( CORBA::is_nil( theIDSource ))
   {
     if ( error ) *error = IDSource_INVALID;
     return false;
index 1485332b3a83ff81669513bfe1d46dc3f9620669..645bad1e2fe5c47a251bef2637fe531d34f09998 100644 (file)
@@ -62,6 +62,7 @@
 #include <gp_Ax3.hxx>
 
 #include <limits>
+#include <numeric>
 
 using namespace std;
 
@@ -977,7 +978,7 @@ bool StdMeshers_Prism_3D::getWallFaces( Prism_3D::TPrismTopo & thePrism,
     {
       iE = 0;
       ++nbE;
-      int nbQuadPrev = nbQuadsPerWire.empty() ? 0 : nbQuadsPerWire.back();
+      int nbQuadPrev = std::accumulate( nbQuadsPerWire.begin(), nbQuadsPerWire.end(), 0 );
       nbQuadsPerWire.push_back( thePrism.myWallQuads.size() - nbQuadPrev );
     }
   }
@@ -1291,6 +1292,7 @@ bool StdMeshers_Prism_3D::compute(const Prism_3D::TPrismTopo& thePrism)
   if ( !smDS ) return toSM( error(COMPERR_BAD_INPUT_MESH, "Null submesh"));
 
   // loop on bottom mesh faces
+  vector< const TNodeColumn* > columns;
   SMDS_ElemIteratorPtr faceIt = smDS->GetElements();
   while ( faceIt->more() )
   {
@@ -1300,7 +1302,7 @@ bool StdMeshers_Prism_3D::compute(const Prism_3D::TPrismTopo& thePrism)
 
     // find node columns for each node
     int nbNodes = face->NbCornerNodes();
-    vector< const TNodeColumn* > columns( nbNodes );
+    columns.resize( nbNodes );
     for ( int i = 0; i < nbNodes; ++i )
     {
       const SMDS_MeshNode* n = face->GetNode( i );
@@ -1317,7 +1319,8 @@ bool StdMeshers_Prism_3D::compute(const Prism_3D::TPrismTopo& thePrism)
       }
     }
     // create prisms
-    AddPrisms( columns, myHelper );
+    if ( !AddPrisms( columns, myHelper ))
+      return toSM( error("Different 'vertical' discretization"));
 
   } // loop on bottom mesh faces
 
@@ -1815,17 +1818,20 @@ bool StdMeshers_Prism_3D::Evaluate(SMESH_Mesh&         theMesh,
  */
 //================================================================================
 
-void StdMeshers_Prism_3D::AddPrisms( vector<const TNodeColumn*> & columns,
+bool StdMeshers_Prism_3D::AddPrisms( vector<const TNodeColumn*> & columns,
                                      SMESH_MesherHelper*          helper)
 {
-  int nbNodes = columns.size();
-  int nbZ     = columns[0]->size();
-  if ( nbZ < 2 ) return;
+  size_t nbNodes = columns.size();
+  size_t nbZ     = columns[0]->size();
+  if ( nbZ < 2 ) return false;
+  for ( size_t i = 1; i < nbNodes; ++i )
+    if ( columns[i]->size() != nbZ )
+      return false;
 
   // find out orientation
   bool isForward = true;
   SMDS_VolumeTool vTool;
-  int z = 1;
+  size_t z = 1;
   switch ( nbNodes ) {
   case 3: {
     SMDS_VolumeOfNodes tmpPenta ( (*columns[0])[z-1], // bottom
@@ -1911,7 +1917,7 @@ void StdMeshers_Prism_3D::AddPrisms( vector<const TNodeColumn*> & columns,
     vector<const SMDS_MeshNode*> nodes( 2*nbNodes + 4*nbNodes);
     for ( z = 1; z < nbZ; ++z )
     {
-      for ( int i = 0; i < nbNodes; ++i ) {
+      for ( size_t i = 0; i < nbNodes; ++i ) {
         nodes[ i             ] = (*columns[ i ])[z+iBase1]; // bottom or top
         nodes[ 2*nbNodes-i-1 ] = (*columns[ i ])[z+iBase2]; // top or bottom
         // side
@@ -1925,6 +1931,8 @@ void StdMeshers_Prism_3D::AddPrisms( vector<const TNodeColumn*> & columns,
     }
 
   } // switch ( nbNodes )
+
+  return true;
 }
 
 //================================================================================
index 708a92ca9226a29e0cf948b76d84266fc737812c..0a24284016a0edcaacbe61ab6ff731be9cbb679e 100644 (file)
@@ -472,7 +472,7 @@ public:
     * \param nodeColumns - columns of nodes generated from nodes of a mesh face
     * \param helper - helper initialized by mesh and shape to add prisms to
    */
-  static void AddPrisms( std::vector<const TNodeColumn*> & nodeColumns,
+  static bool AddPrisms( std::vector<const TNodeColumn*> & nodeColumns,
                          SMESH_MesherHelper*               helper);
 
   static bool IsApplicable(const TopoDS_Shape & aShape, bool toCheckAll);