3#ifdef VIENNALS_VTK_RENDERING
9#include <vtkAutoInit.h>
11#include <vtkCellArray.h>
12#include <vtkCellData.h>
13#include <vtkDataSetMapper.h>
14#include <vtkInteractorStyleImage.h>
15#include <vtkLookupTable.h>
17#include <vtkPolyData.h>
18#include <vtkPolyDataMapper.h>
19#include <vtkProperty.h>
20#include <vtkRenderWindow.h>
21#include <vtkRenderWindowInteractor.h>
22#include <vtkRenderer.h>
23#include <vtkUnstructuredGrid.h>
25#ifndef VIENNALS_VTK_MODULE_INIT
26VTK_MODULE_INIT(vtkRenderingOpenGL2);
27VTK_MODULE_INIT(vtkInteractionStyle);
28VTK_MODULE_INIT(vtkRenderingFreeType);
29VTK_MODULE_INIT(vtkRenderingUI);
32class ImagePanInteractorStyle :
public vtkInteractorStyleImage {
34 static ImagePanInteractorStyle *New();
35 vtkTypeMacro(ImagePanInteractorStyle, vtkInteractorStyleImage);
37 void OnLeftButtonDown()
override { this->StartPan(); }
39 void OnLeftButtonUp()
override { this->EndPan(); }
42vtkStandardNewMacro(ImagePanInteractorStyle);
46template <
typename T>
class VTKRenderWindow {
48 VTKRenderWindow() { initialize(); }
49 VTKRenderWindow(SmartPointer<Mesh<T>> passedMesh) {
56 interactor->SetRenderWindow(
nullptr);
59 renderWindow->RemoveRenderer(renderer);
63 void setMesh(SmartPointer<Mesh<T>> passedMesh) {
65 auto matIds =
mesh->getCellData().getScalarData(
"MaterialIds",
false);
67 materialIds = *matIds;
72 void setVolumeMesh(vtkSmartPointer<vtkUnstructuredGrid> volumeVTK) {
73 volumeMesh = volumeVTK;
77 void setMaterialIds(
const std::vector<T> &ids) { materialIds = ids; }
79 auto setBackgroundColor(
const std::array<double, 3> &color) {
80 backgroundColor = color;
82 renderer->SetBackground(backgroundColor.data());
89 if (mesh ==
nullptr && volumeMesh ==
nullptr) {
90 VIENNACORE_LOG_WARNING(
"No mesh set for rendering.");
96 if (
auto style = interactor->GetInteractorStyle()) {
97 style->SetDefaultRenderer(renderer);
98 style->SetCurrentRenderer(renderer);
100 interactor->SetRenderWindow(renderWindow);
101 renderWindow->AddRenderer(renderer);
103 renderWindow->Render();
104 interactor->Initialize();
108 vtkSmartPointer<vtkRenderWindow> getRenderWindow() {
return renderWindow; }
110 auto enable2DMode() {
111 assert(renderer &&
"Renderer not initialized");
112 vtkCamera *cam = renderer->GetActiveCamera();
113 cam->ParallelProjectionOn();
115 auto style = vtkSmartPointer<ImagePanInteractorStyle>::New();
117 style->SetDefaultRenderer(renderer);
118 style->SetCurrentRenderer(renderer);
119 interactor->SetInteractorStyle(style);
125 renderer = vtkSmartPointer<vtkRenderer>::New();
126 renderer->SetBackground(backgroundColor.data());
128 renderWindow = vtkSmartPointer<vtkRenderWindow>::New();
129 renderWindow->SetWindowName(
"ViennaLS Render Window");
130 renderWindow->SetSize(windowSize.data());
133 interactor = vtkSmartPointer<vtkRenderWindowInteractor>::New();
134 interactor->SetRenderWindow(renderWindow);
137 auto style = vtkSmartPointer<vtkInteractorStyleTrackballCamera>::New();
138 interactor->SetInteractorStyle(style);
140 renderWindow->SetInteractor(interactor);
143 void updatePolyData() {
144 if (mesh ==
nullptr) {
148 polyData = vtkSmartPointer<vtkPolyData>::New();
151 vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
152 for (
const auto &node :
mesh->getNodes()) {
153 points->InsertNextPoint(node[0], node[1], node[2]);
155 polyData->SetPoints(points);
159 polyData->GetBounds(bounds);
162 if (!
mesh->vertices.empty()) {
163 vtkSmartPointer<vtkCellArray> verts =
164 vtkSmartPointer<vtkCellArray>::New();
165 for (
const auto &vertex :
mesh->vertices) {
166 verts->InsertNextCell(1);
167 verts->InsertCellPoint(vertex[0]);
169 polyData->SetVerts(verts);
173 if (!
mesh->lines.empty()) {
174 vtkSmartPointer<vtkCellArray> lines =
175 vtkSmartPointer<vtkCellArray>::New();
176 for (
const auto &line :
mesh->lines) {
177 lines->InsertNextCell(2);
178 lines->InsertCellPoint(line[0]);
179 lines->InsertCellPoint(line[1]);
181 polyData->SetLines(lines);
187 if (!
mesh->triangles.empty()) {
188 vtkSmartPointer<vtkCellArray> polys =
189 vtkSmartPointer<vtkCellArray>::New();
190 for (
const auto &triangle :
mesh->triangles) {
191 polys->InsertNextCell(3);
192 polys->InsertCellPoint(triangle[0]);
193 polys->InsertCellPoint(triangle[1]);
194 polys->InsertCellPoint(triangle[2]);
196 polyData->SetPolys(polys);
200 bool useMaterialIds =
201 !materialIds.empty() &&
202 (materialIds.size() ==
mesh->lines.size() +
mesh->triangles.size());
203 int minId = std::numeric_limits<int>::max();
204 int maxId = std::numeric_limits<int>::min();
205 if (useMaterialIds) {
206 vtkSmartPointer<vtkIntArray> matIdArray =
207 vtkSmartPointer<vtkIntArray>::New();
208 matIdArray->SetName(
"MaterialIds");
209 for (
const auto &
id : materialIds) {
210 int mId =
static_cast<int>(id);
211 matIdArray->InsertNextValue(mId);
212 minId = std::min(minId, mId);
213 maxId = std::max(maxId, mId);
215 polyData->GetCellData()->AddArray(matIdArray);
216 polyData->GetCellData()->SetActiveScalars(
"MaterialIds");
217 VIENNACORE_LOG_DEBUG(
"Added MaterialIds array to cell data.");
220 auto mapper = vtkSmartPointer<vtkPolyDataMapper>::New();
221 mapper->SetInputData(polyData);
223 if (useMaterialIds) {
224 mapper->SetScalarModeToUseCellData();
225 mapper->ScalarVisibilityOn();
226 mapper->SelectColorArray(
"MaterialIds");
228 vtkSmartPointer<vtkLookupTable> lut =
229 vtkSmartPointer<vtkLookupTable>::New();
231 lut->SetNumberOfTableValues(256);
232 lut->SetHueRange(0.667, 0.0);
233 lut->SetSaturationRange(1.0, 1.0);
234 lut->SetValueRange(1.0, 1.0);
237 mapper->SetLookupTable(lut);
238 mapper->SetScalarRange(minId, maxId);
241 auto actor = vtkSmartPointer<vtkActor>::New();
242 actor->SetMapper(mapper);
243 actor->GetProperty()->SetLineWidth(3.0);
245 renderer->AddActor(actor);
246 renderer->ResetCamera();
249 void updateVolumeMesh() {
250 auto mapper = vtkSmartPointer<vtkDataSetMapper>::New();
251 mapper->SetInputData(volumeMesh);
253 auto actor = vtkSmartPointer<vtkActor>::New();
254 actor->SetMapper(mapper);
256 renderer->AddActor(actor);
257 renderer->ResetCamera();
261 SmartPointer<Mesh<T>>
mesh =
nullptr;
262 std::vector<T> materialIds;
264 vtkSmartPointer<vtkUnstructuredGrid> volumeMesh =
nullptr;
266 vtkSmartPointer<vtkRenderer> renderer;
267 vtkSmartPointer<vtkRenderWindow> renderWindow;
268 vtkSmartPointer<vtkRenderWindowInteractor> interactor;
269 vtkSmartPointer<vtkPolyData> polyData;
271 std::array<double, 3> backgroundColor = {84.0 / 255, 89.0 / 255, 109.0 / 255};
272 std::array<int, 2> windowSize = {800, 600};
tuple bounds
Definition AirGapDeposition.py:23
mesh
Definition AirGapDeposition.py:36
Definition lsAdvect.hpp:37