TTK
Loading...
Searching...
No Matches
ttkGridLayout.cpp
Go to the documentation of this file.
1#include <ttkGridLayout.h>
2
3#include <vtkObjectFactory.h>
4
5#include <vtkInformation.h>
6#include <vtkInformationVector.h>
7
8#include <vtkSmartPointer.h>
9
10#include <vtkDataSet.h>
11#include <vtkImageData.h>
12#include <vtkMultiBlockDataSet.h>
13#include <vtkPointSet.h>
14
16
18 this->setDebugMsgPrefix("GridLayout");
19
20 this->SetNumberOfInputPorts(1);
21 this->SetNumberOfOutputPorts(1);
22}
23
25
26int ttkGridLayout::FillInputPortInformation(int port, vtkInformation *info) {
27 if(port == 0)
28 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkMultiBlockDataSet");
29 else
30 return 0;
31 return 1;
32}
33
34int ttkGridLayout::FillOutputPortInformation(int port, vtkInformation *info) {
35 if(port == 0)
37 else
38 return 0;
39 return 1;
40}
41
42int ttkGridLayout::CopyObject(vtkDataObject *output, vtkDataObject *input) {
43 if(input->IsA("vtkImageData")) {
44 output->ShallowCopy(input);
45 } else if(input->IsA("vtkPointSet")) {
46 auto outputAsPS = (vtkPointSet *)output;
47 outputAsPS->ShallowCopy(input);
48
49 auto points = vtkSmartPointer<vtkPoints>::New();
50 points->DeepCopy(outputAsPS->GetPoints());
51
52 outputAsPS->SetPoints(points);
53 } else {
54 output->DeepCopy(input);
55 }
56
57 return 1;
58}
59
60int ttkGridLayout::TranslateObject(vtkDataObject *input,
61 const size_t &colAxis,
62 const size_t &rowAxis,
63 const double &dw,
64 const double &dh) {
65 auto inputMB = vtkMultiBlockDataSet::SafeDownCast(input);
66 auto inputAsPS = vtkPointSet::SafeDownCast(input);
67 auto inputAsID = vtkImageData::SafeDownCast(input);
68
69 if(inputMB) {
70 size_t const n = inputMB->GetNumberOfBlocks();
71 for(size_t i = 0; i < n; i++)
72 if(!this->TranslateObject(inputMB->GetBlock(i), colAxis, rowAxis, dw, dh))
73 return 0;
74
75 } else if(inputAsPS) {
76 size_t const nCoords = inputAsPS->GetNumberOfPoints() * 3;
77 auto points = inputAsPS->GetPoints();
78 auto pointCoords = (float *)points->GetVoidPointer(0);
79 for(size_t i = 0; i < nCoords; i += 3) {
80 pointCoords[i + colAxis] += dw;
81 pointCoords[i + rowAxis] += dh;
82 }
83
84 } else if(inputAsID) {
85 double origin[3] = {0, 0, 0};
86 origin[colAxis] = dw;
87 origin[rowAxis] = dh;
88 inputAsID->SetOrigin(origin);
89 } else {
90 return 0;
91 }
92
93 return 1;
94}
95
96int ttkGridLayout::RequestData(vtkInformation *ttkNotUsed(request),
97 vtkInformationVector **inputVector,
98 vtkInformationVector *outputVector) {
99 ttk::Timer globalTimer;
100
101 // Get Input and Output
102 auto inputMB = vtkMultiBlockDataSet::GetData(inputVector[0]);
103 auto outputMB = vtkMultiBlockDataSet::GetData(outputVector);
104
105 // Determine Grid Axis
106 int const colAxis = this->GetColAxis();
107 int const rowAxis = this->GetRowAxis();
108
109 // Compute width and height of grid cells
110 double bounds[6];
111 double maxWidth = 0;
112 double maxHeight = 0;
113
114 size_t const nBlocks = inputMB->GetNumberOfBlocks();
115
116 this->printMsg("Translating " + std::to_string(nBlocks) + " object(s)", 0,
118
119 for(size_t i = 0; i < nBlocks; i++) {
120 auto block = inputMB->GetBlock(i);
121 if(block->IsA("vtkMultiBlockDataSet")) {
122 auto blockAsMB = vtkMultiBlockDataSet::SafeDownCast(block);
123 blockAsMB->GetBounds(bounds);
124 } else if(block->IsA("vtkDataSet")) {
125 auto blockAsDS = vtkDataSet::SafeDownCast(block);
126 blockAsDS->GetBounds(bounds);
127 } else {
128 this->printErr("Unable to determine bounding box of block #"
129 + std::to_string(i) + " with type "
130 + std::string(block->GetClassName()) + ".");
131 return 0;
132 }
133
134 double const blockWidth = bounds[colAxis * 2 + 1] - bounds[colAxis * 2];
135 double const blockHeight = bounds[rowAxis * 2 + 1] - bounds[rowAxis * 2];
136 if(maxWidth < blockWidth)
137 maxWidth = blockWidth;
138 if(maxHeight < blockHeight)
139 maxHeight = blockHeight;
140 }
141
142 // apply gap
143 maxWidth += maxWidth * this->GetColGap() / 100.;
144 maxHeight += maxHeight * this->GetRowGap() / 100.;
145
146 // Determine Grid Structure
147 const size_t nRows
148 = this->GetNumberOfRows() < 1 ? 0 : (size_t)this->GetNumberOfRows();
149 const size_t nColumns
150 = nRows == 0 ? std::ceil(std::sqrt(nBlocks)) : std::ceil(nBlocks / nRows);
151
152 for(size_t i = 0; i < nBlocks; i++) {
153 // get block
154 auto block = inputMB->GetBlock(i);
155
156 auto outBlock = vtkSmartPointer<vtkDataObject>::Take(block->NewInstance());
157 this->CopyObject(outBlock, block);
158
159 const size_t row = std::floor(i / nColumns);
160 const size_t col = i % nColumns;
161
162 if(!this->TranslateObject(
163 outBlock, colAxis, rowAxis, col * maxWidth, row * maxHeight)) {
164 this->printErr("Unable to translate block #" + std::to_string(i)
165 + " of type '" + std::string(outBlock->GetClassName())
166 + "'.");
167 return 0;
168 }
169
170 outputMB->SetBlock(i, outBlock);
171 }
172
173 this->printMsg("Translating " + std::to_string(nBlocks) + " object(s)", 1,
174 globalTimer.getElapsedTime());
175
176 return 1;
177}
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
static vtkInformationIntegerKey * SAME_DATA_TYPE_AS_INPUT_PORT()
TTK VTK-filter that arranges vtkDataSets on a grid.
virtual int GetColAxis()
virtual int GetRowAxis()
int CopyObject(vtkDataObject *output, vtkDataObject *input)
~ttkGridLayout() override
int FillInputPortInformation(int port, vtkInformation *info) override
int TranslateObject(vtkDataObject *input, const size_t &colAxis, const size_t &rowAxis, const double &dw, const double &dh)
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillOutputPortInformation(int port, vtkInformation *info) override
virtual double GetRowGap()
virtual double GetColGap()
virtual int GetNumberOfRows()
void setDebugMsgPrefix(const std::string &prefix)
Definition Debug.h:364
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
double getElapsedTime()
Definition Timer.h:15
vtkStandardNewMacro(ttkGridLayout)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)