8#include <TPolyLine3D.h>
9#include <TPolyMarker3D.h>
14#include <TVirtualViewer3D.h>
22 m_driftLines.reserve(1000);
23 m_tracks.reserve(100);
40 m_markerSizeCluster = size;
42 std::cerr <<
m_className <<
"::SetClusterMarkerSize: Size must be > 0.\n";
48 m_markerSizeCollision = size;
50 std::cerr <<
m_className <<
"::SetCollisionMarkerSize: Size must be > 0.\n";
55 std::vector<std::array<float, 3> >& driftLine,
56 bool& electron)
const {
58 if (i >= m_driftLines.size())
return;
59 std::copy(m_driftLines[i].first.begin(), m_driftLines[i].first.end(),
60 std::back_inserter(driftLine));
61 if (m_driftLines[i].second == Particle::Electron) {
69 const float x0,
const float y0,
71 std::lock_guard<std::mutex> guard(m_mutex);
73 std::array<float, 3> p = {x0, y0, z0};
74 constexpr size_t minSize = 1;
75 std::vector<std::array<float, 3> > dl(std::max(minSize, np), p);
76 m_driftLines.push_back(std::make_pair(std::move(dl), Particle::Electron));
78 id = m_driftLines.size() - 1;
82 const float y0,
const float z0) {
83 std::lock_guard<std::mutex> guard(m_mutex);
84 std::array<float, 3> p = {x0, y0, z0};
85 constexpr size_t minSize = 1;
86 std::vector<std::array<float, 3> > dl(std::max(minSize, np), p);
87 m_driftLines.push_back(std::make_pair(std::move(dl), Particle::Hole));
89 id = m_driftLines.size() - 1;
93 const float y0,
const float z0) {
94 std::lock_guard<std::mutex> guard(m_mutex);
95 std::array<float, 3> p = {x0, y0, z0};
96 constexpr size_t minSize = 1;
97 std::vector<std::array<float, 3> > dl(std::max(minSize, np), p);
98 m_driftLines.push_back(std::make_pair(std::move(dl), Particle::Ion));
100 id = m_driftLines.size() - 1;
104 const float x1,
const float y1,
const float z1) {
105 std::lock_guard<std::mutex> guard(m_mutex);
106 std::array<float, 3> p0 = {x0, y0, z0};
107 std::array<float, 3> p1 = {x1, y1, z1};
108 m_photons.push_back({p0, p1});
112 const float x0,
const float y0,
114 std::lock_guard<std::mutex> guard(m_mutex);
116 constexpr size_t minSize = 1;
117 std::vector<std::array<float, 3> > track(std::max(minSize, np));
118 track[0] = {x0, y0, z0};
119 m_tracks.push_back(std::move(track));
121 id = m_tracks.size() - 1;
125 const float x,
const float y,
127 std::lock_guard<std::mutex> guard(m_mutex);
128 if (iL >= m_driftLines.size() || iP >= m_driftLines[iL].first.size()) {
129 std::cerr <<
m_className <<
"::SetDriftLinePoint: Index out of range.\n";
132 m_driftLines[iL].first[iP] = {x, y, z};
136 const float y,
const float z) {
137 std::lock_guard<std::mutex> guard(m_mutex);
138 if (iL >= m_driftLines.size()) {
139 std::cerr <<
m_className <<
"::AddDriftLinePoint: Index out of range.\n";
142 std::array<float, 3> p = {x, y, z};
143 m_driftLines[iL].first.push_back(std::move(p));
147 const float x,
const float y,
const float z) {
148 std::lock_guard<std::mutex> guard(m_mutex);
149 if (iL >= m_tracks.size() || iP >= m_tracks[iL].size()) {
150 std::cerr <<
m_className <<
"::SetTrackPoint: Index out of range.\n";
153 m_tracks[iL][iP] = {x, y, z};
157 const float y,
const float z) {
158 std::lock_guard<std::mutex> guard(m_mutex);
159 if (iL >= m_tracks.size()) {
160 std::cerr <<
m_className <<
"::AddTrackPoint: Index out of range.\n";
163 std::array<float, 3> p = {x, y, z};
164 m_tracks[iL].push_back(std::move(p));
168 std::lock_guard<std::mutex> guard(m_mutex);
169 std::array<float, 3> p = {x, y, z};
170 m_exc.push_back(std::move(p));
174 std::lock_guard<std::mutex> guard(m_mutex);
175 std::array<float, 3> p = {x, y, z};
176 m_ion.push_back(std::move(p));
180 std::lock_guard<std::mutex> guard(m_mutex);
181 std::array<float, 3> p = {x, y, z};
182 m_att.push_back(std::move(p));
196 pad->SetTitle(
"Drift lines");
198 const bool rangeSet =
RangeSet(pad);
199 if (axis || !rangeSet) {
201 if (!SetPlotLimits2d()) {
203 <<
" Could not determine the plot limits.\n";
210 frame->GetXaxis()->SetTitle(
LabelX().c_str());
211 frame->GetYaxis()->SetTitle(
LabelY().c_str());
212 }
else if (!rangeSet) {
216 for (
const auto& driftLine : m_driftLines) {
218 if (driftLine.second == Particle::Electron) {
219 DrawLine(driftLine.first, m_colElectron, lw);
220 }
else if (driftLine.second == Particle::Hole) {
221 DrawLine(driftLine.first, m_colHole, lw);
223 DrawLine(driftLine.first, m_colIon, lw);
229 gr.SetMarkerColor(m_colTrack);
230 gr.SetMarkerSize(m_markerSizeCluster);
231 for (
const auto& track : m_tracks) {
233 if (!m_drawClusters)
continue;
234 std::vector<float> xgr;
235 std::vector<float> ygr;
236 for (
const auto& p : track) {
237 if (!
InBox(p))
continue;
238 float xp = 0., yp = 0.;
239 ToPlane(p[0], p[1], p[2], xp, yp);
243 gr.DrawGraph(xgr.size(), xgr.data(), ygr.data(),
"Psame");
246 gr.SetLineColor(m_colPhoton);
248 for (
const auto& photon : m_photons) {
249 float xp0 = 0., yp0 = 0.;
250 float xp1 = 0., yp1 = 0.;
251 ToPlane(photon[0][0], photon[0][1], photon[0][2], xp0, yp0);
252 ToPlane(photon[1][0], photon[1][1], photon[1][2], xp1, yp1);
253 std::vector<float> xgr = {xp0, xp1};
254 std::vector<float> ygr = {yp0, yp1};
255 gr.DrawGraph(2, xgr.data(), ygr.data(),
"Lsame");
258 gr.SetMarkerSize(m_markerSizeCollision);
259 gr.SetMarkerStyle(20);
260 if (!m_exc.empty()) {
261 gr.SetMarkerColor(m_colExcitation);
262 std::vector<float> xgr;
263 std::vector<float> ygr;
264 for (
const auto& p : m_exc) {
265 if (!
InBox(p))
continue;
266 float xp = 0., yp = 0.;
267 ToPlane(p[0], p[1], p[2], xp, yp);
272 gr.DrawGraph(xgr.size(), xgr.data(), ygr.data(),
"Psame");
275 if (!m_ion.empty()) {
276 gr.SetMarkerColor(m_colIonisation);
277 std::vector<float> xgr;
278 std::vector<float> ygr;
279 for (
const auto& p : m_ion) {
280 if (!
InBox(p))
continue;
281 float xp = 0., yp = 0.;
282 ToPlane(p[0], p[1], p[2], xp, yp);
287 gr.DrawGraph(xgr.size(), xgr.data(), ygr.data(),
"Psame");
290 if (!m_att.empty()) {
291 gr.SetMarkerColor(m_colAttachment);
292 std::vector<float> xgr;
293 std::vector<float> ygr;
294 for (
const auto& p : m_att) {
295 if (!
InBox(p))
continue;
296 float xp = 0., yp = 0.;
297 ToPlane(p[0], p[1], p[2], xp, yp);
302 gr.DrawGraph(xgr.size(), xgr.data(), ygr.data(),
"Psame");
312 pad->SetTitle(
"Drift lines");
313 if (!pad->GetView()) {
314 if (!SetPlotLimits3d()) {
316 <<
" Could not determine the plot limits.\n";
318 auto view = TView::CreateView(1, 0, 0);
321 if (axis) view->ShowAxis();
323 if (ogl) pad->GetViewer3D(
"ogl");
325 for (
const auto& driftLine : m_driftLines) {
326 std::vector<float> points;
327 for (
const auto& p : driftLine.first) {
328 points.push_back(p[0]);
329 points.push_back(p[1]);
330 points.push_back(p[2]);
332 const int nP = driftLine.first.size();
333 TPolyLine3D pl(nP, points.data());
334 if (driftLine.second == Particle::Electron) {
335 pl.SetLineColor(m_colElectron);
336 }
else if (driftLine.second == Particle::Hole) {
337 pl.SetLineColor(m_colHole);
339 pl.SetLineColor(m_colIon);
342 pl.DrawPolyLine(nP, points.data(),
"same");
345 for (
const auto& track : m_tracks) {
346 std::vector<float> points;
347 for (
const auto& p : track) {
348 points.push_back(p[0]);
349 points.push_back(p[1]);
350 points.push_back(p[2]);
352 const int nP = track.size();
353 TPolyLine3D pl(nP, points.data());
354 pl.SetLineColor(m_colTrack);
356 pl.DrawPolyLine(nP, points.data(),
"same");
357 if (!m_drawClusters)
continue;
358 TPolyMarker3D pm(nP, points.data());
359 pm.SetMarkerColor(m_colTrack);
360 pm.SetMarkerSize(m_markerSizeCluster);
361 pm.DrawPolyMarker(nP, points.data(), 20,
"same");
364 if (!m_exc.empty()) {
365 const size_t nP = m_exc.size();
366 std::vector<float> points;
367 for (
size_t i = 0; i < nP; ++i) {
368 points.push_back(m_exc[i][0]);
369 points.push_back(m_exc[i][1]);
370 points.push_back(m_exc[i][2]);
372 TPolyMarker3D pm(nP, points.data());
373 pm.SetMarkerColor(m_colExcitation);
374 pm.SetMarkerSize(m_markerSizeCollision);
375 pm.DrawPolyMarker(nP, points.data(), 20,
"same");
378 if (!m_ion.empty()) {
379 const size_t nP = m_ion.size();
380 std::vector<float> points;
381 for (
size_t i = 0; i < nP; ++i) {
382 points.push_back(m_ion[i][0]);
383 points.push_back(m_ion[i][1]);
384 points.push_back(m_ion[i][2]);
386 TPolyMarker3D pm(nP, points.data());
387 pm.SetMarkerColor(m_colIonisation);
388 pm.SetMarkerSize(m_markerSizeCollision);
389 pm.DrawPolyMarker(nP, points.data(), 20,
"same");
392 if (!m_att.empty()) {
393 const size_t nP = m_att.size();
394 std::vector<float> points;
395 for (
size_t i = 0; i < nP; ++i) {
396 points.push_back(m_att[i][0]);
397 points.push_back(m_att[i][1]);
398 points.push_back(m_att[i][2]);
400 TPolyMarker3D pm(nP, points.data());
401 pm.SetMarkerColor(m_colAttachment);
402 pm.SetMarkerSize(m_markerSizeCollision);
403 pm.DrawPolyMarker(nP, points.data(), 20,
"same");
409 TAxis3D* ax3d = TAxis3D::GetPadAxis();
411 ax3d->SetLabelColor(kGray + 2);
412 ax3d->SetAxisColor(kGray + 2);
413 ax3d->SetXTitle(
"x");
414 ax3d->SetYTitle(
"y");
415 ax3d->SetZTitle(
"z");
421bool ViewDrift::SetPlotLimits2d() {
424 double xmin = 0., ymin = 0., xmax = 0., ymax = 0.;
436 std::array<double, 3> bbmin;
437 std::array<double, 3> bbmax;
438 bbmin.fill(std::numeric_limits<double>::max());
439 bbmax.fill(-std::numeric_limits<double>::max());
440 for (
const auto& driftLine : m_driftLines) {
441 for (
const auto& p : driftLine.first) {
442 for (
unsigned int i = 0; i < 3; ++i) {
443 bbmin[i] = std::min(bbmin[i],
double(p[i]));
444 bbmax[i] = std::max(bbmax[i],
double(p[i]));
448 for (
const auto& track : m_tracks) {
449 for (
const auto& p : track) {
450 for (
unsigned int i = 0; i < 3; ++i) {
451 bbmin[i] = std::min(bbmin[i],
double(p[i]));
452 bbmax[i] = std::max(bbmax[i],
double(p[i]));
460bool ViewDrift::SetPlotLimits3d() {
463 if (m_driftLines.empty() && m_tracks.empty())
return false;
465 std::array<double, 3> bbmin;
466 std::array<double, 3> bbmax;
467 bbmin.fill(std::numeric_limits<double>::max());
468 bbmax.fill(-std::numeric_limits<double>::max());
469 for (
const auto& driftLine : m_driftLines) {
470 for (
const auto& p : driftLine.first) {
471 for (
unsigned int i = 0; i < 3; ++i) {
472 bbmin[i] = std::min(bbmin[i],
double(p[i]));
473 bbmax[i] = std::max(bbmax[i],
double(p[i]));
477 for (
const auto& track : m_tracks) {
478 for (
const auto& p : track) {
479 for (
unsigned int i = 0; i < 3; ++i) {
480 bbmin[i] = std::min(bbmin[i],
double(p[i]));
481 bbmax[i] = std::max(bbmax[i],
double(p[i]));
486 for (
size_t i = 0; i < 3; ++i) r +=
fabs(bbmax[i] - bbmin[i]);
Base class for visualization classes.
bool InBox(const std::array< T, 3 > &x) const
void DrawLine(const std::vector< std::array< float, 3 > > &xl, const short col, const short lw)
TPad * GetCanvas()
Retrieve the canvas.
void SetRange(TPad *pad, const double x0, const double y0, const double x1, const double y1)
bool PlotLimits(Sensor *sensor, double &xmin, double &ymin, double &xmax, double &ymax) const
void ToPlane(const T x, const T y, const T z, T &xp, T &yp) const
bool PlotLimitsFromUserBox(double &xmin, double &ymin, double &xmax, double &ymax) const
void AddIonisation(const float x, const float y, const float z)
void NewHoleDriftLine(const size_t np, size_t &id, const float x0, const float y0, const float z0)
void AddDriftLinePoint(const size_t iL, const float x, const float y, const float z)
void Plot(const bool twod=false, const bool axis=true)
Draw the drift lines.
void AddTrackPoint(const size_t iL, const float x, const float y, const float z)
void NewChargedParticleTrack(const size_t np, size_t &id, const float x0, const float y0, const float z0)
void NewIonDriftLine(const size_t np, size_t &id, const float x0, const float y0, const float z0)
void SetTrackPoint(const size_t iL, const size_t iP, const float x, const float y, const float z)
void SetClusterMarkerSize(const double size)
Set the size of the cluster markers (see TAttMarker).
void SetDriftLinePoint(const size_t iL, const size_t iP, const float x, const float y, const float z)
void Plot3d(const bool axis, const bool ogl)
Make a 3D plot of the drift lines.
void Plot2d(const bool axis)
Make a 2D plot of the drift lines in the current viewing plane.
void AddAttachment(const float x, const float y, const float z)
void GetDriftLine(const size_t i, std::vector< std::array< float, 3 > > &driftLine, bool &electron) const
Retrieve the coordinates of a given drift line.
void Clear()
Delete existing drift lines, tracks and markers.
void AddPhoton(const float x0, const float y0, const float z0, const float x1, const float y1, const float z1)
void NewElectronDriftLine(const size_t np, size_t &id, const float x0, const float y0, const float z0)
void AddExcitation(const float x, const float y, const float z)
void SetCollisionMarkerSize(const double size)
Set the size of the collision markers (see TAttMarker).
DoubleAc fabs(const DoubleAc &f)