]> SALOME platform Git repositories - modules/smesh.git/commitdiff
Salome HOME
#17351 [CEA] Mesh with Polyhedron eap/17351_polyhedron
authoreap <eap@opencascade.com>
Thu, 24 Oct 2019 13:34:16 +0000 (16:34 +0300)
committereap <eap@opencascade.com>
Thu, 24 Oct 2019 13:34:16 +0000 (16:34 +0300)
12 files changed:
doc/salome/gui/SMESH/input/basic_meshing_algos.rst
idl/SMESH_BasicHypothesis.idl
resources/StdMeshers.xml.in
src/SMESH_SWIG/StdMeshersBuilder.py
src/StdMeshers/CMakeLists.txt
src/StdMeshers/StdMeshers_PolygonPerFace_2D.cxx
src/StdMeshers/StdMeshers_PolyhedronPerSolid_3D.cxx [new file with mode: 0644]
src/StdMeshers/StdMeshers_PolyhedronPerSolid_3D.hxx [new file with mode: 0644]
src/StdMeshers_I/CMakeLists.txt
src/StdMeshers_I/StdMeshers_PolyhedronPerSolid_3D_i.cxx [new file with mode: 0644]
src/StdMeshers_I/StdMeshers_PolyhedronPerSolid_3D_i.hxx [new file with mode: 0644]
src/StdMeshers_I/StdMeshers_i.cxx

index 692232fd13eff17bfbc9a6f6c4cae80fbda0063f..d636b15f9e100c5506d1e1686350113e5dc7b8ce 100644 (file)
@@ -53,6 +53,7 @@ An algorithm represents either an implementation of a certain meshing technique
    * :ref:`Extrusion 3D <prism_3d_algo_page>` - for meshing prismatic 3D shapes with hexahedra and prisms.
    * :ref:`Quadrangle: Medial Axis Projection <quad_from_ma_algo_page>` - for quadrangle meshing of faces with sinuous borders and rings.
    * **Polygon per Face** meshing algorithm - generates one mesh face (either a triangle, a quadrangle or a polygon) per a geometrical face using all nodes from the face boundary.
+   * **Polyhedron per Solid** meshing algorithm - generates one mesh volume (of a classical type or a polyhedron) per a geometrical solid using all faces of the solid boundary. It does not require that 2D mesh is generated on geometrical faces. It creates one mesh edge per geometrical edges and applies **Polygon per Face** to faces if they are not meshed by optional algorithms of lower dimensions.
    * :ref:`Projection algorithms <projection_algos_page>` - for meshing by projection of another mesh.
    * :ref:`Import algorithms <import_algos_page>` - for meshing by importing elements from another mesh.
    * :ref:`Radial Prism <radial_prism_algo_page>` - for meshing 3D geometrical objects with cavities with hexahedra and prisms.
index a4bb746a64ef4ad059e0c35039c5e9afcef0c6b9..1f9439f0126056287999a2c6c436d00d79b5e995 100644 (file)
@@ -1114,6 +1114,13 @@ module StdMeshers
   {
   };
 
+  /*!
+   * StdMeshers_PolyhedronPerSolid_3D: interface of "Polyhedron Per Solid" 3D algorithm
+   */
+  interface StdMeshers_PolyhedronPerSolid_3D : SMESH::SMESH_3D_Algo
+  {
+  };
+
   /*!
    * StdMeshers_Hexa_3D: interface of "Hexahedron (i,j,k)" algorithm
    */
index 14ed93a0c9b0abefee7f9e86d5df842ed55b678a..fc46a1350bcb7973d406568c478c33f727d95a19 100644 (file)
       </python-wrap>
     </algorithm>
 
+    <algorithm type     ="PolyhedronPerSolid_3D"
+               label-id ="Polyhedron per Solid"
+               icon-id  ="mesh_algo_polygon.png"
+               opt-hypos="ViscousLayers"
+               input    ="POLYGON,QUAD,TRIA,EDGE"
+               dim      ="3">
+      <python-wrap>
+        <algo>PolyhedronPerSolid_3D=Polyhedron()</algo>
+        <hypo>ViscousLayers=ViscousLayers(SetTotalThickness(),SetNumberLayers(),SetStretchFactor(),SetIgnoreEdges())</hypo>
+      </python-wrap>
+    </algorithm>
+
     <algorithm type     ="Hexa_3D"
                label-id ="Hexahedron (i,j,k)"
                icon-id  ="mesh_algo_hexa.png"
index 858820238ad520dc15547e89a0b6bbd581b030d3..735e8382682e0e9d8fb2b7325aad8fcea8a6e861 100644 (file)
@@ -75,6 +75,11 @@ POLYGON     = "PolygonPerFace_2D"
 Algorithm type: Polygon Per Face 2D algorithm, see :class:`~StdMeshersBuilder.StdMeshersBuilder_PolygonPerFace`
 """
 
+POLYHEDRON = "PolyhedronPerSolid_3D"
+"""
+Algorithm type: Polyhedron Per Solid 3D algorithm, see :class:`~StdMeshersBuilder.StdMeshersBuilder_PolyhedronPerSolid`
+"""
+
 # import items of enums
 for e in StdMeshers.QuadType._items: exec('%s = StdMeshers.%s'%(e,e))
 for e in StdMeshers.VLExtrusionMethod._items: exec('%s = StdMeshers.%s'%(e,e))
@@ -1670,6 +1675,44 @@ class StdMeshersBuilder_PolygonPerFace(Mesh_Algorithm):
 
     pass
 
+class StdMeshersBuilder_PolyhedronPerSolid(Mesh_Algorithm):
+    """ Defines a Polyhedron Per Solid 3D algorithm.
+        It is created by calling smeshBuilder.Mesh.Polyhedron(geom=0)
+    """
+
+    meshMethod = "Polyhedron"
+    """
+    name of the dynamic method in smeshBuilder.Mesh class
+    """
+    algoType   = POLYHEDRON
+    """
+    type of algorithm used with helper function in smeshBuilder.Mesh class
+    """
+    isDefault  = True
+    """
+    flag pointing whether this algorithm should be used by default in dynamic method
+        of smeshBuilder.Mesh class
+    """
+    docHelper  = "Create polyhedron 3D algorithm for solids"
+    """
+    doc string of the method
+    """
+
+    def __init__(self, mesh, geom=0):
+        """
+        Private constructor.
+
+        Parameters:
+            mesh: parent mesh object algorithm is assigned to
+            geom: geometry (shape/sub-shape) algorithm is assigned to;
+                if it is :code:`0` (default), the algorithm is assigned to the main shape
+        """
+        Mesh_Algorithm.__init__(self)
+        self.Create(mesh, geom, self.algoType)
+        pass
+
+    pass
+
 class StdMeshersBuilder_UseExistingElements_1D(Mesh_Algorithm):
     """ Defines a Use Existing Elements 1D algorithm.
 
index dfad28a865eb5ac0332f018592be1ca607c280c3..df1d952432f7598f66c2873971d91fa0162107b4 100644 (file)
@@ -130,6 +130,7 @@ SET(StdMeshers_HEADERS
   StdMeshers_Cartesian_3D.hxx
   StdMeshers_QuadFromMedialAxis_1D2D.hxx
   StdMeshers_PolygonPerFace_2D.hxx
+  StdMeshers_PolyhedronPerSolid_3D.hxx
 )
 
 IF(SALOME_SMESH_ENABLE_MEFISTO)
@@ -195,6 +196,7 @@ SET(StdMeshers_SOURCES
   StdMeshers_Adaptive1D.cxx
   StdMeshers_QuadFromMedialAxis_1D2D.cxx
   StdMeshers_PolygonPerFace_2D.cxx
+  StdMeshers_PolyhedronPerSolid_3D.cxx
 )
 
 IF(SALOME_SMESH_ENABLE_MEFISTO)
index c8d203d12eb54c9444e8d2da5d287b3af0a4ee29..76afc6a0c6b8c3b3678a26e0a29e9b0ac3a576b9 100644 (file)
 
 #include <TopExp_Explorer.hxx>
 #include <TopoDS_Face.hxx>
+#include <TopoDS.hxx>
 
 #include <vector>
-#include <TopoDS.hxx>
 
 using namespace std;
 
 //=======================================================================
 //function : StdMeshers_PolygonPerFace_2D
-//purpose  : 
+//purpose  :
 //=======================================================================
 
 StdMeshers_PolygonPerFace_2D::StdMeshers_PolygonPerFace_2D(int        hypId,
diff --git a/src/StdMeshers/StdMeshers_PolyhedronPerSolid_3D.cxx b/src/StdMeshers/StdMeshers_PolyhedronPerSolid_3D.cxx
new file mode 100644 (file)
index 0000000..6695eea
--- /dev/null
@@ -0,0 +1,618 @@
+// Copyright (C) 2007-2019  CEA/DEN, EDF R&D, OPEN CASCADE
+//
+// Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
+// CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
+//
+// 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, or (at your option) any later version.
+//
+// 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
+//
+// See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
+//
+
+// File      : StdMeshers_PolyhedronPerSolid_3D.cxx
+// Module    : SMESH
+// Created   : Fri Oct 20 11:37:07 2006
+// Author    : Edward AGAPOV (eap)
+//
+#include "StdMeshers_PolyhedronPerSolid_3D.hxx"
+
+#include "SMESHDS_Mesh.hxx"
+#include "SMESH_ControlsDef.hxx"
+#include "SMESH_Gen.hxx"
+#include "SMESH_Mesh.hxx"
+#include "SMESH_MeshAlgos.hxx"
+#include "SMESH_MeshEditor.hxx"
+#include "SMESH_MesherHelper.hxx"
+#include "SMESH_ProxyMesh.hxx"
+#include "SMESH_subMesh.hxx"
+#include "StdMeshers_PolygonPerFace_2D.hxx"
+#include "StdMeshers_Regular_1D.hxx"
+#include "StdMeshers_ViscousLayers.hxx"
+
+#include <TopExp_Explorer.hxx>
+
+#include <vector>
+#include <set>
+
+namespace
+{
+  struct _EdgeMesher : public StdMeshers_Regular_1D
+  {
+    _EdgeMesher( int hypId, SMESH_Gen* gen )
+      : StdMeshers_Regular_1D( hypId, gen )
+    {
+      _hypType = NB_SEGMENTS;
+      _ivalue[ NB_SEGMENTS_IND ] = 1;
+    }
+  };
+
+  //=======================================================================
+  //function : addHexa
+  //purpose  :
+  //=======================================================================
+
+  const SMDS_MeshElement* addHexa( std::vector< const SMDS_MeshElement* >& faces,
+                                   const std::vector< int > &              quantities,
+                                   SMESH_MesherHelper &                    helper )
+  {
+    const SMDS_MeshElement* newHexa = 0;
+
+    // check nb of nodes in faces
+    for ( size_t i = 0; i < quantities.size(); ++i )
+      if ( quantities[ i ] != 4 )
+        return newHexa;
+
+    // look for a top face
+    const SMDS_MeshElement* topFace = 0;
+    const SMDS_MeshElement* botFace = faces[0];
+    std::vector< const SMDS_MeshNode* > nodes( 16 ); // last 8 is a working buffer
+    nodes.assign( botFace->begin_nodes(), botFace->end_nodes() );
+    for ( size_t iF = 1; iF < faces.size() &&  !topFace; ++iF )
+    {
+      bool hasCommonNode = false;
+      for ( int iN = 0; iN < quantities[ 0 ] &&  !hasCommonNode; ++iN )
+        hasCommonNode = ( faces[ iF ]->GetNodeIndex( nodes[ iN ]) >= 0 );
+
+      if ( !hasCommonNode )
+        topFace = faces[ iF ];
+    }
+
+    nodes.resize( 8 ); // set top nodes after hexa nodes - [8-11]
+    nodes.insert( nodes.end(), topFace->begin_nodes(), topFace->end_nodes() );
+    nodes.resize( 12 );
+    nodes.insert( nodes.end(), nodes.begin() + 8, nodes.begin() + 12 );
+
+    // find corresponding nodes of top and bottom by finding a side face including 2 node of each
+    SMESHDS_Mesh* mesh = helper.GetMeshDS();
+    const SMDS_MeshElement* sideFace = 0;
+    size_t i;
+    for ( i = 8; i < nodes.size()-1 &&  !sideFace; ++i )
+    {
+      sideFace = mesh->FindFace( nodes[0], nodes[1], nodes[ i ], nodes[ i + 1 ]);
+    }
+    if ( !sideFace )
+      return newHexa;
+
+    --i; // restore after ++i in the loop
+    bool botOriRight = SMESH_MeshAlgos::IsRightOrder( sideFace, nodes[ 0 ], nodes[ 1 ] );
+    bool topOriRight = SMESH_MeshAlgos::IsRightOrder( sideFace, nodes[ i ], nodes[ i + 1 ] );
+    if ( botOriRight == topOriRight )
+    {
+      nodes[ 4 ] = nodes[ i + 1 ];
+      nodes[ 5 ] = nodes[ i + 0 ];
+      nodes[ 6 ] = nodes[ i + 3 ];
+      nodes[ 7 ] = nodes[ i + 2 ];
+    }
+    else
+    {
+      nodes[ 4 ] = nodes[ i + 0 ];
+      nodes[ 5 ] = nodes[ i + 1 ];
+      nodes[ 6 ] = nodes[ i + 2 ];
+      nodes[ 7 ] = nodes[ i + 3 ];
+    }
+
+    newHexa = helper.AddVolume( nodes[ 0 ], nodes[ 1 ], nodes[ 2 ], nodes[ 3 ],
+                                nodes[ 4 ], nodes[ 5 ], nodes[ 6 ], nodes[ 7 ]);
+
+    return newHexa;
+  }
+
+  //=======================================================================
+  //function : addTetra
+  //purpose  :
+  //=======================================================================
+
+  const SMDS_MeshElement* addTetra( std::vector< const SMDS_MeshElement* >& faces,
+                                    const std::vector< int > &              quantities,
+                                    SMESH_MesherHelper &                    helper )
+  {
+    const SMDS_MeshElement* newTetra = 0;
+
+    // check nb of nodes in faces
+    for ( size_t i = 0; i < quantities.size(); ++i )
+      if ( quantities[ i ] != 3 )
+        return newTetra;
+
+    const SMDS_MeshElement* botFace = faces[0];
+
+    std::vector< const SMDS_MeshNode* > nodes( 6 );
+    nodes.assign( botFace->begin_nodes(), botFace->end_nodes() );
+    nodes.resize( 3 );
+
+    const SMDS_MeshNode* topNode = 0;
+    for ( size_t i = 0; i < 3 &&  !topNode; ++i )
+    {
+      topNode = faces[ 1 ]->GetNode( i );
+      if ( botFace->GetNodeIndex( topNode ) >= 0 )
+        topNode = 0;
+    }
+    if ( !topNode )
+      return newTetra;
+
+    nodes.push_back( topNode );
+
+    newTetra = helper.AddVolume( nodes[ 0 ], nodes[ 1 ], nodes[ 2 ], nodes[ 3 ]);
+
+    return newTetra;
+  }
+
+  //=======================================================================
+  //function : addPenta
+  //purpose  :
+  //=======================================================================
+
+  const SMDS_MeshElement* addPenta( std::vector< const SMDS_MeshElement* >& faces,
+                                    const std::vector< int > &              quantities,
+                                    SMESH_MesherHelper &                    helper )
+  {
+    const SMDS_MeshElement* newPenta = 0;
+
+    // check nb of nodes in faces and find triangle faces
+    int trias[2] = { -1, -1 };
+    for ( size_t i = 0; i < quantities.size(); ++i )
+      if ( quantities[ i ] != 4 )
+      {
+        if ( quantities[ i ] != 3 )
+          return newPenta;
+        int iTria = ( trias[0] != -1 );
+        if ( trias[ iTria ] != -1 )
+          return newPenta;
+        trias[ iTria ] = i;
+      }
+    if ( trias[1] == -1 )
+      return newPenta;
+
+    int iSide = trias[0] + 1;
+    if ( iSide == trias[1] )
+      ++iSide;
+
+    const SMDS_MeshElement* botFace  = faces[ trias[0]];
+    const SMDS_MeshElement* topFace  = faces[ trias[1]];
+    const SMDS_MeshElement* sideFace = faces[ iSide ];
+    const SMDS_MeshNode* nodes[ 6 ] = { 0,0,0,0,0,0 };
+    for ( int i = 0 ; i < 3; ++i )
+    {
+      const SMDS_MeshNode* botNode = botFace->GetNode( i );
+      if ( sideFace->GetNodeIndex( botNode ) < 0 )
+        nodes[2] = botNode;
+      else
+        nodes[ bool( nodes[0] )] = botNode;
+
+      const SMDS_MeshNode* topNode = topFace->GetNode( i );
+      if ( sideFace->GetNodeIndex( topNode ) < 0 )
+        nodes[5] = topNode;
+      else
+        nodes[ 3 + bool( nodes[3]) ] = topNode;
+    }
+    bool botOriRight = SMESH_MeshAlgos::IsRightOrder( sideFace, nodes[ 0 ], nodes[ 1 ]);
+    bool topOriRight = SMESH_MeshAlgos::IsRightOrder( sideFace, nodes[ 3 ], nodes[ 4 ]);
+    if ( botOriRight == topOriRight )
+      std::swap( nodes[ 3 ], nodes[ 4 ]);
+
+    newPenta = helper.AddVolume( nodes[ 0 ], nodes[ 1 ], nodes[ 2 ],
+                                 nodes[ 3 ], nodes[ 4 ], nodes[ 5 ]);
+
+    return newPenta;
+  }
+
+  //=======================================================================
+  //function : addPyra
+  //purpose  :
+  //=======================================================================
+
+  const SMDS_MeshElement* addPyra( std::vector< const SMDS_MeshElement* >& faces,
+                                   const std::vector< int > &              quantities,
+                                   SMESH_MesherHelper &                    helper )
+  {
+    const SMDS_MeshElement* newPyra = 0;
+
+    // check nb of nodes in faces
+    int iBot = -1;
+    for ( size_t i = 0; i < quantities.size(); ++i )
+      if ( quantities[ i ] != 3 )
+      {
+        if ( quantities[ i ] != 4 || iBot != -1 )
+          return newPyra;
+        iBot = i;
+      }
+
+    const SMDS_MeshElement* botFace = faces[ iBot ];
+
+    std::vector< const SMDS_MeshNode* > nodes( 8 );
+    nodes.assign( botFace->begin_nodes(), botFace->end_nodes() );
+    nodes.resize( 4 );
+
+    const SMDS_MeshNode* topNode = 0;
+    for ( size_t i = 0; i < 4 &&  !topNode; ++i )
+    {
+      topNode = faces[ 1 ]->GetNode( i );
+      if ( botFace->GetNodeIndex( topNode ) >= 0 )
+        topNode = 0;
+    }
+    if ( !topNode )
+      return newPyra;
+
+    nodes.push_back( topNode );
+
+    newPyra = helper.AddVolume( nodes[ 0 ], nodes[ 1 ], nodes[ 2 ], nodes[ 3 ], nodes[4] );
+
+    return newPyra;
+  }
+
+  //=======================================================================
+  //function : addHPrism
+  //purpose  : add hexagonal prism
+  //=======================================================================
+
+  const SMDS_MeshElement* addHPrism( std::vector< const SMDS_MeshElement* >& faces,
+                                     const std::vector< int > &              quantities,
+                                     SMESH_MesherHelper &                    helper )
+  {
+    const SMDS_MeshElement* newHexPrism = 0;
+
+    // check nb of nodes in faces and find hexagons
+    int hexa[2] = { -1, -1 };
+    for ( size_t i = 0; i < quantities.size(); ++i )
+      if ( quantities[ i ] != 4 )
+      {
+        if ( quantities[ i ] != 6 )
+          return newHexPrism;
+        int iHex = ( hexa[0] != -1 );
+        if ( hexa[ iHex ] != -1 )
+          return newHexPrism;
+        hexa[ iHex ] = i;
+      }
+    if ( hexa[1] == -1 )
+      return newHexPrism;
+
+    int iSide = hexa[0] + 1;
+    if ( iSide == hexa[1] )
+      ++iSide;
+
+    const SMDS_MeshElement* botFace = faces[ hexa[ 0 ]];
+    const SMDS_MeshElement* topFace = faces[ hexa[ 1 ]];
+    std::vector< const SMDS_MeshNode* > nodes( 24 ); // last 12 is a working buffer
+
+    nodes.assign( botFace->begin_nodes(), botFace->end_nodes() );
+    nodes.resize( 12 ); // set top nodes after hexa nodes - [12-17]
+    nodes.insert( nodes.end(), topFace->begin_nodes(), topFace->end_nodes() );
+    nodes.resize( 18 );
+    nodes.insert( nodes.end(), nodes.begin() + 12, nodes.begin() + 18 );
+
+    // find corresponding nodes of top and bottom by finding a side face including 2 node of each
+    SMESHDS_Mesh* mesh = helper.GetMeshDS();
+    const SMDS_MeshElement* sideFace = 0;
+    size_t i;
+    for ( i = 12; i < nodes.size()-1 &&  !sideFace; ++i )
+    {
+      sideFace = mesh->FindFace( nodes[0], nodes[1], nodes[ i ], nodes[ i + 1 ]);
+    }
+    if ( !sideFace )
+      return newHexPrism;
+
+    --i; // restore after ++i in the loop
+    bool botOriRight = SMESH_MeshAlgos::IsRightOrder( sideFace, nodes[ 0 ], nodes[ 1 ] );
+    bool topOriRight = SMESH_MeshAlgos::IsRightOrder( sideFace, nodes[ i ], nodes[ i + 1 ] );
+    if ( botOriRight == topOriRight )
+    {
+      nodes[ 6  ] = nodes[ i + 1 ];
+      nodes[ 7  ] = nodes[ i + 0 ];
+      nodes[ 8  ] = nodes[ i + 5 ];
+      nodes[ 9  ] = nodes[ i + 4 ];
+      nodes[ 10 ] = nodes[ i + 3 ];
+      nodes[ 11 ] = nodes[ i + 2 ];
+    }
+    else
+    {
+      nodes[ 6  ] = nodes[ i + 0 ];
+      nodes[ 7  ] = nodes[ i + 1 ];
+      nodes[ 8  ] = nodes[ i + 2 ];
+      nodes[ 9  ] = nodes[ i + 3 ];
+      nodes[ 10 ] = nodes[ i + 4 ];
+      nodes[ 11 ] = nodes[ i + 5 ];
+    }
+
+    newHexPrism = helper.AddVolume( nodes[ 0 ], nodes[ 1 ], nodes[ 2 ],
+                                    nodes[ 3 ], nodes[ 4 ], nodes[ 5 ],
+                                    nodes[ 6 ], nodes[ 7 ], nodes[ 8 ],
+                                    nodes[ 9 ], nodes[10 ], nodes[11 ]);
+
+    return newHexPrism;
+  }
+
+  //=======================================================================
+  //function : addPoly
+  //purpose  :
+  //=======================================================================
+
+  const SMDS_MeshElement* addPoly( std::vector< const SMDS_MeshElement* >& faces,
+                                   const std::vector< int > &              quantities,
+                                   SMESH_MesherHelper &                    helper )
+  {
+    const SMDS_MeshElement* newPoly = 0;
+
+    std::vector< const SMDS_MeshNode* > nodes;
+    for ( size_t iF = 0; iF < faces.size(); ++iF )
+      nodes.insert( nodes.end(), faces[iF]->begin_nodes(), faces[iF]->end_nodes() );
+
+    newPoly = helper.AddPolyhedralVolume( nodes, quantities );
+
+    return newPoly;
+  }
+
+} // namespace
+
+//=======================================================================
+//function : StdMeshers_PolyhedronPerSolid_3D
+//purpose  :
+//=======================================================================
+
+StdMeshers_PolyhedronPerSolid_3D::StdMeshers_PolyhedronPerSolid_3D(int        hypId,
+                                                                   SMESH_Gen* gen)
+  :SMESH_3D_Algo(hypId, gen),
+   myEdgeMesher( new _EdgeMesher( gen->GetANewId(), gen )),
+   myFaceMesher( new StdMeshers_PolygonPerFace_2D( gen->GetANewId(), gen ))
+{
+  _name = "PolyhedronPerSolid_3D";
+  _requireDiscreteBoundary = false;
+  _supportSubmeshes = true;
+  _compatibleHypothesis.push_back("ViscousLayers");
+  _neededLowerHyps[0] = _neededLowerHyps[1] = _neededLowerHyps[2] = true;
+}
+
+//=======================================================================
+//function : ~StdMeshers_PolyhedronPerSolid_3D
+//purpose  :
+//=======================================================================
+
+StdMeshers_PolyhedronPerSolid_3D::~StdMeshers_PolyhedronPerSolid_3D()
+{
+  delete myEdgeMesher;
+  delete myFaceMesher;
+}
+
+//=======================================================================
+//function : CheckHypothesis
+//purpose  :
+//=======================================================================
+
+bool StdMeshers_PolyhedronPerSolid_3D::CheckHypothesis(SMESH_Mesh&         theMesh,
+                                                       const TopoDS_Shape& theShape,
+                                                       Hypothesis_Status&  theStatus)
+{
+  myViscousLayersHyp = NULL;
+
+  const std::list<const SMESHDS_Hypothesis*>& hyps =
+    GetUsedHypothesis( theMesh, theShape, /*ignoreAuxiliary=*/false);
+  std::list <const SMESHDS_Hypothesis* >::const_iterator h = hyps.begin();
+  if ( h == hyps.end())
+  {
+    theStatus = SMESH_Hypothesis::HYP_OK;
+    return true;
+  }
+
+  // only StdMeshers_ViscousLayers can be used
+  theStatus = HYP_OK;
+  for ( ; h != hyps.end(); ++h )
+  {
+    if ( !(myViscousLayersHyp = dynamic_cast< const StdMeshers_ViscousLayers*> ( *h )))
+      break;
+  }
+  if ( !myViscousLayersHyp )
+    theStatus = HYP_INCOMPATIBLE;
+  else
+    error( myViscousLayersHyp->CheckHypothesis( theMesh, theShape, theStatus ));
+
+  return theStatus == HYP_OK;
+}
+
+//=======================================================================
+//function : Compute
+//purpose  :
+//=======================================================================
+
+bool StdMeshers_PolyhedronPerSolid_3D::Compute(SMESH_Mesh&         theMesh,
+                                               const TopoDS_Shape& theShape)
+{
+  const SMDS_MeshElement* newVolume = 0;
+
+  SMESH_subMesh* sm = theMesh.GetSubMesh( theShape );
+  SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator( /*includeSelf=*/true,
+                                                            /*complexFirst=*/false);
+  while ( smIt->more() )
+  {
+    sm = smIt->next();
+    if ( !sm->IsEmpty() )
+      continue;
+
+    const TopoDS_Shape & shape = sm->GetSubShape();
+    switch ( shape.ShapeType() )
+    {
+    case TopAbs_VERTEX:
+      sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
+      break;
+
+    case TopAbs_EDGE:
+      sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
+      if ( sm->IsEmpty() )
+        myEdgeMesher->Compute( theMesh, shape );
+      break;
+
+    case TopAbs_FACE:
+      sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
+      if ( sm->IsEmpty() && !myFaceMesher->Compute( theMesh, shape ))
+      {
+        sm->GetComputeError() = myFaceMesher->GetComputeError();
+        sm->GetComputeError()->myAlgo = myFaceMesher;
+        return false;
+      }
+      break;
+
+    case TopAbs_SOLID:
+    {
+      SMESH_MesherHelper helper( theMesh );
+      helper.SetElementsOnShape( true );
+      _quadraticMesh = helper.IsQuadraticSubMesh( shape );
+
+      SMESH_ProxyMesh::Ptr proxymesh( new SMESH_ProxyMesh( theMesh ));
+      if ( myViscousLayersHyp )
+      {
+        proxymesh = myViscousLayersHyp->Compute( theMesh, theShape );
+        if ( !proxymesh )
+          return false;
+      }
+
+      std::vector< const SMDS_MeshElement* > faces;
+      faces.reserve( 20 );
+
+      for ( TopExp_Explorer faceEx( shape, TopAbs_FACE ); faceEx.More(); faceEx.Next() )
+      {
+        const SMESHDS_SubMesh* smDS = proxymesh->GetSubMesh( faceEx.Current() );
+        for ( SMDS_ElemIteratorPtr faceIt = smDS->GetElements(); faceIt->more(); )
+          faces.push_back( faceIt->next() );
+      }
+
+      bool useMediumNodes = false;
+      if ( !_quadraticMesh && theMesh.GetMeshDS()->GetMeshInfo().NbFaces( ORDER_QUADRATIC ))
+        for ( size_t i = 0; i < faces.size() &&  !useMediumNodes ; ++i )
+          useMediumNodes = faces[ i ]->IsQuadratic();
+
+      std::vector< int > quantities( faces.size() );
+      std::set< const SMDS_MeshNode* > nodes;
+      for ( size_t i = 0; i < faces.size(); ++i )
+      {
+        quantities[ i ] = useMediumNodes ? faces[ i ]->NbNodes() : faces[ i ]->NbCornerNodes();
+        for ( int iN = 0; iN < quantities[ i ]; ++iN )
+          nodes.insert( faces[ i ]->GetNode( iN ));
+      }
+
+      const size_t nbNodes = nodes.size(), nbFaces = faces.size();
+      if (      nbNodes == 8  && nbFaces == 6 ) newVolume = addHexa  ( faces, quantities, helper );
+      else if ( nbNodes == 4  && nbFaces == 4 ) newVolume = addTetra ( faces, quantities, helper );
+      else if ( nbNodes == 6  && nbFaces == 5 ) newVolume = addPenta ( faces, quantities, helper );
+      else if ( nbNodes == 5  && nbFaces == 5 ) newVolume = addPyra  ( faces, quantities, helper );
+      else if ( nbNodes == 12 && nbFaces == 8 ) newVolume = addHPrism( faces, quantities, helper );
+      if ( !newVolume )
+        newVolume = addPoly ( faces, quantities, helper );
+
+      if ( newVolume )
+      {
+        SMESH::Controls::BadOrientedVolume checker;
+        checker.SetMesh( theMesh.GetMeshDS() );
+        if ( checker.IsSatisfy( newVolume->GetID() ))
+        {
+          SMESH_MeshEditor editor( &theMesh );
+          editor.Reorient( newVolume );
+        }
+      }
+    }
+    default:;
+
+    } // switch ( shape.ShapeType() )
+  } // loop on sub-meshes
+
+  return newVolume;
+}
+
+//=======================================================================
+//function : Evaluate
+//purpose  :
+//=======================================================================
+
+bool StdMeshers_PolyhedronPerSolid_3D::Evaluate(SMESH_Mesh&         theMesh,
+                                                const TopoDS_Shape& theShape,
+                                                MapShapeNbElems&    theResMap)
+{
+  _quadraticMesh = false;
+
+  SMESH_subMesh* sm = theMesh.GetSubMesh( theShape );
+  SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator( /*includeSelf=*/true,
+                                                            /*complexFirst=*/false);
+  while ( smIt->more() )
+  {
+    sm = smIt->next();
+
+    MapShapeNbElems::iterator sm2vec = theResMap.find( sm );
+    if ( sm2vec != theResMap.end() && !sm2vec->second.empty() )
+      continue;
+
+    const TopoDS_Shape & shape = sm->GetSubShape();
+    switch ( shape.ShapeType() )
+    {
+    case TopAbs_EDGE:
+      myEdgeMesher->Evaluate( theMesh, shape, theResMap );
+      break;
+
+    case TopAbs_FACE:
+    {
+      myFaceMesher->Evaluate( theMesh, shape, theResMap );
+      std::vector<int> & quantities = theResMap[ sm ];
+      _quadraticMesh = ( !quantities.empty() &&
+                         ( quantities[ SMDSEntity_Quad_Triangle   ] +
+                           quantities[ SMDSEntity_Quad_Quadrangle ] +
+                           quantities[ SMDSEntity_Quad_Polygon    ]));
+      break;
+    }
+
+    case TopAbs_SOLID:
+    {
+      std::vector<int> & quantities = theResMap[ sm ];
+      quantities.resize( SMDSEntity_Last, 0 );
+
+      SMESH_MesherHelper helper( theMesh );
+      const int nbNodes = helper.Count( shape, TopAbs_VERTEX, /*ignoreSame=*/true );
+      const int nbFaces = helper.Count( shape, TopAbs_FACE,   /*ignoreSame=*/false );
+
+      if (      nbNodes == 8 && nbFaces == 6 )
+        quantities[ _quadraticMesh ? SMDSEntity_Quad_Hexa : SMDSEntity_Hexa ] = 1;
+      else if ( nbNodes == 4 && nbFaces == 4 )
+        quantities[ _quadraticMesh ? SMDSEntity_Quad_Tetra : SMDSEntity_Tetra ] = 1;
+      else if ( nbNodes == 6 && nbFaces == 5 )
+        quantities[ _quadraticMesh ? SMDSEntity_Quad_Penta : SMDSEntity_Penta ] = 1;
+      else if ( nbNodes == 5 && nbFaces == 5 )
+        quantities[ _quadraticMesh ? SMDSEntity_Quad_Pyramid : SMDSEntity_Pyramid ] = 1;
+      else if ( nbNodes == 12 && nbFaces == 8 )
+        quantities[ /*_quadraticMesh ? SMDSEntity_Quad_Pyramid :*/ SMDSEntity_Hexagonal_Prism ] = 1;
+      else
+        quantities[ /*_quadraticMesh ? SMDSEntity_Quad_Polyhedra : */SMDSEntity_Polyhedra ] = 1;
+
+      return true;
+    }
+    default:;
+
+    } // switch ( shape.ShapeType() )
+  } // loop on sub-meshes
+
+  return false;
+}
diff --git a/src/StdMeshers/StdMeshers_PolyhedronPerSolid_3D.hxx b/src/StdMeshers/StdMeshers_PolyhedronPerSolid_3D.hxx
new file mode 100644 (file)
index 0000000..f6296aa
--- /dev/null
@@ -0,0 +1,58 @@
+// Copyright (C) 2007-2019  CEA/DEN, EDF R&D, OPEN CASCADE
+//
+// Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
+// CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
+//
+// 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, or (at your option) any later version.
+//
+// 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
+//
+// See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
+//
+
+//  File   : StdMeshers_PolyhedronPerSolid_3D.hxx
+//  Module : SMESH
+//
+#ifndef _SMESH_PolyhedronPerSolid_3D_HXX_
+#define _SMESH_PolyhedronPerSolid_3D_HXX_
+
+#include "SMESH_StdMeshers.hxx"
+#include "SMESH_Algo.hxx"
+
+class StdMeshers_Regular_1D;
+class StdMeshers_PolygonPerFace_2D;
+class StdMeshers_ViscousLayers;
+
+class STDMESHERS_EXPORT StdMeshers_PolyhedronPerSolid_3D: public SMESH_3D_Algo
+{
+ public:
+  StdMeshers_PolyhedronPerSolid_3D(int hypId, SMESH_Gen* gen);
+  ~StdMeshers_PolyhedronPerSolid_3D();
+
+  virtual bool CheckHypothesis(SMESH_Mesh&                          aMesh,
+                               const TopoDS_Shape&                  aShape,
+                               SMESH_Hypothesis::Hypothesis_Status& aStatus);
+
+  virtual bool Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape);
+
+  virtual bool Evaluate(SMESH_Mesh & aMesh, const TopoDS_Shape & aShape,
+                        MapShapeNbElems& aResMap);
+
+ private:
+
+  StdMeshers_Regular_1D*          myEdgeMesher;
+  StdMeshers_PolygonPerFace_2D*   myFaceMesher;
+  const StdMeshers_ViscousLayers* myViscousLayersHyp;
+};
+
+#endif
index ceed790b4cb1e5195e0af9cfccddf506c1e22c7b..4b1c390bad06c391ab4d9d11ff38b9611157fb72 100644 (file)
@@ -119,6 +119,7 @@ SET(StdMeshersEngine_HEADERS
   StdMeshers_CartesianParameters3D_i.hxx
   StdMeshers_Cartesian_3D_i.hxx
   StdMeshers_PolygonPerFace_2D_i.hxx
+  StdMeshers_PolyhedronPerSolid_3D_i.hxx
 )
 IF(SALOME_SMESH_ENABLE_MEFISTO)
   SET(StdMeshersEngine_HEADERS ${StdMeshersEngine_HEADERS} StdMeshers_MEFISTO_2D_i.hxx)
@@ -174,6 +175,7 @@ SET(StdMeshersEngine_SOURCES
   StdMeshers_Cartesian_3D_i.cxx
   StdMeshers_Adaptive1D_i.cxx 
   StdMeshers_PolygonPerFace_2D_i.cxx
+  StdMeshers_PolyhedronPerSolid_3D_i.cxx
 )
 
 IF(SALOME_SMESH_ENABLE_MEFISTO)
diff --git a/src/StdMeshers_I/StdMeshers_PolyhedronPerSolid_3D_i.cxx b/src/StdMeshers_I/StdMeshers_PolyhedronPerSolid_3D_i.cxx
new file mode 100644 (file)
index 0000000..94d5648
--- /dev/null
@@ -0,0 +1,57 @@
+// Copyright (C) 2007-2019  CEA/DEN, EDF R&D, OPEN CASCADE
+//
+// Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
+// CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
+//
+// 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, or (at your option) any later version.
+//
+// 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
+//
+// See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
+//
+
+//  File   : StdMeshers_PolyhedronPerSolid_3D_i.cxx
+//  Module : SMESH
+//
+
+#include "StdMeshers_PolyhedronPerSolid_3D_i.hxx"
+
+#include "SMESH_Gen.hxx"
+#include "StdMeshers_PolyhedronPerSolid_3D.hxx"
+
+//=============================================================================
+/*!
+ *  Constructor
+ */
+//=============================================================================
+
+StdMeshers_PolyhedronPerSolid_3D_i::StdMeshers_PolyhedronPerSolid_3D_i( PortableServer::POA_ptr thePOA,
+                                                                        ::SMESH_Gen*            theGenImpl )
+  : SALOME::GenericObj_i( thePOA ),
+    SMESH_Hypothesis_i( thePOA ),
+    SMESH_Algo_i( thePOA ),
+    SMESH_3D_Algo_i( thePOA )
+{
+  myBaseImpl = new ::StdMeshers_PolyhedronPerSolid_3D( theGenImpl->GetANewId(),
+                                                       theGenImpl );
+}
+
+//=============================================================================
+/*!
+ *  Destructor
+ */
+//=============================================================================
+
+StdMeshers_PolyhedronPerSolid_3D_i::~StdMeshers_PolyhedronPerSolid_3D_i()
+{
+}
diff --git a/src/StdMeshers_I/StdMeshers_PolyhedronPerSolid_3D_i.hxx b/src/StdMeshers_I/StdMeshers_PolyhedronPerSolid_3D_i.hxx
new file mode 100644 (file)
index 0000000..26bdffc
--- /dev/null
@@ -0,0 +1,53 @@
+// Copyright (C) 2007-2019  CEA/DEN, EDF R&D, OPEN CASCADE
+//
+// Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
+// CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
+//
+// 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, or (at your option) any later version.
+//
+// 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
+//
+// See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
+//
+
+//  File   : StdMeshers_PolyhedronPerSolid_3D_i.hxx
+//  Module : SMESH
+//
+#ifndef _SMESH_PolyhedronPerSolid_3D_I_HXX_
+#define _SMESH_PolyhedronPerSolid_3D_I_HXX_
+
+#include "SMESH_StdMeshers_I.hxx"
+
+#include <SALOMEconfig.h>
+#include CORBA_SERVER_HEADER(SMESH_BasicHypothesis)
+
+#include "SMESH_3D_Algo_i.hxx"
+
+class SMESH_Gen;
+
+// ======================================================
+// Polyhedron Per Solid 3D algorithm
+// ======================================================
+class STDMESHERS_I_EXPORT StdMeshers_PolyhedronPerSolid_3D_i:
+  public virtual POA_StdMeshers::StdMeshers_PolyhedronPerSolid_3D,
+  public virtual SMESH_3D_Algo_i
+{
+ public:
+  // Constructor
+  StdMeshers_PolyhedronPerSolid_3D_i( PortableServer::POA_ptr thePOA,
+                                      ::SMESH_Gen*            theGenImpl );
+  // Destructor
+  virtual ~StdMeshers_PolyhedronPerSolid_3D_i();
+};
+
+#endif
index 602b00ae2fb18923ec64f2f3fcf1b4d1b2444431..676de466d35cb1a9a6073c495f219d00e8e64650 100644 (file)
@@ -57,6 +57,7 @@
 #include "StdMeshers_NumberOfLayers_i.hxx"
 #include "StdMeshers_NumberOfSegments_i.hxx"
 #include "StdMeshers_PolygonPerFace_2D_i.hxx"
+#include "StdMeshers_PolyhedronPerSolid_3D_i.hxx"
 #include "StdMeshers_Prism_3D_i.hxx"
 #include "StdMeshers_ProjectionSource1D_i.hxx"
 #include "StdMeshers_ProjectionSource2D_i.hxx"
@@ -250,6 +251,8 @@ STDMESHERS_I_EXPORT
       aCreator = new StdHypothesisCreator_i<StdMeshers_Cartesian_3D_i>;
     else if (strcmp(aHypName, "PolygonPerFace_2D") == 0)
       aCreator = new StdHypothesisCreator_i<StdMeshers_PolygonPerFace_2D_i>;
+    else if (strcmp(aHypName, "PolyhedronPerSolid_3D") == 0)
+      aCreator = new StdHypothesisCreator_i<StdMeshers_PolyhedronPerSolid_3D_i>;
     else ;
 
     return aCreator;