X-Git-Url: http://git.salome-platform.org/gitweb/?a=blobdiff_plain;f=src%2FSMDS%2FSMDS_UnstructuredGrid.cxx;h=1eb5c445d60e5e13ed1af797359a7f258ceb53d5;hb=d627437f5d8c20e0750973df76c62d55fc19a5e7;hp=099b8cbe85235b332d2cdfb89df257bde51bfde0;hpb=d4a710ce52f6e76786a7b3845e2f7975dc9a00b1;p=modules%2Fsmesh.git diff --git a/src/SMDS/SMDS_UnstructuredGrid.cxx b/src/SMDS/SMDS_UnstructuredGrid.cxx index 099b8cbe8..1eb5c445d 100644 --- a/src/SMDS/SMDS_UnstructuredGrid.cxx +++ b/src/SMDS/SMDS_UnstructuredGrid.cxx @@ -27,7 +27,9 @@ #include "utilities.h" #include +#include #include +#include #include #include @@ -218,10 +220,12 @@ void SMDS_UnstructuredGrid::compactGrid(std::vector& idNodesOldToNew, int n ++i; int endBloc = i; if ( endBloc > startBloc ) - copyBloc(newTypes, idCellsOldToNew, idNodesOldToNew, newConnectivity, newLocations, pointsCell, - alreadyCopied, startBloc, endBloc); + copyBloc(newTypes, + idCellsOldToNew, idNodesOldToNew, + newConnectivity, newLocations, + pointsCell, alreadyCopied, + startBloc, endBloc); } - newConnectivity->Squeeze(); if (1/*newNodeSize*/) @@ -231,6 +235,19 @@ void SMDS_UnstructuredGrid::compactGrid(std::vector& idNodesOldToNew, int n 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) { vtkIdTypeArray *newFaceLocations = vtkIdTypeArray::New(); @@ -275,7 +292,9 @@ void SMDS_UnstructuredGrid::compactGrid(std::vector& idNodesOldToNew, int n newFaces->Delete(); } else + { this->SetCells(newTypes, newLocations, newConnectivity, FaceLocations, Faces); + } newPoints->Delete(); newTypes->Delete(); @@ -299,10 +318,15 @@ void SMDS_UnstructuredGrid::copyNodes(vtkPoints *newPoints, std::vector& id } } -void SMDS_UnstructuredGrid::copyBloc(vtkUnsignedCharArray *newTypes, std::vector& idCellsOldToNew, - std::vector& idNodesOldToNew, vtkCellArray* newConnectivity, - vtkIdTypeArray* newLocations, vtkIdType* pointsCell, int& alreadyCopied, - int start, int end) +void SMDS_UnstructuredGrid::copyBloc(vtkUnsignedCharArray *newTypes, + std::vector& idCellsOldToNew, + std::vector& idNodesOldToNew, + vtkCellArray* newConnectivity, + vtkIdTypeArray* newLocations, + vtkIdType* pointsCell, + int& alreadyCopied, + int start, + int end) { MESSAGE("copyBloc " << alreadyCopied << " " << start << " " << end << " size: " << end - start << " total: " << alreadyCopied + end - start); for (int j = start; j < end; j++) @@ -739,7 +763,7 @@ void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool withEdges) * @param vtkId the vtk id of the cell * @return number of neighbors */ -int SMDS_UnstructuredGrid::GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId) +int SMDS_UnstructuredGrid::GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId, bool getSkin) { int vtkType = this->GetCellType(vtkId); int cellDim = SMDS_Downward::getCellDimension(vtkType); @@ -774,9 +798,27 @@ int SMDS_UnstructuredGrid::GetNeighbors(int* neighborsVtkIds, int* downIds, unsi downIds[nb] = downId; downTypes[nb] = cellType; nb++; + if (nb >= NBMAXNEIGHBORS) + { + INFOS("SMDS_UnstructuredGrid::GetNeighbors problem: NBMAXNEIGHBORS=" <getVtkCellId(downId); // OK if skin present + downIds[nb] = downId; + downTypes[nb] = cellType; + nb++; + if (nb >= NBMAXNEIGHBORS) + { + INFOS("SMDS_UnstructuredGrid::GetNeighbors problem: NBMAXNEIGHBORS=" <= NBMAXNEIGHBORS) - assert(0); } return nb; } @@ -1072,3 +1114,50 @@ SMDS_MeshCell* SMDS_UnstructuredGrid::extrudeVolumeFromFace(int vtkVolId, // TODO update sub-shape list of elements and nodes 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; +} +