0021459: EDF 1495 SMESH: Manipulation of discrete elements with attributes

+  void AllocateDiameters( vtkIdType maxVtkID );
+  void SetBallDiameter( vtkIdType vtkID, double diameter );
+  double GetBallDiameter( vtkIdType vtkID ) const;
This commit is contained in:
eap 2012-07-19 13:14:28 +00:00
parent f738cc36e8
commit 9f74998b92
2 changed files with 88 additions and 11 deletions

View File

@ -27,7 +27,9 @@
#include "utilities.h" #include "utilities.h"
#include <vtkCellArray.h> #include <vtkCellArray.h>
#include <vtkCellData.h>
#include <vtkCellLinks.h> #include <vtkCellLinks.h>
#include <vtkDoubleArray.h>
#include <vtkIdTypeArray.h> #include <vtkIdTypeArray.h>
#include <vtkUnsignedCharArray.h> #include <vtkUnsignedCharArray.h>
@ -218,10 +220,12 @@ void SMDS_UnstructuredGrid::compactGrid(std::vector<int>& idNodesOldToNew, int n
++i; ++i;
int endBloc = i; int endBloc = i;
if ( endBloc > startBloc ) if ( endBloc > startBloc )
copyBloc(newTypes, idCellsOldToNew, idNodesOldToNew, newConnectivity, newLocations, pointsCell, copyBloc(newTypes,
alreadyCopied, startBloc, endBloc); idCellsOldToNew, idNodesOldToNew,
newConnectivity, newLocations,
pointsCell, alreadyCopied,
startBloc, endBloc);
} }
newConnectivity->Squeeze(); newConnectivity->Squeeze();
if (1/*newNodeSize*/) if (1/*newNodeSize*/)
@ -231,6 +235,19 @@ void SMDS_UnstructuredGrid::compactGrid(std::vector<int>& idNodesOldToNew, int n
MESSAGE("NumberOfPoints: " << this->GetNumberOfPoints()); MESSAGE("NumberOfPoints: " << this->GetNumberOfPoints());
} }
if (vtkDoubleArray* diameters =
vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() )) // Balls
{
for (int oldCellID = 0; oldCellID < oldCellSize; oldCellID++)
{
if (this->Types->GetValue(oldCellID) == VTK_EMPTY_CELL)
continue;
int newCellId = idCellsOldToNew[ oldCellID ];
if (newTypes->GetValue(newCellId) == VTK_POLY_VERTEX)
diameters->SetValue( newCellId, diameters->GetValue( oldCellID ));
}
}
if (this->FaceLocations) if (this->FaceLocations)
{ {
vtkIdTypeArray *newFaceLocations = vtkIdTypeArray::New(); vtkIdTypeArray *newFaceLocations = vtkIdTypeArray::New();
@ -275,7 +292,9 @@ void SMDS_UnstructuredGrid::compactGrid(std::vector<int>& idNodesOldToNew, int n
newFaces->Delete(); newFaces->Delete();
} }
else else
{
this->SetCells(newTypes, newLocations, newConnectivity, FaceLocations, Faces); this->SetCells(newTypes, newLocations, newConnectivity, FaceLocations, Faces);
}
newPoints->Delete(); newPoints->Delete();
newTypes->Delete(); newTypes->Delete();
@ -299,10 +318,15 @@ void SMDS_UnstructuredGrid::copyNodes(vtkPoints *newPoints, std::vector<int>& id
} }
} }
void SMDS_UnstructuredGrid::copyBloc(vtkUnsignedCharArray *newTypes, std::vector<int>& idCellsOldToNew, void SMDS_UnstructuredGrid::copyBloc(vtkUnsignedCharArray *newTypes,
std::vector<int>& idNodesOldToNew, vtkCellArray* newConnectivity, std::vector<int>& idCellsOldToNew,
vtkIdTypeArray* newLocations, vtkIdType* pointsCell, int& alreadyCopied, std::vector<int>& idNodesOldToNew,
int start, int end) vtkCellArray* newConnectivity,
vtkIdTypeArray* newLocations,
vtkIdType* pointsCell,
int& alreadyCopied,
int start,
int end)
{ {
MESSAGE("copyBloc " << alreadyCopied << " " << start << " " << end << " size: " << end - start << " total: " << alreadyCopied + end - start); MESSAGE("copyBloc " << alreadyCopied << " " << start << " " << end << " size: " << end - start << " total: " << alreadyCopied + end - start);
for (int j = start; j < end; j++) for (int j = start; j < end; j++)
@ -1072,3 +1096,50 @@ SMDS_MeshCell* SMDS_UnstructuredGrid::extrudeVolumeFromFace(int vtkVolId,
// TODO update sub-shape list of elements and nodes // TODO update sub-shape list of elements and nodes
return 0; return 0;
} }
//================================================================================
/*!
* \brief Allocates data array for ball diameters
* \param MaxVtkID - max ID of a ball element
*/
//================================================================================
void SMDS_UnstructuredGrid::AllocateDiameters( vtkIdType MaxVtkID )
{
SetBallDiameter( MaxVtkID, 0 );
}
//================================================================================
/*!
* \brief Sets diameter of a ball element
* \param vtkID - vtk id of the ball element
* \param diameter - diameter of the ball element
*/
//================================================================================
void SMDS_UnstructuredGrid::SetBallDiameter( vtkIdType vtkID, double diameter )
{
vtkDoubleArray* array = vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() );
if ( !array )
{
array = vtkDoubleArray::New();
array->SetNumberOfComponents(1);
vtkDataSet::CellData->SetScalars( array );
}
array->InsertValue( vtkID, diameter );
}
//================================================================================
/*!
* \brief Returns diameter of a ball element
* \param vtkID - vtk id of the ball element
*/
//================================================================================
double SMDS_UnstructuredGrid::GetBallDiameter( vtkIdType vtkID ) const
{
if ( vtkDataSet::CellData )
return vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() )->GetValue( vtkID );
return 0;
}

View File

@ -66,9 +66,10 @@ class SMDS_EXPORT SMDS_UnstructuredGrid: public vtkUnstructuredGrid
{ {
public: public:
void setSMDS_mesh(SMDS_Mesh *mesh); void setSMDS_mesh(SMDS_Mesh *mesh);
void compactGrid(std::vector<int>& idNodesOldToNew, int newNodeSize, std::vector<int>& idCellsOldToNew, void compactGrid(std::vector<int>& idNodesOldToNew,
int newCellSize); int newNodeSize,
std::vector<int>& idCellsOldToNew,
int newCellSize);
virtual unsigned long GetMTime(); virtual unsigned long GetMTime();
virtual void Update(); virtual void Update();
virtual void UpdateInformation(); virtual void UpdateInformation();
@ -89,7 +90,8 @@ public:
void ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds); void ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds);
int getOrderedNodesOfFace(int vtkVolId, int& dim, std::vector<vtkIdType>& orderedNodes); int getOrderedNodesOfFace(int vtkVolId, int& dim, std::vector<vtkIdType>& orderedNodes);
void BuildLinks(); void BuildLinks();
SMDS_MeshCell* extrudeVolumeFromFace(int vtkVolId, int domain1, int domain2, std::set<int>& originalNodes, SMDS_MeshCell* extrudeVolumeFromFace(int vtkVolId, int domain1, int domain2,
std::set<int>& originalNodes,
std::map<int, std::map<int, int> >& nodeDomains, std::map<int, std::map<int, int> >& nodeDomains,
std::map<int, std::map<long,int> >& nodeQuadDomains); std::map<int, std::map<long,int> >& nodeQuadDomains);
vtkCellLinks* GetLinks() vtkCellLinks* GetLinks()
@ -100,6 +102,10 @@ public:
{ {
return _downArray[vtkType]; return _downArray[vtkType];
} }
void AllocateDiameters( vtkIdType maxVtkID );
void SetBallDiameter( vtkIdType vtkID, double diameter );
double GetBallDiameter( vtkIdType vtkID ) const;
static SMDS_UnstructuredGrid* New(); static SMDS_UnstructuredGrid* New();
SMDS_Mesh *_mesh; SMDS_Mesh *_mesh;
protected: protected: