1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410
|
/*=========================================================================
Program: Visualization Toolkit
Module: vtkGenericCutter.cxx
Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
All rights reserved.
See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notice for more information.
=========================================================================*/
#include "vtkGenericCutter.h"
#include "vtkCellArray.h"
#include "vtkCellData.h"
#include "vtkContourValues.h"
#include "vtkDataSet.h"
#include "vtkDoubleArray.h"
#include "vtkGenericCell.h"
#include "vtkImplicitFunction.h"
#include "vtkMergePoints.h"
#include "vtkInformation.h"
#include "vtkInformationVector.h"
#include "vtkObjectFactory.h"
#include "vtkPointData.h"
#include "vtkPolyData.h"
#include "vtkUnstructuredGrid.h"
#include "vtkPoints.h"
#include "vtkGenericCellIterator.h"
#include "vtkGenericAdaptorCell.h"
#include "vtkGenericPointIterator.h"
#include "vtkGenericDataSet.h"
#include "vtkGenericAttributeCollection.h"
#include "vtkGenericAttribute.h"
#include "vtkGenericCellTessellator.h"
#include "vtkIncrementalPointLocator.h"
#include <cmath>
vtkStandardNewMacro(vtkGenericCutter);
vtkCxxSetObjectMacro(vtkGenericCutter,CutFunction,vtkImplicitFunction);
vtkCxxSetObjectMacro(vtkGenericCutter,Locator,vtkIncrementalPointLocator);
//----------------------------------------------------------------------------
// Construct with user-specified implicit function; initial value of 0.0; and
// generating cut scalars turned off.
//
vtkGenericCutter::vtkGenericCutter(vtkImplicitFunction *cf)
{
this->ContourValues = vtkContourValues::New();
this->CutFunction = cf;
this->GenerateCutScalars = 0;
this->Locator = NULL;
this->InternalPD = vtkPointData::New();
this->SecondaryPD = vtkPointData::New();
this->SecondaryCD = vtkCellData::New();
}
//----------------------------------------------------------------------------
vtkGenericCutter::~vtkGenericCutter()
{
this->ContourValues->Delete();
this->SetCutFunction(NULL);
if ( this->Locator )
{
this->Locator->UnRegister(this);
this->Locator = NULL;
}
this->InternalPD->Delete();
this->SecondaryPD->Delete();
this->SecondaryCD->Delete();
}
//----------------------------------------------------------------------------
// Description:
// Set a particular contour value at contour number i. The index i ranges
// between 0<=i<NumberOfContours.
void vtkGenericCutter::SetValue(int i, double value)
{
this->ContourValues->SetValue(i,value);
}
//-----------------------------------------------------------------------------
// Description:
// Get the ith contour value.
double vtkGenericCutter::GetValue(int i)
{
return this->ContourValues->GetValue(i);
}
//-----------------------------------------------------------------------------
// Description:
// Get a pointer to an array of contour values. There will be
// GetNumberOfContours() values in the list.
double *vtkGenericCutter::GetValues()
{
return this->ContourValues->GetValues();
}
//-----------------------------------------------------------------------------
// Description:
// Fill a supplied list with contour values. There will be
// GetNumberOfContours() values in the list. Make sure you allocate
// enough memory to hold the list.
void vtkGenericCutter::GetValues(double *contourValues)
{
this->ContourValues->GetValues(contourValues);
}
//-----------------------------------------------------------------------------
// Description:
// Set the number of contours to place into the list. You only really
// need to use this method to reduce list size. The method SetValue()
// will automatically increase list size as needed.
void vtkGenericCutter::SetNumberOfContours(int number)
{
this->ContourValues->SetNumberOfContours(number);
}
//-----------------------------------------------------------------------------
// Description:
// Get the number of contours in the list of contour values.
int vtkGenericCutter::GetNumberOfContours()
{
return this->ContourValues->GetNumberOfContours();
}
//-----------------------------------------------------------------------------
// Description:
// Generate numContours equally spaced contour values between specified
// range. Contour values will include min/max range values.
void vtkGenericCutter::GenerateValues(int numContours, double range[2])
{
this->ContourValues->GenerateValues(numContours, range);
}
// Description:
// Generate numContours equally spaced contour values between specified
// range. Contour values will include min/max range values.
void vtkGenericCutter::GenerateValues(int numContours, double rangeStart,
double rangeEnd)
{
this->ContourValues->GenerateValues(numContours, rangeStart, rangeEnd);
}
//----------------------------------------------------------------------------
// Overload standard modified time function. If cut functions is modified,
// or contour values modified, then this object is modified as well.
//
vtkMTimeType vtkGenericCutter::GetMTime()
{
vtkMTimeType mTime = this->Superclass::GetMTime();
vtkMTimeType contourValuesMTime = this->ContourValues->GetMTime();
vtkMTimeType time;
mTime = ( contourValuesMTime > mTime ? contourValuesMTime : mTime );
if ( this->CutFunction != NULL )
{
time = this->CutFunction->GetMTime();
mTime = ( time > mTime ? time : mTime );
}
if ( this->Locator != NULL )
{
time = this->Locator->GetMTime();
mTime = ( time > mTime ? time : mTime );
}
return mTime;
}
//----------------------------------------------------------------------------
// Cut through data generating surface.
//
int vtkGenericCutter::RequestData(
vtkInformation *vtkNotUsed(request),
vtkInformationVector **inputVector,
vtkInformationVector *outputVector)
{
// get the info objects
vtkInformation *inInfo = inputVector[0]->GetInformationObject(0);
vtkInformation *outInfo = outputVector->GetInformationObject(0);
// get the input and output
vtkGenericDataSet *input = vtkGenericDataSet::SafeDownCast(
inInfo->Get(vtkDataObject::DATA_OBJECT()));
vtkPolyData *output = vtkPolyData::SafeDownCast(
outInfo->Get(vtkDataObject::DATA_OBJECT()));
vtkDebugMacro(<< "Executing cutter");
if (input==0)
{
vtkErrorMacro("No input specified");
return 1;
}
if (this->CutFunction==0)
{
vtkErrorMacro("No cut function specified");
return 1;
}
if ( input->GetNumberOfPoints()<1 )
{
vtkErrorMacro("Input data set is empty");
return 1;
}
vtkPointData *outPd = output->GetPointData();
vtkCellData *outCd = output->GetCellData();
// Create objects to hold output of contour operation
//
vtkIdType numCells=input->GetNumberOfCells();
int numContours = this->ContourValues->GetNumberOfContours();
vtkIdType estimatedSize = static_cast<vtkIdType>(
pow(static_cast<double>(numCells), .75)) * numContours;
estimatedSize = estimatedSize / 1024 * 1024; //multiple of 1024
if (estimatedSize < 1024)
{
estimatedSize = 1024;
}
vtkPoints *newPts = vtkPoints::New();
newPts->Allocate(estimatedSize,estimatedSize);
vtkCellArray *newVerts = vtkCellArray::New();
newVerts->Allocate(estimatedSize,estimatedSize);
vtkCellArray *newLines = vtkCellArray::New();
newLines->Allocate(estimatedSize,estimatedSize);
vtkCellArray *newPolys = vtkCellArray::New();
newPolys->Allocate(estimatedSize,estimatedSize);
output->Allocate(numCells);
// locator used to merge potentially duplicate points
if(this->Locator==0)
{
this->CreateDefaultLocator();
}
this->Locator->InitPointInsertion(newPts,input->GetBounds(),estimatedSize);
// prepare the output attributes
vtkGenericAttributeCollection *attributes=input->GetAttributes();
vtkGenericAttribute *attribute;
vtkDataArray *attributeArray;
int c=attributes->GetNumberOfAttributes();
vtkDataSetAttributes *secondaryAttributes;
int attributeType;
for(vtkIdType i = 0; i<c; ++i)
{
attribute = attributes->GetAttribute(i);
attributeType = attribute->GetType();
if(attribute->GetCentering() == vtkPointCentered)
{
secondaryAttributes = this->SecondaryPD;
attributeArray=vtkDataArray::CreateDataArray(attribute->GetComponentType());
attributeArray->SetNumberOfComponents(attribute->GetNumberOfComponents());
attributeArray->SetName(attribute->GetName());
this->InternalPD->AddArray(attributeArray);
attributeArray->Delete();
if(this->InternalPD->GetAttribute(attributeType)==0)
{
this->InternalPD->SetActiveAttribute(
this->InternalPD->GetNumberOfArrays()-1,attributeType);
}
}
else // vtkCellCentered
{
secondaryAttributes = this->SecondaryCD;
}
attributeArray = vtkDataArray::CreateDataArray(attribute->GetComponentType());
attributeArray->SetNumberOfComponents(attribute->GetNumberOfComponents());
attributeArray->SetName(attribute->GetName());
secondaryAttributes->AddArray(attributeArray);
attributeArray->Delete();
if(secondaryAttributes->GetAttribute(attributeType)==0)
{
secondaryAttributes->SetActiveAttribute(secondaryAttributes->GetNumberOfArrays()-1,
attributeType);
}
}
outPd->InterpolateAllocate(this->SecondaryPD,estimatedSize,estimatedSize);
outCd->CopyAllocate(this->SecondaryCD,estimatedSize,estimatedSize);
vtkGenericAdaptorCell *cell;
//----------- Begin of contouring algorithm --------------------//
vtkGenericCellIterator *cellIt = input->NewCellIterator();
vtkIdType updateCount = numCells/20 + 1; // update roughly every 5%
vtkIdType count = 0;
int abortExecute=0;
input->GetTessellator()->InitErrorMetrics(input);
for(cellIt->Begin(); !cellIt->IsAtEnd() && !abortExecute; cellIt->Next())
{
if ( !(count % updateCount) )
{
this->UpdateProgress(static_cast<double>(count) / numCells);
abortExecute = this->GetAbortExecute();
}
cell = cellIt->GetCell();
cell->Contour(this->ContourValues, this->CutFunction, input->GetAttributes(),
input->GetTessellator(),
this->Locator, newVerts, newLines, newPolys, outPd, outCd,
this->InternalPD,this->SecondaryPD,this->SecondaryCD);
++count;
} // for each cell
cellIt->Delete();
vtkDebugMacro(<<"Created: "
<< newPts->GetNumberOfPoints() << " points, "
<< newVerts->GetNumberOfCells() << " verts, "
<< newLines->GetNumberOfCells() << " lines, "
<< newPolys->GetNumberOfCells() << " triangles");
//----------- End of contouring algorithm ----------------------//
// Update ourselves. Because we don't know up front how many verts, lines,
// polys we've created, take care to reclaim memory.
//
output->SetPoints(newPts);
newPts->Delete();
if (newVerts->GetNumberOfCells()>0)
{
output->SetVerts(newVerts);
}
newVerts->Delete();
if (newLines->GetNumberOfCells()>0)
{
output->SetLines(newLines);
}
newLines->Delete();
if (newPolys->GetNumberOfCells()>0)
{
output->SetPolys(newPolys);
}
newPolys->Delete();
this->Locator->Initialize();//releases leftover memory
output->Squeeze();
return 1;
}
//----------------------------------------------------------------------------
// Specify a spatial locator for merging points. By default,
// an instance of vtkMergePoints is used.
void vtkGenericCutter::CreateDefaultLocator()
{
if ( this->Locator == NULL )
{
this->Locator = vtkMergePoints::New();
this->Locator->Register(this);
this->Locator->Delete();
}
}
//----------------------------------------------------------------------------
void vtkGenericCutter::PrintSelf(ostream& os, vtkIndent indent)
{
this->Superclass::PrintSelf(os,indent);
os << indent << "Cut Function: " << this->CutFunction << "\n";
if ( this->Locator )
{
os << indent << "Locator: " << this->Locator << "\n";
}
else
{
os << indent << "Locator: (none)\n";
}
this->ContourValues->PrintSelf(os,indent.GetNextIndent());
os << indent << "Generate Cut Scalars: "
<< (this->GenerateCutScalars ? "On\n" : "Off\n");
}
//----------------------------------------------------------------------------
int vtkGenericCutter::FillInputPortInformation(int port,
vtkInformation* info)
{
if(!this->Superclass::FillInputPortInformation(port, info))
{
return 0;
}
info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkGenericDataSet");
return 1;
}
|