Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
vtkTensorGlyphColor Class Reference

scale and orient glyph(s) according to eigenvalues and eigenvectors of symmetrical part of tensor More...

#include <vtkTensorGlyphColor.h>

+ Inheritance diagram for vtkTensorGlyphColor:

Public Types

enum  { COLOR_BY_SCALARS , COLOR_BY_EIGENVALUES }
 

Public Member Functions

 vtkTypeMacro (vtkTensorGlyphColor, vtkPolyDataAlgorithm) void PrintSelf(ostream &os
 
void SetSourceData (vtkPolyData *source)
 
vtkPolyData * GetSource ()
 
void SetSourceConnection (int id, vtkAlgorithmOutput *algOutput)
 
void SetSourceConnection (vtkAlgorithmOutput *algOutput)
 
 vtkSetMacro (Scaling, vtkTypeBool)
 
 vtkGetMacro (Scaling, vtkTypeBool)
 
 vtkBooleanMacro (Scaling, vtkTypeBool)
 
 vtkSetMacro (ScaleFactor, double)
 
 vtkGetMacro (ScaleFactor, double)
 
 vtkSetMacro (ThreeGlyphs, vtkTypeBool)
 
 vtkGetMacro (ThreeGlyphs, vtkTypeBool)
 
 vtkBooleanMacro (ThreeGlyphs, vtkTypeBool)
 
 vtkSetMacro (Symmetric, vtkTypeBool)
 
 vtkGetMacro (Symmetric, vtkTypeBool)
 
 vtkBooleanMacro (Symmetric, vtkTypeBool)
 
 vtkSetMacro (Length, double)
 
 vtkGetMacro (Length, double)
 
 vtkSetMacro (ExtractEigenvalues, vtkTypeBool)
 
 vtkBooleanMacro (ExtractEigenvalues, vtkTypeBool)
 
 vtkGetMacro (ExtractEigenvalues, vtkTypeBool)
 
 vtkSetMacro (ColorGlyphs, vtkTypeBool)
 
 vtkGetMacro (ColorGlyphs, vtkTypeBool)
 
 vtkBooleanMacro (ColorGlyphs, vtkTypeBool)
 
 vtkSetClampMacro (ColorMode, int, COLOR_BY_SCALARS, COLOR_BY_EIGENVALUES)
 
 vtkGetMacro (ColorMode, int)
 
void SetColorModeToScalars ()
 
void SetColorModeToEigenvalues ()
 
 vtkSetMacro (ClampScaling, vtkTypeBool)
 
 vtkGetMacro (ClampScaling, vtkTypeBool)
 
 vtkBooleanMacro (ClampScaling, vtkTypeBool)
 
 vtkSetMacro (MaxScaleFactor, double)
 
 vtkGetMacro (MaxScaleFactor, double)
 

Static Public Member Functions

static vtkTensorGlyphColorNew ()
 

Public Attributes

vtkIndent indent override
 

Protected Member Functions

 vtkTensorGlyphColor ()
 
 ~vtkTensorGlyphColor () override
 
int RequestUpdateExtent (vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
 
int RequestData (vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
 
int FillInputPortInformation (int port, vtkInformation *info) override
 

Protected Attributes

vtkTypeBool Scaling
 
double ScaleFactor
 
vtkTypeBool ExtractEigenvalues
 
vtkTypeBool ColorGlyphs
 
int ColorMode
 
vtkTypeBool ClampScaling
 
double MaxScaleFactor
 
vtkTypeBool ThreeGlyphs
 
vtkTypeBool Symmetric
 
double Length
 

Detailed Description

scale and orient glyph(s) according to eigenvalues and eigenvectors of symmetrical part of tensor

vtkTensorGlyphColor is a filter that copies a geometric representation (specified as polygonal data) to every input point. The geometric representation, or glyph, can be scaled and/or rotated according to the tensor at the input point. Scaling and rotation is controlled by the eigenvalues/eigenvectors of the symmetrical part of the tensor as follows: For each tensor, the eigenvalues (and associated eigenvectors) are sorted to determine the major, medium, and minor eigenvalues/eigenvectors. The eigenvalue decomposition only makes sense for symmetric tensors, hence the need to only consider the symmetric part of the tensor, which is 1/2 (T + T.transposed()).

If the boolean variable ThreeGlyphs is not set the major eigenvalue scales the glyph in the x-direction, the medium in the y-direction, and the minor in the z-direction. Then, the glyph is rotated so that the glyph's local x-axis lies along the major eigenvector, y-axis along the medium eigenvector, and z-axis along the minor.

If the boolean variable ThreeGlyphs is set three glyphs are produced, each of them oriented along an eigenvector and scaled according to the corresponding eigenvector.

If the boolean variable Symmetric is set each glyph is mirrored (2 or 6 glyphs will be produced)

The x-axis of the source glyph will correspond to the eigenvector on output. Point (0,0,0) in the source will be placed in the data point. Variable Length will normally correspond to the distance from the origin to the tip of the source glyph along the x-axis, but can be changed to produce other results when Symmetric is on, e.g. glyphs that do not touch or that overlap.

Please note that when Symmetric is false it will generally be better to place the source glyph from (-0.5,0,0) to (0.5,0,0), i.e. centred at the origin. When symmetric is true the placement from (0,0,0) to (1,0,0) will generally be more convenient.

A scale factor is provided to control the amount of scaling. Also, you can turn off scaling completely if desired. The boolean variable ClampScaling controls the maximum scaling (in conjunction with MaxScaleFactor.) This is useful in certain applications where singularities or large order of magnitude differences exist in the eigenvalues.

If the boolean variable ColorGlyphs is set to true the glyphs are colored. The glyphs can be colored using the input scalars (SetColorModeToScalars), which is the default, or colored using the eigenvalues (SetColorModeToEigenvalues).

Another instance variable, ExtractEigenvalues, has been provided to control extraction of eigenvalues/eigenvectors. If this boolean is false, then eigenvalues/eigenvectors are not extracted, and the columns of the tensor are taken as the eigenvectors (the norm of column, always positive, is the eigenvalue). This allows additional capability over the vtkGlyph3D object. That is, the glyph can be oriented in three directions instead of one.

Thanks:
Thanks to Jose Paulo Moitinho de Almeida for enhancements.
See also
vtkGlyph3D vtkPointLoad vtkHyperStreamline

Definition at line 89 of file vtkTensorGlyphColor.h.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
Enumerator
COLOR_BY_SCALARS 
COLOR_BY_EIGENVALUES 

Definition at line 191 of file vtkTensorGlyphColor.h.

Constructor & Destructor Documentation

◆ vtkTensorGlyphColor()

vtkTensorGlyphColor::vtkTensorGlyphColor ( )
protected

◆ ~vtkTensorGlyphColor()

vtkTensorGlyphColor::~vtkTensorGlyphColor ( )
overrideprotecteddefault

Member Function Documentation

◆ FillInputPortInformation()

int vtkTensorGlyphColor::FillInputPortInformation ( int  port,
vtkInformation *  info 
)
overrideprotected

Definition at line 538 of file vtkTensorGlyphColor.cxx.

539{
540 if (port == 1)
541 {
542 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkPolyData");
543 return 1;
544 }
545 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
546 return 1;
547}

◆ GetSource()

vtkPolyData * vtkTensorGlyphColor::GetSource ( )

Definition at line 528 of file vtkTensorGlyphColor.cxx.

529{
530 if (this->GetNumberOfInputConnections(1) < 1)
531 {
532 return nullptr;
533 }
534 return vtkPolyData::SafeDownCast(this->GetExecutive()->GetInputData(1, 0));
535}

Referenced by RequestData().

◆ New()

static vtkTensorGlyphColor * vtkTensorGlyphColor::New ( )
static

Construct object with scaling on and scale factor 1.0. Eigenvalues are extracted, glyphs are colored with input scalar data, and logarithmic scaling is turned off.

◆ RequestData()

int vtkTensorGlyphColor::RequestData ( vtkInformation *  ,
vtkInformationVector **  ,
vtkInformationVector *   
)
overrideprotected

Definition at line 97 of file vtkTensorGlyphColor.cxx.

101{
102 // get the info objects
103 vtkInformation *inInfo = inputVector[0]->GetInformationObject(0);
104 vtkInformation *sourceInfo = inputVector[1]->GetInformationObject(0);
105 vtkInformation *outInfo = outputVector->GetInformationObject(0);
106
107 // get the input and output
108 vtkDataSet *input = vtkDataSet::SafeDownCast(
109 inInfo->Get(vtkDataObject::DATA_OBJECT()));
110 vtkPolyData *source = vtkPolyData::SafeDownCast(
111 sourceInfo->Get(vtkDataObject::DATA_OBJECT()));
112 vtkPolyData *output = vtkPolyData::SafeDownCast(
113 outInfo->Get(vtkDataObject::DATA_OBJECT()));
114
115 vtkDataArray *inTensors;
116 double tensor[9];
117 vtkDataArray *inScalars;
118 vtkIdType numPts, numSourcePts, numSourceCells, inPtId, i;
119 int j;
120 vtkPoints *sourcePts;
121 vtkDataArray *sourceNormals;
122 vtkCellArray *sourceCells, *cells;
123 vtkPoints *newPts;
124 vtkFloatArray *newScalars=nullptr;
125 vtkFloatArray *newNormals=nullptr;
126 double x[3], s;
127 vtkTransform *trans;
128 vtkCell *cell;
129 vtkIdList *cellPts;
130 int npts;
131 vtkIdType *pts;
132 vtkIdType ptIncr, cellId;
133 vtkIdType subIncr;
134 int numDirs, dir, eigen_dir, symmetric_dir;
135 vtkMatrix4x4 *matrix;
136 double *m[3], w[3], *v[3];
137 double m0[3], m1[3], m2[3];
138 double v0[3], v1[3], v2[3];
139 double xv[3], yv[3], zv[3];
140 double maxScale;
141
142 numDirs = (this->ThreeGlyphs?3:1)*(this->Symmetric+1);
143
144 // set up working matrices
145 m[0] = m0; m[1] = m1; m[2] = m2;
146 v[0] = v0; v[1] = v1; v[2] = v2;
147
148 vtkDebugMacro(<<"Generating tensor glyphs");
149
150 vtkPointData *outPD = output->GetPointData();
151 inTensors = this->GetInputArrayToProcess(0, inputVector);
152 inScalars = this->GetInputArrayToProcess(1, inputVector);
153 numPts = input->GetNumberOfPoints();
154
155 if ( !inTensors || numPts < 1 )
156 {
157 vtkErrorMacro(<<"No data to glyph!");
158 return 1;
159 }
160
161 pts = new vtkIdType[source->GetMaxCellSize()];
162 trans = vtkTransform::New();
163 matrix = vtkMatrix4x4::New();
164
165 //
166 // Allocate storage for output PolyData
167 //
168 sourcePts = source->GetPoints();
169 numSourcePts = sourcePts->GetNumberOfPoints();
170 numSourceCells = source->GetNumberOfCells();
171
172 newPts = vtkPoints::New();
173 newPts->Allocate(numDirs*numPts*numSourcePts);
174
175 // Setting up for calls to PolyData::InsertNextCell()
176 if ( (sourceCells=source->GetVerts())->GetNumberOfCells() > 0 )
177 {
178 cells = vtkCellArray::New();
179 cells->Allocate(numDirs*numPts*sourceCells->GetSize());
180 output->SetVerts(cells);
181 cells->Delete();
182 }
183 if ( (sourceCells=this->GetSource()->GetLines())->GetNumberOfCells() > 0 )
184 {
185 cells = vtkCellArray::New();
186 cells->Allocate(numDirs*numPts*sourceCells->GetSize());
187 output->SetLines(cells);
188 cells->Delete();
189 }
190 if ( (sourceCells=this->GetSource()->GetPolys())->GetNumberOfCells() > 0 )
191 {
192 cells = vtkCellArray::New();
193 cells->Allocate(numDirs*numPts*sourceCells->GetSize());
194 output->SetPolys(cells);
195 cells->Delete();
196 }
197 if ( (sourceCells=this->GetSource()->GetStrips())->GetNumberOfCells() > 0 )
198 {
199 cells = vtkCellArray::New();
200 cells->Allocate(numDirs*numPts*sourceCells->GetSize());
201 output->SetStrips(cells);
202 cells->Delete();
203 }
204
205 // only copy scalar data through
206 vtkPointData *pd = this->GetSource()->GetPointData();
207 // generate scalars if eigenvalues are chosen or if scalars exist.
208 if (this->ColorGlyphs &&
209 ((this->ColorMode == COLOR_BY_EIGENVALUES) ||
210 (inScalars && (this->ColorMode == COLOR_BY_SCALARS)) ) )
211 {
212 newScalars = vtkFloatArray::New();
213 newScalars->SetNumberOfComponents(4);
214 newScalars->Allocate(numDirs*numPts*numSourcePts);
215 if (this->ColorMode == COLOR_BY_EIGENVALUES)
216 {
217 newScalars->SetName("MaxEigenvalue");
218 }
219 else
220 {
221 newScalars->SetName(inScalars->GetName());
222 }
223 }
224 else
225 {
226 outPD->CopyAllOff();
227 outPD->CopyScalarsOn();
228 outPD->CopyAllocate(pd,numDirs*numPts*numSourcePts);
229 }
230 if ( (sourceNormals = pd->GetNormals()) )
231 {
232 newNormals = vtkFloatArray::New();
233 newNormals->SetNumberOfComponents(3);
234 newNormals->SetName("Normals");
235 newNormals->Allocate(numDirs*3*numPts*numSourcePts);
236 }
237 //
238 // First copy all topology (transformation independent)
239 //
240 for (inPtId=0; inPtId < numPts; inPtId++)
241 {
242 ptIncr = numDirs * inPtId * numSourcePts;
243 for (cellId=0; cellId < numSourceCells; cellId++)
244 {
245 cell = this->GetSource()->GetCell(cellId);
246 cellPts = cell->GetPointIds();
247 npts = (int)cellPts->GetNumberOfIds();
248 for (dir=0; dir < numDirs; dir++)
249 {
250 // This variable may be removed, but that
251 // will not improve readability
252 subIncr = ptIncr + dir*numSourcePts;
253 for (i=0; i < npts; i++)
254 {
255 pts[i] = cellPts->GetId(i) + subIncr;
256 }
257 output->InsertNextCell(cell->GetCellType(),npts,pts);
258 }
259 }
260 }
261 //
262 // Traverse all Input points, transforming glyph at Source points
263 //
264 trans->PreMultiply();
265
266 for (inPtId=0; inPtId < numPts; inPtId++)
267 {
268 ptIncr = numDirs * inPtId * numSourcePts;
269
270 // Translation is postponed
271 // Symmetric tensor support
272 inTensors->GetTuple(inPtId, tensor);
273 if (inTensors->GetNumberOfComponents() == 6)
274 {
275 vtkMath::TensorFromSymmetricTensor(tensor);
276 }
277
278 // compute orientation vectors and scale factors from tensor
279 if ( this->ExtractEigenvalues ) // extract appropriate eigenfunctions
280 {
281 // We are interested in the symmetrical part of the tensor only, since
282 // eigenvalues are real if and only if the matrice of reals is symmetrical
283 for (j=0; j<3; j++)
284 {
285 for (i=0; i<3; i++)
286 {
287 m[i][j] = 0.5 * (tensor[i + 3 * j] + tensor[j + 3 * i]);
288 }
289 }
290 vtkMath::Jacobi(m, w, v);
291
292 //copy eigenvectors
293 xv[0] = v[0][0]; xv[1] = v[1][0]; xv[2] = v[2][0];
294 yv[0] = v[0][1]; yv[1] = v[1][1]; yv[2] = v[2][1];
295 zv[0] = v[0][2]; zv[1] = v[1][2]; zv[2] = v[2][2];
296 }
297 else //use tensor columns as eigenvectors
298 {
299 for (i=0; i<3; i++)
300 {
301 xv[i] = tensor[i];
302 yv[i] = tensor[i+3];
303 zv[i] = tensor[i+6];
304 }
305 w[0] = vtkMath::Normalize(xv);
306 w[1] = vtkMath::Normalize(yv);
307 w[2] = vtkMath::Normalize(zv);
308 }
309
310 // compute scale factors
311 w[0] *= this->ScaleFactor;
312 w[1] *= this->ScaleFactor;
313 w[2] *= this->ScaleFactor;
314
315 if ( this->ClampScaling )
316 {
317 for (maxScale=0.0, i=0; i<3; i++)
318 {
319 if ( maxScale < fabs(w[i]) )
320 {
321 maxScale = fabs(w[i]);
322 }
323 }
324 if ( maxScale > this->MaxScaleFactor )
325 {
326 maxScale = this->MaxScaleFactor / maxScale;
327 for (i=0; i<3; i++)
328 {
329 w[i] *= maxScale; //preserve overall shape of glyph
330 }
331 }
332 }
333
334 // normalization is postponed
335
336 // make sure scale is okay (non-zero) and scale data
337 for (maxScale=0.0, i=0; i<3; i++)
338 {
339 if ( w[i] > maxScale )
340 {
341 maxScale = w[i];
342 }
343 }
344 if ( maxScale == 0.0 )
345 {
346 maxScale = 1.0;
347 }
348 for (i=0; i<3; i++)
349 {
350 if ( w[i] == 0.0 )
351 {
352 w[i] = maxScale * 1.0e-06;
353 }
354 }
355
356 // Now do the real work for each "direction"
357
358 for (dir=0; dir < numDirs; dir++)
359 {
360 eigen_dir = dir%(this->ThreeGlyphs?3:1);
361 symmetric_dir = dir/(this->ThreeGlyphs?3:1);
362
363 // Remove previous scales ...
364 trans->Identity();
365
366 // translate Source to Input point
367 input->GetPoint(inPtId, x);
368 trans->Translate(x[0], x[1], x[2]);
369
370 // normalized eigenvectors rotate object for eigen direction 0
371 matrix->Element[0][0] = xv[0];
372 matrix->Element[0][1] = yv[0];
373 matrix->Element[0][2] = zv[0];
374 matrix->Element[1][0] = xv[1];
375 matrix->Element[1][1] = yv[1];
376 matrix->Element[1][2] = zv[1];
377 matrix->Element[2][0] = xv[2];
378 matrix->Element[2][1] = yv[2];
379 matrix->Element[2][2] = zv[2];
380 trans->Concatenate(matrix);
381
382 if (eigen_dir == 1)
383 {
384 trans->RotateZ(90.0);
385 }
386
387 if (eigen_dir == 2)
388 {
389 trans->RotateY(-90.0);
390 }
391
392 if (this->ThreeGlyphs)
393 {
394 trans->Scale(w[eigen_dir], this->ScaleFactor, this->ScaleFactor);
395 }
396 else
397 {
398 trans->Scale(w[0], w[1], w[2]);
399 }
400
401 // Mirror second set to the symmetric position
402 if (symmetric_dir == 1)
403 {
404 trans->Scale(-1.,1.,1.);
405 }
406
407 // if the eigenvalue is negative, shift to reverse direction.
408 // The && is there to ensure that we do not change the
409 // old behaviour of vtkTensorGlyphColors (which only used one dir),
410 // in case there is an oriented glyph, e.g. an arrow.
411 if (w[eigen_dir] < 0 && numDirs > 1)
412 {
413 trans->Translate(-this->Length, 0., 0.);
414 }
415
416 // multiply points (and normals if available) by resulting
417 // matrix
418 trans->TransformPoints(sourcePts,newPts);
419
420 // Apply the transformation to a series of points,
421 // and append the results to outPts.
422 if ( newNormals )
423 {
424 // a negative determinant means the transform turns the
425 // glyph surface inside out, and its surface normals all
426 // point inward. The following scale corrects the surface
427 // normals to point outward.
428 if (trans->GetMatrix()->Determinant() < 0)
429 {
430 trans->Scale(-1.0,-1.0,-1.0);
431 }
432 trans->TransformNormals(sourceNormals,newNormals);
433 }
434
435 // Copy point data from source
436 if ( this->ColorGlyphs && inScalars &&
437 (this->ColorMode == COLOR_BY_SCALARS) )
438 {
439 for (i=0; i < numSourcePts; i++)
440 {
441 auto st = inScalars->GetTuple(inPtId);
442 newScalars->InsertTuple4(ptIncr+i,st[0],st[1],st[2],st[3]);
443 }
444 }
445 else if (this->ColorGlyphs &&
447 {
448 // If ThreeGlyphs is false we use the first (largest)
449 // eigenvalue as scalar.
450 s = w[eigen_dir];
451 for (i=0; i < numSourcePts; i++)
452 {
453 newScalars->InsertTuple(ptIncr+i, &s);
454 }
455 }
456 else
457 {
458 for (i=0; i < numSourcePts; i++)
459 {
460 outPD->CopyData(pd,i,ptIncr+i);
461 }
462 }
463 ptIncr += numSourcePts;
464 }
465 }
466 vtkDebugMacro(<<"Generated " << numPts <<" tensor glyphs");
467 //
468 // Update output and release memory
469 //
470 delete [] pts;
471
472 output->SetPoints(newPts);
473 newPts->Delete();
474
475 if ( newScalars )
476 {
477 int idx = outPD->AddArray(newScalars);
478 outPD->SetActiveAttribute(idx, vtkDataSetAttributes::SCALARS);
479 newScalars->Delete();
480 }
481
482 if ( newNormals )
483 {
484 outPD->SetNormals(newNormals);
485 newNormals->Delete();
486 }
487
488 output->Squeeze();
489 trans->Delete();
490 matrix->Delete();
491
492 return 1;
493}

◆ RequestUpdateExtent()

int vtkTensorGlyphColor::RequestUpdateExtent ( vtkInformation *  ,
vtkInformationVector **  ,
vtkInformationVector *   
)
overrideprotected

Definition at line 65 of file vtkTensorGlyphColor.cxx.

69{
70 // get the info objects
71 vtkInformation *inInfo = inputVector[0]->GetInformationObject(0);
72 vtkInformation *sourceInfo = inputVector[1]->GetInformationObject(0);
73 vtkInformation *outInfo = outputVector->GetInformationObject(0);
74
75 if (sourceInfo)
76 {
77 sourceInfo->Set(vtkStreamingDemandDrivenPipeline::UPDATE_PIECE_NUMBER(),
78 0);
79 sourceInfo->Set(vtkStreamingDemandDrivenPipeline::UPDATE_NUMBER_OF_PIECES(),
80 1);
81 sourceInfo->Set(vtkStreamingDemandDrivenPipeline::UPDATE_NUMBER_OF_GHOST_LEVELS(),
82 0);
83 }
84
85 inInfo->Set(vtkStreamingDemandDrivenPipeline::UPDATE_PIECE_NUMBER(),
86 outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_PIECE_NUMBER()));
87 inInfo->Set(vtkStreamingDemandDrivenPipeline::UPDATE_NUMBER_OF_PIECES(),
88 outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_NUMBER_OF_PIECES()));
89 inInfo->Set(vtkStreamingDemandDrivenPipeline::UPDATE_NUMBER_OF_GHOST_LEVELS(),
90 outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_NUMBER_OF_GHOST_LEVELS()));
91 inInfo->Set(vtkStreamingDemandDrivenPipeline::EXACT_EXTENT(), 1);
92
93 return 1;
94}

◆ SetColorModeToEigenvalues()

void vtkTensorGlyphColor::SetColorModeToEigenvalues ( )
inline

Definition at line 212 of file vtkTensorGlyphColor.h.

213 {this->SetColorMode(COLOR_BY_EIGENVALUES);};

◆ SetColorModeToScalars()

void vtkTensorGlyphColor::SetColorModeToScalars ( )
inline

Definition at line 210 of file vtkTensorGlyphColor.h.

211 {this->SetColorMode(COLOR_BY_SCALARS);};

◆ SetSourceConnection() [1/2]

void vtkTensorGlyphColor::SetSourceConnection ( int  id,
vtkAlgorithmOutput *  algOutput 
)

Specify a source object at a specified table location. New style. Source connection is stored in port 1. This method is equivalent to SetInputConnection(1, id, outputPort).

Definition at line 496 of file vtkTensorGlyphColor.cxx.

497{
498 if (id < 0)
499 {
500 vtkErrorMacro("Bad index " << id << " for source.");
501 return;
502 }
503
504 int numConnections = this->GetNumberOfInputConnections(1);
505 if (id < numConnections)
506 {
507 this->SetNthInputConnection(1, id, algOutput);
508 }
509 else if (id == numConnections && algOutput)
510 {
511 this->AddInputConnection(1, algOutput);
512 }
513 else if (algOutput)
514 {
515 vtkWarningMacro("The source id provided is larger than the maximum "
516 "source id, using " << numConnections << " instead.");
517 this->AddInputConnection(1, algOutput);
518 }
519}

◆ SetSourceConnection() [2/2]

void vtkTensorGlyphColor::SetSourceConnection ( vtkAlgorithmOutput *  algOutput)
inline

Definition at line 120 of file vtkTensorGlyphColor.h.

121 {
122 this->SetSourceConnection(0, algOutput);
123 }
void SetSourceConnection(int id, vtkAlgorithmOutput *algOutput)

◆ SetSourceData()

void vtkTensorGlyphColor::SetSourceData ( vtkPolyData *  source)

Specify the geometry to copy to each point. Note that this method does not connect the pipeline. The algorithm will work on the input data as it is without updating the producer of the data. See SetSourceConnection for connecting the pipeline.

Definition at line 522 of file vtkTensorGlyphColor.cxx.

523{
524 this->SetInputData(1, source);
525}

◆ vtkBooleanMacro() [1/6]

vtkTensorGlyphColor::vtkBooleanMacro ( ClampScaling  ,
vtkTypeBool   
)

◆ vtkBooleanMacro() [2/6]

vtkTensorGlyphColor::vtkBooleanMacro ( ColorGlyphs  ,
vtkTypeBool   
)

◆ vtkBooleanMacro() [3/6]

vtkTensorGlyphColor::vtkBooleanMacro ( ExtractEigenvalues  ,
vtkTypeBool   
)

◆ vtkBooleanMacro() [4/6]

vtkTensorGlyphColor::vtkBooleanMacro ( Scaling  ,
vtkTypeBool   
)

◆ vtkBooleanMacro() [5/6]

vtkTensorGlyphColor::vtkBooleanMacro ( Symmetric  ,
vtkTypeBool   
)

◆ vtkBooleanMacro() [6/6]

vtkTensorGlyphColor::vtkBooleanMacro ( ThreeGlyphs  ,
vtkTypeBool   
)

◆ vtkGetMacro() [1/10]

vtkTensorGlyphColor::vtkGetMacro ( ClampScaling  ,
vtkTypeBool   
)

◆ vtkGetMacro() [2/10]

vtkTensorGlyphColor::vtkGetMacro ( ColorGlyphs  ,
vtkTypeBool   
)

◆ vtkGetMacro() [3/10]

vtkTensorGlyphColor::vtkGetMacro ( ColorMode  ,
int   
)

◆ vtkGetMacro() [4/10]

vtkTensorGlyphColor::vtkGetMacro ( ExtractEigenvalues  ,
vtkTypeBool   
)

◆ vtkGetMacro() [5/10]

vtkTensorGlyphColor::vtkGetMacro ( Length  ,
double   
)

◆ vtkGetMacro() [6/10]

vtkTensorGlyphColor::vtkGetMacro ( MaxScaleFactor  ,
double   
)

◆ vtkGetMacro() [7/10]

vtkTensorGlyphColor::vtkGetMacro ( ScaleFactor  ,
double   
)

◆ vtkGetMacro() [8/10]

vtkTensorGlyphColor::vtkGetMacro ( Scaling  ,
vtkTypeBool   
)

◆ vtkGetMacro() [9/10]

vtkTensorGlyphColor::vtkGetMacro ( Symmetric  ,
vtkTypeBool   
)

◆ vtkGetMacro() [10/10]

vtkTensorGlyphColor::vtkGetMacro ( ThreeGlyphs  ,
vtkTypeBool   
)

◆ vtkSetClampMacro()

vtkTensorGlyphColor::vtkSetClampMacro ( ColorMode  ,
int  ,
COLOR_BY_SCALARS  ,
COLOR_BY_EIGENVALUES   
)

Set the color mode to be used for the glyphs. This can be set to use the input scalars (default) or to use the eigenvalues at the point. If ThreeGlyphs is set and the eigenvalues are chosen for coloring then each glyph is colored by the corresponding eigenvalue and if not set the color corresponding to the largest eigenvalue is chosen. The recognized values are: COLOR_BY_SCALARS = 0 (default) COLOR_BY_EIGENVALUES = 1

◆ vtkSetMacro() [1/9]

vtkTensorGlyphColor::vtkSetMacro ( ClampScaling  ,
vtkTypeBool   
)

Turn on/off scalar clamping. If scalar clamping is on, the ivar MaxScaleFactor is used to control the maximum scale factor. (This is useful to prevent uncontrolled scaling near singularities.)

◆ vtkSetMacro() [2/9]

vtkTensorGlyphColor::vtkSetMacro ( ColorGlyphs  ,
vtkTypeBool   
)

Turn on/off coloring of glyph with input scalar data or eigenvalues. If false, or input scalar data not present, then the scalars from the source object are passed through the filter.

◆ vtkSetMacro() [3/9]

vtkTensorGlyphColor::vtkSetMacro ( ExtractEigenvalues  ,
vtkTypeBool   
)

Turn on/off extraction of eigenvalues from tensor.

◆ vtkSetMacro() [4/9]

vtkTensorGlyphColor::vtkSetMacro ( Length  ,
double   
)

Set/Get the distance, along x, from the origin to the end of the source glyph. It is used to draw the symmetric glyphs.

◆ vtkSetMacro() [5/9]

vtkTensorGlyphColor::vtkSetMacro ( MaxScaleFactor  ,
double   
)

Set/Get the maximum allowable scale factor. This value is compared to the combination of the scale factor times the eigenvalue. If less, the scale factor is reset to the MaxScaleFactor. The boolean ClampScaling has to be "on" for this to work.

◆ vtkSetMacro() [6/9]

vtkTensorGlyphColor::vtkSetMacro ( ScaleFactor  ,
double   
)

Specify scale factor to scale object by. (Scale factor always affects output even if scaling is off.)

◆ vtkSetMacro() [7/9]

vtkTensorGlyphColor::vtkSetMacro ( Scaling  ,
vtkTypeBool   
)

Turn on/off scaling of glyph with eigenvalues.

◆ vtkSetMacro() [8/9]

vtkTensorGlyphColor::vtkSetMacro ( Symmetric  ,
vtkTypeBool   
)

Turn on/off drawing a mirror of each glyph

◆ vtkSetMacro() [9/9]

vtkTensorGlyphColor::vtkSetMacro ( ThreeGlyphs  ,
vtkTypeBool   
)

Turn on/off drawing three glyphs

◆ vtkTypeMacro()

vtkTensorGlyphColor::vtkTypeMacro ( vtkTensorGlyphColor  ,
vtkPolyDataAlgorithm   
) &

Member Data Documentation

◆ ClampScaling

vtkTypeBool vtkTensorGlyphColor::ClampScaling
protected

Definition at line 251 of file vtkTensorGlyphColor.h.

Referenced by RequestData(), and vtkStandardNewMacro().

◆ ColorGlyphs

vtkTypeBool vtkTensorGlyphColor::ColorGlyphs
protected

Definition at line 249 of file vtkTensorGlyphColor.h.

Referenced by RequestData(), and vtkStandardNewMacro().

◆ ColorMode

int vtkTensorGlyphColor::ColorMode
protected

Definition at line 250 of file vtkTensorGlyphColor.h.

Referenced by RequestData(), and vtkStandardNewMacro().

◆ ExtractEigenvalues

vtkTypeBool vtkTensorGlyphColor::ExtractEigenvalues
protected

Definition at line 248 of file vtkTensorGlyphColor.h.

Referenced by RequestData(), and vtkStandardNewMacro().

◆ Length

double vtkTensorGlyphColor::Length
protected

Definition at line 255 of file vtkTensorGlyphColor.h.

Referenced by RequestData(), and vtkStandardNewMacro().

◆ MaxScaleFactor

double vtkTensorGlyphColor::MaxScaleFactor
protected

Definition at line 252 of file vtkTensorGlyphColor.h.

Referenced by RequestData(), and vtkStandardNewMacro().

◆ override

vtkIndent indent vtkTensorGlyphColor::override

Definition at line 93 of file vtkTensorGlyphColor.h.

◆ ScaleFactor

double vtkTensorGlyphColor::ScaleFactor
protected

Definition at line 247 of file vtkTensorGlyphColor.h.

Referenced by RequestData(), and vtkStandardNewMacro().

◆ Scaling

vtkTypeBool vtkTensorGlyphColor::Scaling
protected

Definition at line 246 of file vtkTensorGlyphColor.h.

Referenced by vtkStandardNewMacro().

◆ Symmetric

vtkTypeBool vtkTensorGlyphColor::Symmetric
protected

Definition at line 254 of file vtkTensorGlyphColor.h.

Referenced by RequestData(), and vtkStandardNewMacro().

◆ ThreeGlyphs

vtkTypeBool vtkTensorGlyphColor::ThreeGlyphs
protected

Definition at line 253 of file vtkTensorGlyphColor.h.

Referenced by RequestData(), and vtkStandardNewMacro().


The documentation for this class was generated from the following files: