#include "SMESH_ObjectDef.h"#include "SMESH_ActorUtils.h"#include "SMDS_Mesh.hxx"#include "SMDS_PolyhedralVolumeOfNodes.hxx"#include "SMESH_Actor.h"#include "SMESH_ControlsDef.hxx"#include "SalomeApp_Application.h"#include "VTKViewer_ExtractUnstructuredGrid.h"#include "VTKViewer_CellLocationsArray.h"#include <CORBA_SERVER_HEADER(SMESH_Gen)>#include <CORBA_SERVER_HEADER(SALOME_Exception)>#include <vtkCell.h>#include <vtkIdList.h>#include <vtkCellArray.h>#include <vtkUnsignedCharArray.h>#include <vtkUnstructuredGrid.h>#include <memory>#include <sstream>#include <stdexcept>#include <set>#include "utilities.h"
Go to the source code of this file.
Defines | |
| #define | EXCEPTION(TYPE, MSG) |
Functions | |
| static vtkIdType | getCellType (const SMDSAbs_ElementType theType, const bool thePoly, const int theNbNodes) |
| static int | getNodesFromElems (SMESH::long_array_var &theElemIds, const SMDS_Mesh *theMesh, std::list< const SMDS_MeshElement * > &theResList) |
| static int | getPointers (const SMDSAbs_ElementType theRequestType, SMESH::long_array_var &theElemIds, const SMDS_Mesh *theMesh, std::list< const SMDS_MeshElement * > &theResList) |
Variables | |
| static int | MYDEBUG = 0 |
| static int | MYDEBUGWITHFILES = 0 |
| #define EXCEPTION | ( | TYPE, | |
| MSG | |||
| ) |
{\
std::ostringstream aStream;\
aStream<<__FILE__<<"["<<__LINE__<<"]::"<<MSG;\
throw TYPE(aStream.str());\
}
Definition at line 59 of file SMESH_Object.cxx.
| static vtkIdType getCellType | ( | const SMDSAbs_ElementType | theType, |
| const bool | thePoly, | ||
| const int | theNbNodes | ||
| ) | [static] |
Definition at line 84 of file SMESH_Object.cxx.
References SMDSAbs_0DElement, SMDSAbs_Edge, SMDSAbs_Face, and SMDSAbs_Volume.
Referenced by SMESH_VisualObjDef.buildElemPrs().
{
switch( theType )
{
case SMDSAbs_0DElement:
return VTK_VERTEX;
case SMDSAbs_Edge:
if( theNbNodes == 2 ) return VTK_LINE;
else if ( theNbNodes == 3 ) return VTK_QUADRATIC_EDGE;
else return VTK_EMPTY_CELL;
case SMDSAbs_Face :
if (thePoly && theNbNodes>2 ) return VTK_POLYGON;
else if ( theNbNodes == 3 ) return VTK_TRIANGLE;
else if ( theNbNodes == 4 ) return VTK_QUAD;
else if ( theNbNodes == 6 ) return VTK_QUADRATIC_TRIANGLE;
else if ( theNbNodes == 8 ) return VTK_QUADRATIC_QUAD;
else return VTK_EMPTY_CELL;
case SMDSAbs_Volume:
if (thePoly && theNbNodes>3 ) return VTK_POLYHEDRON; //VTK_CONVEX_POINT_SET;
else if ( theNbNodes == 4 ) return VTK_TETRA;
else if ( theNbNodes == 5 ) return VTK_PYRAMID;
else if ( theNbNodes == 6 ) return VTK_WEDGE;
else if ( theNbNodes == 8 ) return VTK_HEXAHEDRON;
else if ( theNbNodes == 10 ) {
return VTK_QUADRATIC_TETRA;
}
else if ( theNbNodes == 20 ) {
return VTK_QUADRATIC_HEXAHEDRON;
}
else if ( theNbNodes == 15 ) {
return VTK_QUADRATIC_WEDGE;
}
else if ( theNbNodes==13 ) {
return VTK_QUADRATIC_PYRAMID; //VTK_CONVEX_POINT_SET;
}
else return VTK_EMPTY_CELL;
default: return VTK_EMPTY_CELL;
}
}
| static int getNodesFromElems | ( | SMESH::long_array_var & | theElemIds, |
| const SMDS_Mesh * | theMesh, | ||
| std::list< const SMDS_MeshElement * > & | theResList | ||
| ) | [static] |
Definition at line 871 of file SMESH_Object.cxx.
References SMDS_Mesh.FindElement(), and SMDS_MeshElement.nodesIterator().
Referenced by SMESH_subMeshObj.GetEntities(), and SMESH_GroupObj.GetEntities().
{
set<const SMDS_MeshElement*> aNodeSet;
for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
{
const SMDS_MeshElement* anElem = theMesh->FindElement( theElemIds[ i ] );
if ( anElem != 0 )
{
SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
while ( anIter->more() )
{
const SMDS_MeshElement* aNode = anIter->next();
if ( aNode != 0 )
aNodeSet.insert( aNode );
}
}
}
set<const SMDS_MeshElement*>::const_iterator anIter;
for ( anIter = aNodeSet.begin(); anIter != aNodeSet.end(); ++anIter )
theResList.push_back( *anIter );
return theResList.size();
}
| static int getPointers | ( | const SMDSAbs_ElementType | theRequestType, |
| SMESH::long_array_var & | theElemIds, | ||
| const SMDS_Mesh * | theMesh, | ||
| std::list< const SMDS_MeshElement * > & | theResList | ||
| ) | [static] |
Definition at line 903 of file SMESH_Object.cxx.
References SMDS_Mesh.FindElement(), SMDS_Mesh.FindNode(), and SMDSAbs_Node.
Referenced by SMESH_subMeshObj.GetEntities(), and SMESH_GroupObj.GetEntities().
{
for ( CORBA::Long i = 0, n = theElemIds->length(); i < n; i++ )
{
const SMDS_MeshElement* anElem = theRequestType == SMDSAbs_Node
? theMesh->FindNode( theElemIds[ i ] ) : theMesh->FindElement( theElemIds[ i ] );
if ( anElem != 0 )
theResList.push_back( anElem );
}
return theResList.size();
}
Definition at line 70 of file SMESH_Object.cxx.
Referenced by SMESH_VisualObjDef.buildElemPrs(), SMESH_GroupObj.SMESH_GroupObj(), SMESH_MeshObj.SMESH_MeshObj(), SMESH_subMeshObj.SMESH_subMeshObj(), SMESH_SubMeshObj.SMESH_SubMeshObj(), SMESH_GroupObj.~SMESH_GroupObj(), SMESH_MeshObj.~SMESH_MeshObj(), and SMESH_subMeshObj.~SMESH_subMeshObj().
int MYDEBUGWITHFILES = 0 [static] |
Definition at line 71 of file SMESH_Object.cxx.
Referenced by SMESH_VisualObjDef.buildPrs().