oILAB
Loading...
Searching...
No Matches
BicrystalActor.cpp
Go to the documentation of this file.
1/* This file is part of MODEL, the Mechanics Of Defect Evolution Library.
2 *
3 * Copyright (C) 2011 by Giacomo Po <gpo@ucla.edu>.
4 *
5 * model is distributed without any warranty under the
6 * GNU General Public License (GPL) v2 <http://www.gnu.org/licenses/>.
7 */
8
9#ifndef model_BicrystalActor_cpp_
10#define model_BicrystalActor_cpp_
11
12#include <iostream>
13#include <deque>
14#include <string>
15
16
17#include <vtkVersion.h>
18#include <vtkSmartPointer.h>
19#include <vtkPolyData.h>
20#include <vtkPolyDataMapper.h>
21#include <vtkActor.h>
22#include <vtkMath.h>
23#include <vtkProperty.h>
24#include <vtkTubeFilter.h>
25#include <vtkPolyLine.h>
26#include <vtkSphereSource.h>
27#include <vtkArrowSource.h>
28#include <vtkGlyph3D.h>
29#include <vtkDoubleArray.h>
30#include <vtkPointData.h>
31#include <vtkLabeledDataMapper.h>
32#include <vtkFloatArray.h>
33#include <vtkTextProperty.h>
34#include <vtkActor2D.h>
35#include <vtkProperty2D.h>
36#include <vtkRenderer.h>
37
38#include <TerminalColors.h>
39#include <BicrystalActor.h>
40
41namespace gbLAB
42{
43
44
45 /**********************************************************************/
46 BicrystalActor::BicrystalActor(vtkGenericOpenGLRenderWindow* const renWin,vtkRenderer* const rndr) :
47 /* init */ renderWindow(renWin)
48 /* init */,renderer(rndr)
49 /* init */,mainLayout(new QGridLayout(this))
50 /* init */,showA(new QCheckBox(this))
51 /* init */,showB(new QCheckBox(this))
52
53// /* init */,showNodeLabels(new QCheckBox(this))
54// /* init */,showVelocities(new QCheckBox(this))
55// /* init */,velocityScaleEdit(new QLineEdit("1"))
56 /* init */,aPolyData(vtkSmartPointer<vtkPolyData>::New())
57 /* init */,aGlyphs(vtkSmartPointer<vtkGlyph3D>::New())
58 /* init */,aMapper(vtkSmartPointer<vtkPolyDataMapper>::New())
59 /* init */,aActor(vtkSmartPointer<vtkActor>::New())
60/* init */,bPolyData(vtkSmartPointer<vtkPolyData>::New())
61/* init */,bGlyphs(vtkSmartPointer<vtkGlyph3D>::New())
62/* init */,bMapper(vtkSmartPointer<vtkPolyDataMapper>::New())
63/* init */,bActor(vtkSmartPointer<vtkActor>::New())
64
65// /* init */,labelPolyData(vtkSmartPointer<vtkPolyData>::New())
66// /* init */,labelMapper(vtkSmartPointer<vtkLabeledDataMapper>::New())
67// /* init */,labelActor(vtkSmartPointer<vtkActor2D>::New())
68// /* init */,velocityPolyData(vtkSmartPointer<vtkPolyData>::New())
69// /* init */,velocityGlyphs(vtkSmartPointer<vtkGlyph3D>::New())
70// /* init */,velocityMapper(vtkSmartPointer<vtkPolyDataMapper>::New())
71// /* init */,velocityActor(vtkSmartPointer<vtkActor>::New())
72// /* init */,singleNodeLabelPolyData(vtkSmartPointer<vtkPolyData>::New())
73// /* init */,singleNodeLabelMapper(vtkSmartPointer<vtkLabeledDataMapper>::New())
74// /* init */,singleNodeLabelActor(vtkSmartPointer<vtkActor2D>::New())
75// /* init */,singleNodeID(0)
76// /* init */,nodeClr{{100,100,100},{0,255,255},{255,0,255},{1,1,1}}
77 {
78
79 showA->setText("lattice A");
80 showA->setChecked(true);
81 aActor->SetVisibility(true);
82 showB->setText("lattice B");
83 showB->setChecked(true);
84 bActor->SetVisibility(true);
85
86// showNodeLabels->setText("node labels");
87// showNodeLabels->setChecked(false);
88// labelActor->SetVisibility(false);
89//
90// showVelocities->setText("velocities");
91// showVelocities->setChecked(false);
92// velocityActor->SetVisibility(false);
93// velocityScaleEdit->setEnabled(false);
94
95
96 mainLayout->addWidget(showA,0,0,1,1);
97 mainLayout->addWidget(showB,1,0,1,1);
98
99// mainLayout->addWidget(showNodeLabels,1,0,1,1);
100// mainLayout->addWidget(showVelocities,2,0,1,1);
101// mainLayout->addWidget(velocityScaleEdit,2,1,1,1);
102 this->setLayout(mainLayout);
103
104 connect(showA,SIGNAL(stateChanged(int)), this, SLOT(modify()));
105 connect(showB,SIGNAL(stateChanged(int)), this, SLOT(modify()));
106
107// connect(showNodeLabels,SIGNAL(stateChanged(int)), this, SLOT(modify()));
108// connect(showVelocities,SIGNAL(stateChanged(int)), this, SLOT(modify()));
109// connect(velocityScaleEdit,SIGNAL(returnPressed()), this, SLOT(modify()));
110
111
112 aGlyphs->SetInputData(aPolyData);
113 aGlyphs->SetSourceConnection(vtkSmartPointer<vtkSphereSource>::New()->GetOutputPort());
114 aGlyphs->ScalingOn();
115 aGlyphs->SetScaleModeToScaleByVector();
116 aGlyphs->SetScaleFactor(0.25);
117 aGlyphs->SetColorModeToColorByScalar();
118 aGlyphs->Update();
119 aMapper->SetInputConnection(aGlyphs->GetOutputPort());
120 aActor->SetMapper(aMapper);
121 aActor->GetProperty()->SetColor(1.0, 0.0, 0.0); //(R,G,B)
122 renderer->AddActor(aActor);
123
124 bGlyphs->SetInputData(bPolyData);
125 bGlyphs->SetSourceConnection(vtkSmartPointer<vtkSphereSource>::New()->GetOutputPort());
126 bGlyphs->ScalingOn();
127 bGlyphs->SetScaleModeToScaleByVector();
128 bGlyphs->SetScaleFactor(0.25);
129 bGlyphs->SetColorModeToColorByScalar();
130 bGlyphs->Update();
131 bMapper->SetInputConnection(bGlyphs->GetOutputPort());
132 bActor->SetMapper(bMapper);
133 bActor->GetProperty()->SetColor(0.0, 0.0, 1.0); //(R,G,B)
134 renderer->AddActor(bActor);
135
136
137 // Labels
138// labelMapper->SetInputData(labelPolyData);
139// labelMapper->SetLabelModeToLabelScalars();
140// labelMapper->SetLabelFormat("%1.0f");
141// labelMapper->GetLabelTextProperty()->SetFontSize(20);
142// labelActor->SetMapper(labelMapper);
143// labelActor->GetProperty()->SetColor(0.0, 0.0, 0.0); //(R,G,B)
144//
145// // Velocities
146// velocityGlyphs->SetInputData(velocityPolyData);
147// velocityGlyphs->SetSourceConnection(vtkSmartPointer<vtkArrowSource>::New()->GetOutputPort());
148// velocityGlyphs->ScalingOn();
149// velocityGlyphs->SetScaleModeToScaleByVector();
150// velocityGlyphs->OrientOn();
151// velocityGlyphs->ClampingOff();
152// velocityGlyphs->SetVectorModeToUseVector();
153// velocityGlyphs->SetIndexModeToOff();
154// velocityMapper->SetInputConnection(velocityGlyphs->GetOutputPort());
155// velocityMapper->ScalarVisibilityOff();
156// velocityActor->SetMapper(velocityMapper);
157// velocityActor->GetProperty()->SetColor(1.0, 0.0, 1.0); //(R,G,B)
158//
159// // Single node Label
160// singleNodeLabelMapper->SetInputData(singleNodeLabelPolyData);
161// singleNodeLabelMapper->SetLabelModeToLabelScalars();
162// singleNodeLabelMapper->SetLabelFormat("%1.0f");
163// singleNodeLabelMapper->GetLabelTextProperty()->SetFontSize(20);
164// singleNodeLabelActor->SetMapper(singleNodeLabelMapper);
165// singleNodeLabelActor->GetProperty()->SetColor(1.0, 0.0, 0.0); //(R,G,B)
166// singleNodeLabelActor->VisibilityOff();
167
168// renderer->AddActor(velocityActor);
169// renderer->AddActor(labelActor);
170// renderer->AddActor(singleNodeLabelActor);
171
172 }
173
174
175
176 /**********************************************************************/
177 void BicrystalActor::updateConfiguration(const std::shared_ptr<BiCrystal<3>>& bc)
178 {// https://stackoverflow.com/questions/6878263/remove-individual-points-from-vtkpoints
179 std::cout<<"Updating BiCrystal..."<<std::flush;
180 const auto t0= std::chrono::system_clock::now();
181//
182 vtkSmartPointer<vtkPoints> aPoints(vtkSmartPointer<vtkPoints>::New());
183 vtkSmartPointer<vtkPoints> bPoints(vtkSmartPointer<vtkPoints>::New());
184
185 const int N=5;
186 for(int i=-N;i<N+1;++i)
187 {
188 for(int j=-N;j<N+1;++j)
189 {
190 for(int k=-N;k<N+1;++k)
191 {
192 const Eigen::Matrix<double,3,1> Pa(bc->A.latticeBasis*(Eigen::Matrix<double,3,1>()<<i,j,k).finished());
193 aPoints->InsertNextPoint(Pa.data());
194 const Eigen::Matrix<double,3,1> Pb(bc->B.latticeBasis*(Eigen::Matrix<double,3,1>()<<i,j,k).finished());
195 bPoints->InsertNextPoint(Pb.data());
196
197 }
198 }
199 }
200
201
202// vtkSmartPointer<vtkUnsignedCharArray> nodeColors(vtkSmartPointer<vtkUnsignedCharArray>::New());
203// nodeColors->SetNumberOfComponents(3);
204//
205// vtkSmartPointer<vtkDoubleArray> nodeLabels(vtkSmartPointer<vtkDoubleArray>::New());
206// nodeLabels->SetNumberOfComponents(1);
207// nodeLabels->SetName("node IDs");
208//
209// vtkSmartPointer<vtkPoints> singleNodePoint(vtkSmartPointer<vtkPoints>::New());
210// vtkSmartPointer<vtkDoubleArray> singleNodenodeLabels(vtkSmartPointer<vtkDoubleArray>::New());
211//
212// vtkSmartPointer<vtkDoubleArray> velocityVectors(vtkSmartPointer<vtkDoubleArray>::New());
213// velocityVectors->SetNumberOfComponents(3);
214// velocityVectors->SetName("nodeVelocity");
215//
216// for(const auto& node : configIO.nodes())
217// {
218// nodePoints->InsertNextPoint(node.P.data());
219// nodeLabels->InsertNextTuple1(node.sID);
220// velocityVectors->InsertNextTuple(node.V.data()); // arrow vector
221// nodeColors->InsertNextTypedTuple(node.meshLocation>2? this->nodeClr[3] : this->nodeClr[node.meshLocation]);
222//
223// // Single node
224// if(node.sID==singleNodeID)
225// {
226// singleNodePoint->InsertNextPoint(node.P.data());
227// singleNodenodeLabels->InsertNextTuple1(node.sID);
228// }
229// }
230//
231 aPolyData->SetPoints(aPoints);
232 aPolyData->Modified();
233 bPolyData->SetPoints(bPoints);
234 bPolyData->Modified();
235
236//
237// labelPolyData->SetPoints(nodePoints);
238// labelPolyData->GetPointData()->SetScalars(nodeLabels);
239// labelPolyData->Modified();
240//
241// singleNodeLabelPolyData->SetPoints(singleNodePoint);
242// singleNodeLabelPolyData->GetPointData()->SetScalars(singleNodenodeLabels);
243// singleNodeLabelPolyData->Modified();
244//
245// velocityPolyData->SetPoints(nodePoints);
246// velocityPolyData->GetPointData()->SetVectors(velocityVectors);
247// velocityPolyData->Modified();
248 renderer->ResetCamera();
249 renderWindow->Render();
250 std::cout<<magentaColor<<" ["<<(std::chrono::duration<double>(std::chrono::system_clock::now()-t0)).count()<<" sec]"<<defaultColor<<std::endl;
251 }
252
253 /**********************************************************************/
255 {
256
257 aActor->SetVisibility(showA->isChecked());
258 bActor->SetVisibility(showB->isChecked());
259
260// labelActor->SetVisibility(showNodeLabels->isChecked());
261// velocityActor->SetVisibility(showVelocities->isChecked());
262// velocityScaleEdit->setEnabled(showVelocities->isChecked());
263// const double vScaling(std::atof(velocityScaleEdit->text() .toStdString().c_str()));
264// velocityGlyphs->SetScaleFactor(vScaling);
265
266// velocityActor->SetScale(vScaling);
267// aGlyphs->SetScaleFactor(2.0*this->tubeRadius*1.2);
268//
269// if(this->showVelocities)
270// {
271// velocityActor->VisibilityOn();
272// }
273// else
274// {
275// velocityActor->VisibilityOff();
276// }
277//
278// if(this->showNodeIDs)
279// {
280// labelActor->VisibilityOn();
281//
282// }
283// else
284// {
285// labelActor->VisibilityOff();
286// }
287//
288// velocityGlyphs->SetScaleFactor(this->velocityFactor);
289//
290//
291// if(this->showSingleNode)
292// {
293// // HERE WE SHOULD CHANGE THE NODE POSITION BASED ON NODE ID
294// // OTHERWISE THE SELECTED NODE WILL BE VISIBLE ONLY UPON LOADING A NEW FRAME
295// std::cout<<"RELOAD FRAME TO SHOW SELECTED NODE"<<std::endl;
296// singleNodeLabelActor->VisibilityOn();
297// }
298// else
299// {
300// singleNodeLabelActor->VisibilityOff();
301// }
302//
303// if(this->showA)
304// {
305// aActor->VisibilityOn();
306// }
307// else
308// {
309// aActor->VisibilityOff();
310// }
311
312 renderWindow->Render();
313 }
314
315
316} // namespace model
317#endif
static std::string defaultColor
static std::string magentaColor
vtkSmartPointer< vtkPolyData > bPolyData
BicrystalActor(vtkGenericOpenGLRenderWindow *const, vtkRenderer *const)
void updateConfiguration(const std::shared_ptr< BiCrystal< 3 > > &bc)
vtkSmartPointer< vtkGlyph3D > aGlyphs
vtkSmartPointer< vtkPolyDataMapper > aMapper
vtkSmartPointer< vtkActor > bActor
vtkSmartPointer< vtkPolyData > aPolyData
QGridLayout * mainLayout
vtkGenericOpenGLRenderWindow *const renderWindow
vtkSmartPointer< vtkPolyDataMapper > bMapper
vtkSmartPointer< vtkGlyph3D > bGlyphs
vtkRenderer *const renderer
vtkSmartPointer< vtkActor > aActor