VTK  9.3.0
vtkIntegrateAttributes.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
19#ifndef vtkIntegrateAttributes_h
20#define vtkIntegrateAttributes_h
21
22#include "vtkFiltersParallelModule.h" // For export macro
24
25VTK_ABI_NAMESPACE_BEGIN
26class vtkDataSet;
27class vtkIdList;
28class vtkInformation;
32
33class VTKFILTERSPARALLEL_EXPORT vtkIntegrateAttributes : public vtkUnstructuredGridAlgorithm
34{
35public:
38 void PrintSelf(ostream& os, vtkIndent indent) override;
39
41
46 vtkGetObjectMacro(Controller, vtkMultiProcessController);
48
50
54 vtkSetMacro(DivideAllCellDataByVolume, bool);
55 vtkGetMacro(DivideAllCellDataByVolume, bool);
57
58protected:
61
63
64 int RequestData(vtkInformation* request, vtkInformationVector** inputVector,
65 vtkInformationVector* outputVector) override;
66
67 // Create a default executive.
69
71
74
75 // The length, area or volume of the data set. Computed by Execute;
76 double Sum;
77 // ToCompute the location of the output point.
78 double SumCenter[3];
79
81
83 vtkDataSet* input, vtkUnstructuredGrid* output, vtkIdType cellId, vtkIdList* cellPtIds);
85 vtkDataSet* input, vtkUnstructuredGrid* output, vtkIdType cellId, vtkIdList* cellPtIds);
87 vtkDataSet* input, vtkUnstructuredGrid* output, vtkIdType cellId, vtkIdList* cellPtIds);
89 vtkIdType pt1Id, vtkIdType pt2Id, vtkIdType pt3Id);
91 vtkIdType pt1Id, vtkIdType pt2Id, vtkIdType pt3Id, vtkIdType pt4Id);
93 vtkDataSet* input, vtkUnstructuredGrid* output, vtkIdType cellId, vtkIdList* cellPtIds);
95 vtkDataSet* input, vtkUnstructuredGrid* output, vtkIdType cellId, vtkIdList* cellPtIds);
97 vtkDataSet* input, vtkUnstructuredGrid* output, vtkIdType cellId, vtkIdList* cellPtIds);
99 vtkDataSet* input, vtkUnstructuredGrid* output, vtkIdType cellId, vtkIdList* cellPtIds);
101 vtkDataSet* input, vtkUnstructuredGrid* output, vtkIdType cellId, vtkIdList* cellPtIds);
106 void ReceivePiece(vtkUnstructuredGrid* mergeTo, int fromId);
107
108 // This function assumes the data is in the format of the output of this filter with one
109 // point/cell having the value computed as its only tuple. It divides each value by sum,
110 // skipping the last data array if requested (so the volume doesn't get divided by itself
111 // and set to 1).
113 vtkDataSetAttributes* data, bool skipLastArray, double sum);
114
115private:
117 void operator=(const vtkIntegrateAttributes&) = delete;
118
119 class vtkFieldList;
120 vtkFieldList* CellFieldList;
121 vtkFieldList* PointFieldList;
122 int FieldListIndex;
123
124 void AllocateAttributes(vtkFieldList& fieldList, vtkDataSetAttributes* outda);
125 void ExecuteBlock(vtkDataSet* input, vtkUnstructuredGrid* output, int fieldset_index,
126 vtkFieldList& pdList, vtkFieldList& cdList);
127
128 void IntegrateData1(vtkDataSetAttributes* inda, vtkDataSetAttributes* outda, vtkIdType pt1Id,
129 double k, vtkFieldList& fieldlist, int fieldlist_index);
130 void IntegrateData2(vtkDataSetAttributes* inda, vtkDataSetAttributes* outda, vtkIdType pt1Id,
131 vtkIdType pt2Id, double k, vtkFieldList& fieldlist, int fieldlist_index);
132 void IntegrateData3(vtkDataSetAttributes* inda, vtkDataSetAttributes* outda, vtkIdType pt1Id,
133 vtkIdType pt2Id, vtkIdType pt3Id, double k, vtkFieldList& fieldlist, int fieldlist_index);
134 void IntegrateData4(vtkDataSetAttributes* inda, vtkDataSetAttributes* outda, vtkIdType pt1Id,
135 vtkIdType pt2Id, vtkIdType pt3Id, vtkIdType pt4Id, double k, vtkFieldList& fieldlist,
136 int fieldlist_index);
137
138public:
140 {
141 IntegrateAttrInfo = 2000,
142 IntegrateAttrData
143 };
144};
145
146VTK_ABI_NAMESPACE_END
147#endif
represent and manipulate attribute data in a dataset
abstract class to specify dataset behavior
Definition vtkDataSet.h:53
Superclass for all pipeline executives in VTK.
list of point or cell ids
Definition vtkIdList.h:23
a simple class to control print indentation
Definition vtkIndent.h:29
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Integrates lines, surfaces and volume.
~vtkIntegrateAttributes() override
void ZeroAttributes(vtkDataSetAttributes *outda)
void SendPiece(vtkUnstructuredGrid *src)
void ReceivePiece(vtkUnstructuredGrid *mergeTo, int fromId)
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
This is called by the superclass.
void IntegrateVoxel(vtkDataSet *input, vtkUnstructuredGrid *output, vtkIdType cellId, vtkIdList *cellPtIds)
void IntegrateGeneral3DCell(vtkDataSet *input, vtkUnstructuredGrid *output, vtkIdType cellId, vtkIdList *cellPtIds)
void IntegratePixel(vtkDataSet *input, vtkUnstructuredGrid *output, vtkIdType cellId, vtkIdList *cellPtIds)
int PieceNodeMinToNode0(vtkUnstructuredGrid *data)
void IntegrateGeneral1DCell(vtkDataSet *input, vtkUnstructuredGrid *output, vtkIdType cellId, vtkIdList *cellPtIds)
void IntegrateTriangleStrip(vtkDataSet *input, vtkUnstructuredGrid *output, vtkIdType cellId, vtkIdList *cellPtIds)
vtkMultiProcessController * Controller
static void DivideDataArraysByConstant(vtkDataSetAttributes *data, bool skipLastArray, double sum)
void IntegratePolygon(vtkDataSet *input, vtkUnstructuredGrid *output, vtkIdType cellId, vtkIdList *cellPtIds)
void IntegrateTriangle(vtkDataSet *input, vtkUnstructuredGrid *output, vtkIdType cellId, vtkIdType pt1Id, vtkIdType pt2Id, vtkIdType pt3Id)
int FillInputPortInformation(int, vtkInformation *) override
Fill the input port information objects for this algorithm.
int CompareIntegrationDimension(vtkDataSet *output, int dim)
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void IntegrateGeneral2DCell(vtkDataSet *input, vtkUnstructuredGrid *output, vtkIdType cellId, vtkIdList *cellPtIds)
void IntegratePolyLine(vtkDataSet *input, vtkUnstructuredGrid *output, vtkIdType cellId, vtkIdList *cellPtIds)
static vtkIntegrateAttributes * New()
void SetController(vtkMultiProcessController *controller)
Get/Set the parallel controller to use.
vtkExecutive * CreateDefaultExecutive() override
Create a default executive.
void IntegrateSatelliteData(vtkDataSetAttributes *inda, vtkDataSetAttributes *outda)
void IntegrateTetrahedron(vtkDataSet *input, vtkUnstructuredGrid *output, vtkIdType cellId, vtkIdType pt1Id, vtkIdType pt2Id, vtkIdType pt3Id, vtkIdType pt4Id)
Multiprocessing communication superclass.
Superclass for algorithms that produce only unstructured grid as output.
dataset represents arbitrary combinations of all possible cell types
int vtkIdType
Definition vtkType.h:315