Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
Garfield::ViewBase Class Reference

Base class for visualization classes. More...

#include <ViewBase.hh>

+ Inheritance diagram for Garfield::ViewBase:

Public Member Functions

 ViewBase ()=delete
 Default constructor.
 
 ViewBase (const std::string &name)
 Constructor.
 
virtual ~ViewBase ()=default
 Destructor.
 
void SetCanvas (TPad *pad)
 Set the canvas to be painted on.
 
void SetCanvas ()
 Unset an external canvas.
 
TPad * GetCanvas ()
 Retrieve the canvas.
 
void SetArea (const double xmin, const double ymin, const double xmax, const double ymax)
 
virtual void SetArea (const double xmin, const double ymin, const double zmin, const double xmax, const double ymax, const double zmax)
 Set a bounding box (if applicable).
 
void SetArea ()
 
virtual void SetPlane (const double fx, const double fy, const double fz, const double x0, const double y0, const double z0)
 
virtual void SetPlane (const double fx, const double fy, const double fz, const double x0, const double y0, const double z0, const double hx, const double hy, const double hz)
 Set the projection plane specifying a hint for the in-plane x axis.
 
void Rotate (const double angle)
 Rotate the viewing plane (angle in radian).
 
void SetPlaneXY ()
 Set the viewing plane to x-y.
 
void SetPlaneXZ ()
 Set the viewing plane to x-z.
 
void SetPlaneYZ ()
 Set the viewing plane to y-z.
 
void EnableDebugging (const bool on=true)
 Switch on/off debugging output.
 

Static Public Member Functions

static std::string FindUnusedFunctionName (const std::string &s)
 Find an unused function name.
 
static std::string FindUnusedHistogramName (const std::string &s)
 Find an unused histogram name.
 
static std::string FindUnusedCanvasName (const std::string &s)
 Find an unused canvas name.
 

Protected Member Functions

void UpdateProjectionMatrix ()
 
template<typename T >
void ToPlane (const T x, const T y, const T z, T &xp, T &yp) const
 
template<typename T >
bool InBox (const std::array< T, 3 > &x) const
 
void Clip (const std::array< float, 3 > &x0, const std::array< float, 3 > &x1, std::array< float, 3 > &xc) const
 
void DrawLine (const std::vector< std::array< float, 3 > > &xl, const short col, const short lw)
 
std::string LabelX ()
 
std::string LabelY ()
 
std::string PlaneDescription ()
 
bool PlotLimits (Sensor *sensor, double &xmin, double &ymin, double &xmax, double &ymax) const
 
bool PlotLimits (Component *cmp, double &xmin, double &ymin, double &xmax, double &ymax) const
 
bool PlotLimitsFromUserBox (double &xmin, double &ymin, double &xmax, double &ymax) const
 
bool PlotLimits (std::array< double, 3 > &bbmin, std::array< double, 3 > &bbmax, double &xmin, double &ymin, double &xmax, double &ymax) const
 
bool RangeSet (TPad *)
 
void SetRange (TPad *pad, const double x0, const double y0, const double x1, const double y1)
 

Protected Attributes

std::string m_className = "ViewBase"
 
bool m_debug = false
 
bool m_userPlotLimits = false
 
double m_xMinPlot = -1.
 
double m_xMaxPlot = 1.
 
double m_yMinPlot = -1.
 
double m_yMaxPlot = 1.
 
bool m_userBox = false
 
double m_xMinBox = -1.
 
double m_xMaxBox = 1.
 
double m_yMinBox = -1.
 
double m_yMaxBox = 1.
 
double m_zMinBox = -1.
 
double m_zMaxBox = 1.
 
std::array< std::array< double, 3 >, 3 > m_proj
 
std::array< double, 4 > m_plane {{0, 0, 1, 0}}
 
std::array< std::array< double, 3 >, 3 > m_prmat
 

Detailed Description

Base class for visualization classes.

Definition at line 18 of file ViewBase.hh.

Constructor & Destructor Documentation

◆ ViewBase() [1/2]

Garfield::ViewBase::ViewBase ( )
delete

Default constructor.

◆ ViewBase() [2/2]

Garfield::ViewBase::ViewBase ( const std::string &  name)

Constructor.

Definition at line 68 of file ViewBase.cc.

68 :
69 m_className(name) {
70
72}
void SetDefaultStyle()
Apply the default Garfield ROOT style.
std::string m_className
Definition: ViewBase.hh:78
PlottingEngine plottingEngine

◆ ~ViewBase()

virtual Garfield::ViewBase::~ViewBase ( )
virtualdefault

Destructor.

Member Function Documentation

◆ Clip()

void Garfield::ViewBase::Clip ( const std::array< float, 3 > &  x0,
const std::array< float, 3 > &  x1,
std::array< float, 3 > &  xc 
) const
protected

Definition at line 347 of file ViewBase.cc.

349 {
350
351 xc.fill(0.);
352 const bool in0 = InBox(x0);
353 const bool in1 = InBox(x1);
354 if (in0 == in1) return;
355 xc = in0 ? x1 : x0;
356 const std::array<float, 3> dx = {x1[0] - x0[0], x1[1] - x0[1],
357 x1[2] - x0[2]};
358 std::array<double, 3> bmin = {m_xMinBox, m_yMinBox, m_zMinBox};
359 std::array<double, 3> bmax = {m_xMaxBox, m_yMaxBox, m_zMaxBox};
360 for (size_t i = 0; i < 3; ++i) {
361 if (dx[i] == 0. || (xc[i] >= bmin[i] && xc[i] <= bmax[i])) continue;
362 const double b = xc[i] < bmin[i] ? bmin[i] : bmax[i];
363 const double s = (b - xc[i]) / dx[i];
364 xc[i] = b;
365 for (size_t j = 0; j < 3; ++j) {
366 if (j != i) xc[j] += dx[j] * s;
367 }
368 }
369}
bool InBox(const std::array< T, 3 > &x) const
Definition: ViewBase.hh:115

Referenced by DrawLine().

◆ DrawLine()

void Garfield::ViewBase::DrawLine ( const std::vector< std::array< float, 3 > > &  xl,
const short  col,
const short  lw 
)
protected

Definition at line 371 of file ViewBase.cc.

372 {
373
374 const size_t nP = xl.size();
375 if (nP < 2) return;
376
377 TGraph gr;
378 gr.SetLineColor(col);
379 gr.SetLineWidth(lw);
380
381 std::vector<float> xgr;
382 std::vector<float> ygr;
383 auto x0 = xl[0];
384 bool in0 = InBox(x0);
385 if (in0) {
386 float xp = 0., yp = 0.;
387 ToPlane(x0[0], x0[1], x0[2], xp, yp);
388 xgr.push_back(xp);
389 ygr.push_back(yp);
390 }
391 for (unsigned int j = 1; j < nP; ++j) {
392 auto x1 = xl[j];
393 bool in1 = InBox(x1);
394 if (in1 != in0) {
395 float xp = 0., yp = 0.;
396 std::array<float, 3> xc;
397 Clip(x0, x1, xc);
398 ToPlane(xc[0], xc[1], xc[2], xp, yp);
399 xgr.push_back(xp);
400 ygr.push_back(yp);
401 }
402 if (in1) {
403 float xp = 0., yp = 0.;
404 ToPlane(x1[0], x1[1], x1[2], xp, yp);
405 xgr.push_back(xp);
406 ygr.push_back(yp);
407 } else if (!xgr.empty()) {
408 gr.DrawGraph(xgr.size(), xgr.data(), ygr.data(), "Lsame");
409 xgr.clear();
410 ygr.clear();
411 }
412 x0 = x1;
413 in0 = in1;
414 }
415 if (!xgr.empty()) {
416 gr.DrawGraph(xgr.size(), xgr.data(), ygr.data(), "Lsame");
417 }
418}
void Clip(const std::array< float, 3 > &x0, const std::array< float, 3 > &x1, std::array< float, 3 > &xc) const
Definition: ViewBase.cc:347
void ToPlane(const T x, const T y, const T z, T &xp, T &yp) const
Definition: ViewBase.hh:109

Referenced by Garfield::ViewDrift::Plot2d(), and Garfield::ViewField::PlotFieldLines().

◆ EnableDebugging()

void Garfield::ViewBase::EnableDebugging ( const bool  on = true)
inline

Switch on/off debugging output.

Definition at line 68 of file ViewBase.hh.

68{ m_debug = on; }

◆ FindUnusedCanvasName()

std::string Garfield::ViewBase::FindUnusedCanvasName ( const std::string &  s)
static

Find an unused canvas name.

Definition at line 310 of file ViewBase.cc.

310 {
311 int idx = 0;
312 std::string hname = s + "_0";
313 while (gROOT->GetListOfCanvases()->FindObject(hname.c_str())) {
314 ++idx;
315 hname = s + "_" + std::to_string(idx);
316 }
317 return hname;
318}

Referenced by GetCanvas(), Garfield::ComponentAnalyticField::MultipoleMoments(), Garfield::TrackSrim::PlotEnergyLoss(), Garfield::TrackSrim::PlotRange(), and Garfield::TrackSrim::PlotStraggling().

◆ FindUnusedFunctionName()

std::string Garfield::ViewBase::FindUnusedFunctionName ( const std::string &  s)
static

Find an unused function name.

Definition at line 290 of file ViewBase.cc.

290 {
291 int idx = 0;
292 std::string fname = s + "_0";
293 while (gROOT->GetListOfFunctions()->FindObject(fname.c_str())) {
294 ++idx;
295 fname = s + "_" + std::to_string(idx);
296 }
297 return fname;
298}

◆ FindUnusedHistogramName()

std::string Garfield::ViewBase::FindUnusedHistogramName ( const std::string &  s)
static

Find an unused histogram name.

Definition at line 300 of file ViewBase.cc.

300 {
301 int idx = 0;
302 std::string hname = s + "_0";
303 while (gDirectory->GetList()->FindObject(hname.c_str())) {
304 ++idx;
305 hname = s + "_" + std::to_string(idx);
306 }
307 return hname;
308}

Referenced by Garfield::ViewSignal::Plot(), and Garfield::ViewSignal::PlotSignal().

◆ GetCanvas()

TPad * Garfield::ViewBase::GetCanvas ( )

Retrieve the canvas.

Definition at line 74 of file ViewBase.cc.

74 {
75 if (!m_pad) {
76 std::string name = FindUnusedCanvasName("c" + m_className);
77 if (!m_canvas) m_canvas.reset(new TCanvas(name.c_str(), ""));
78 m_pad = m_canvas.get();
79 }
80 return m_pad;
81}
static std::string FindUnusedCanvasName(const std::string &s)
Find an unused canvas name.
Definition: ViewBase.cc:310

Referenced by Garfield::ViewFEMesh::Plot(), Garfield::ViewSignal::Plot(), Garfield::ViewGeometry::Plot2d(), Garfield::ViewDrift::Plot2d(), Garfield::ViewGeometry::Plot3d(), Garfield::ViewDrift::Plot3d(), Garfield::ViewField::PlotFieldLines(), Garfield::ViewIsochrons::PlotIsochrons(), and Garfield::ViewSignal::PlotSignal().

◆ InBox()

template<typename T >
bool Garfield::ViewBase::InBox ( const std::array< T, 3 > &  x) const
inlineprotected

Definition at line 115 of file ViewBase.hh.

115 {
116 if (!m_userBox) return true;
117 if (x[0] < m_xMinBox || x[0] > m_xMaxBox ||
118 x[1] < m_yMinBox || x[1] > m_yMaxBox ||
119 x[2] < m_yMinBox || x[2] > m_zMaxBox) return false;
120 return true;
121 }

Referenced by Clip(), DrawLine(), and Garfield::ViewDrift::Plot2d().

◆ LabelX()

std::string Garfield::ViewBase::LabelX ( )
protected

Definition at line 420 of file ViewBase.cc.

420 {
421
422 std::string xLabel = "";
423 constexpr double tol = 1.e-4;
424 // x portion
425 if (fabs(m_proj[0][0] - 1) < tol) {
426 xLabel = "#it{x}";
427 } else if (fabs(m_proj[0][0] + 1) < tol) {
428 xLabel = "#minus#it{x}";
429 } else if (fabs(m_proj[0][0]) > tol) {
430 xLabel = Fmt(m_proj[0][0]) + " #it{x}";
431 }
432
433 // y portion
434 if (!xLabel.empty()) {
435 if (m_proj[0][1] < -tol) {
436 xLabel += " #minus ";
437 } else if (m_proj[0][1] > tol) {
438 xLabel += " #plus ";
439 }
440 if (fabs(m_proj[0][1] - 1) < tol || fabs(m_proj[0][1] + 1) < tol) {
441 xLabel += "#it{y}";
442 } else if (fabs(m_proj[0][1]) > tol) {
443 xLabel += Fmt(fabs(m_proj[0][1])) + " #it{y}";
444 }
445 } else {
446 if (fabs(m_proj[0][1] - 1) < tol) {
447 xLabel = "#it{y}";
448 } else if (fabs(m_proj[0][1] + 1) < tol) {
449 xLabel = "#minus#it{y}";
450 } else if (fabs(m_proj[0][1]) > tol) {
451 xLabel = Fmt(m_proj[0][1]) + " #it{y}";
452 }
453 }
454
455 // z portion
456 if (!xLabel.empty()) {
457 if (m_proj[0][2] < -tol) {
458 xLabel += " #minus ";
459 } else if (m_proj[0][2] > tol) {
460 xLabel += " #plus ";
461 }
462 if (fabs(m_proj[0][2] - 1) < tol || fabs(m_proj[0][2] + 1) < tol) {
463 xLabel += "#it{z}";
464 } else if (fabs(m_proj[0][2]) > tol) {
465 xLabel += Fmt(fabs(m_proj[0][2])) + " #it{z}";
466 }
467 } else {
468 if (fabs(m_proj[0][2] - 1) < tol) {
469 xLabel = "#it{z}";
470 } else if (fabs(m_proj[0][2] + 1) < tol) {
471 xLabel = "#minus#it{z}";
472 } else if (fabs(m_proj[0][2]) > tol) {
473 xLabel = Fmt(m_proj[0][2]) + " #it{z}";
474 }
475 }
476
477 // Unit
478 xLabel += " [cm]";
479 return xLabel;
480
481}
std::array< std::array< double, 3 >, 3 > m_proj
Definition: ViewBase.hh:96
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:615

Referenced by Garfield::ViewFEMesh::CreateDefaultAxes(), Garfield::ViewFEMesh::Plot(), Garfield::ViewDrift::Plot2d(), Garfield::ViewField::PlotFieldLines(), and Garfield::ViewIsochrons::PlotIsochrons().

◆ LabelY()

std::string Garfield::ViewBase::LabelY ( )
protected

Definition at line 483 of file ViewBase.cc.

483 {
484
485 std::string yLabel = "";
486 constexpr double tol = 1.e-4;
487 // x portion
488 if (fabs(m_proj[1][0] - 1) < tol) {
489 yLabel = "#it{x}";
490 } else if (fabs(m_proj[1][0] + 1) < tol) {
491 yLabel = "#minus#it{x}";
492 } else if (fabs(m_proj[1][0]) > tol) {
493 yLabel = Fmt(m_proj[1][0]) + " #it{x}";
494 }
495
496 // y portion
497 if (!yLabel.empty()) {
498 if (m_proj[1][1] < -tol) {
499 yLabel += " #minus ";
500 } else if (m_proj[1][1] > tol) {
501 yLabel += " #plus ";
502 }
503 if (fabs(m_proj[1][1] - 1) < tol || fabs(m_proj[1][1] + 1) < tol) {
504 yLabel += "#it{y}";
505 } else if (fabs(m_proj[1][1]) > tol) {
506 yLabel += Fmt(fabs(m_proj[1][1])) + " #it{y}";
507 }
508 } else {
509 if (fabs(m_proj[1][1] - 1) < tol) {
510 yLabel = "#it{y}";
511 } else if (fabs(m_proj[1][1] + 1) < tol) {
512 yLabel = "#minus#it{y}";
513 } else if (fabs(m_proj[1][1]) > tol) {
514 yLabel = Fmt(m_proj[1][1]) + " #it{y}";
515 }
516 }
517
518 // z portion
519 if (!yLabel.empty()) {
520 if (m_proj[1][2] < -tol) {
521 yLabel += " #minus ";
522 } else if (m_proj[1][2] > tol) {
523 yLabel += " #plus ";
524 }
525 if (fabs(m_proj[1][2] - 1) < tol || fabs(m_proj[1][2] + 1) < tol) {
526 yLabel += "#it{z}";
527 } else if (fabs(m_proj[1][2]) > tol) {
528 yLabel += Fmt(fabs(m_proj[1][2])) + " #it{z}";
529 }
530 } else {
531 if (fabs(m_proj[1][2] - 1) < tol) {
532 yLabel = "#it{z}";
533 } else if (fabs(m_proj[1][2] + 1) < tol) {
534 yLabel = "#minus#it{z}";
535 } else if (fabs(m_proj[1][2]) > tol) {
536 yLabel = Fmt(m_proj[1][2]) + " #it{z}";
537 }
538 }
539
540 // Unit
541 yLabel += " [cm]";
542 return yLabel;
543}

Referenced by Garfield::ViewFEMesh::CreateDefaultAxes(), Garfield::ViewFEMesh::Plot(), Garfield::ViewDrift::Plot2d(), Garfield::ViewField::PlotFieldLines(), and Garfield::ViewIsochrons::PlotIsochrons().

◆ PlaneDescription()

std::string Garfield::ViewBase::PlaneDescription ( )
protected

Definition at line 545 of file ViewBase.cc.

545 {
546
547 std::string description;
548
549 constexpr double tol = 1.e-4;
550 // x portion
551 if (fabs(m_plane[0] - 1) < tol) {
552 description = "x";
553 } else if (fabs(m_plane[0] + 1) < tol) {
554 description = "-x";
555 } else if (fabs(m_plane[0]) > tol) {
556 description = Fmt(m_plane[0]) + " x";
557 }
558
559 // y portion
560 if (!description.empty()) {
561 if (m_plane[1] < -tol) {
562 description += " - ";
563 } else if (m_plane[1] > tol) {
564 description += " + ";
565 }
566 if (fabs(m_plane[1] - 1) < tol || fabs(m_plane[1] + 1) < tol) {
567 description += "y";
568 } else if (fabs(m_plane[1]) > tol) {
569 description += Fmt(fabs(m_plane[1])) + " y";
570 }
571 } else {
572 if (fabs(m_plane[1] - 1) < tol) {
573 description = "y";
574 } else if (fabs(m_plane[1] + 1) < tol) {
575 description = "-y";
576 } else if (fabs(m_plane[1]) > tol) {
577 description = Fmt(m_plane[1]) + " y";
578 }
579 }
580
581 // z portion
582 if (!description.empty()) {
583 if (m_plane[2] < -tol) {
584 description += " - ";
585 } else if (m_plane[2] > tol) {
586 description += " + ";
587 }
588 if (fabs(m_plane[2] - 1) < tol || fabs(m_plane[2] + 1) < tol) {
589 description += "z";
590 } else if (fabs(m_plane[2]) > tol) {
591 description += Fmt(fabs(m_plane[2])) + " z";
592 }
593 } else {
594 if (fabs(m_plane[2] - 1) < tol) {
595 description = "z";
596 } else if (fabs(m_plane[2] + 1) < tol) {
597 description = "-z";
598 } else if (fabs(m_plane[2]) > tol) {
599 description = Fmt(m_plane[2]) + " z";
600 }
601 }
602
603 // Constant
604 description += " = " + Fmt(m_plane[3]);
605 return description;
606}
std::array< double, 4 > m_plane
Definition: ViewBase.hh:99

◆ PlotLimits() [1/3]

bool Garfield::ViewBase::PlotLimits ( Component cmp,
double &  xmin,
double &  ymin,
double &  xmax,
double &  ymax 
) const
protected

Definition at line 626 of file ViewBase.cc.

628 {
629
630 if (!cmp) return false;
631 // Try to get the area/bounding box from the sensor/component.
632 std::array<double, 3> bbmin;
633 std::array<double, 3> bbmax;
634 if (!cmp->GetBoundingBox(bbmin[0], bbmin[1], bbmin[2],
635 bbmax[0], bbmax[1], bbmax[2])) {
636 std::cerr << m_className << "::PlotLimits:\n"
637 << " Bounding box of the component is not defined.\n"
638 << " Please set the plot limits explicitly (SetArea).\n";
639 return false;
640 }
641 if (std::isinf(bbmin[0]) || std::isinf(bbmax[0]) ||
642 std::isinf(bbmin[1]) || std::isinf(bbmax[1]) ||
643 std::isinf(bbmin[2]) || std::isinf(bbmax[2])) {
644 std::array<double, 3> cellmin = {0., 0., 0.};
645 std::array<double, 3> cellmax = {0., 0., 0.};
646 if (!cmp->GetElementaryCell(cellmin[0], cellmin[1], cellmin[2],
647 cellmax[0], cellmax[1], cellmax[2])) {
648 std::cerr << m_className << "::PlotLimits:\n"
649 << " Cell boundaries are not defined.\n"
650 << " Please set the plot limits explicitly (SetArea).\n";
651 return false;
652 }
653 for (size_t i = 0; i < 3; ++i) {
654 if (std::isinf(bbmin[i]) || std::isinf(bbmax[i])) {
655 bbmin[i] = cellmin[i];
656 bbmax[i] = cellmax[i];
657 }
658 }
659 }
660 return PlotLimits(bbmin, bbmax, xmin, ymin, xmax, ymax);
661}
bool PlotLimits(Sensor *sensor, double &xmin, double &ymin, double &xmax, double &ymax) const
Definition: ViewBase.cc:608

◆ PlotLimits() [2/3]

bool Garfield::ViewBase::PlotLimits ( Sensor sensor,
double &  xmin,
double &  ymin,
double &  xmax,
double &  ymax 
) const
protected

Definition at line 608 of file ViewBase.cc.

610 {
611
612 if (!sensor) return false;
613 // Try to get the area/bounding box from the sensor/component.
614 std::array<double, 3> bbmin;
615 std::array<double, 3> bbmax;
616 if (!sensor->GetArea(bbmin[0], bbmin[1], bbmin[2],
617 bbmax[0], bbmax[1], bbmax[2])) {
618 std::cerr << m_className << "::PlotLimits:\n"
619 << " Sensor area is not defined.\n"
620 << " Please set the plot limits explicitly (SetArea).\n";
621 return false;
622 }
623 return PlotLimits(bbmin, bbmax, xmin, ymin, xmax, ymax);
624}

Referenced by Garfield::ViewGeometry::Plot2d(), PlotLimits(), and PlotLimitsFromUserBox().

◆ PlotLimits() [3/3]

bool Garfield::ViewBase::PlotLimits ( std::array< double, 3 > &  bbmin,
std::array< double, 3 > &  bbmax,
double &  xmin,
double &  ymin,
double &  xmax,
double &  ymax 
) const
protected

Definition at line 671 of file ViewBase.cc.

674 {
675 constexpr double tol = 1.e-4;
676 double umin[2] = {-std::numeric_limits<double>::max(),
677 -std::numeric_limits<double>::max()};
678 double umax[2] = {std::numeric_limits<double>::max(),
679 std::numeric_limits<double>::max()};
680 for (unsigned int i = 0; i < 3; ++i) {
681 bbmin[i] -= m_proj[2][i];
682 bbmax[i] -= m_proj[2][i];
683 for (unsigned int j = 0; j < 2; ++j) {
684 if (fabs(m_proj[j][i]) < tol) continue;
685 const double t1 = bbmin[i] / m_proj[j][i];
686 const double t2 = bbmax[i] / m_proj[j][i];
687 const double tmin = std::min(t1, t2);
688 const double tmax = std::max(t1, t2);
689 if (tmin > umin[j] && tmin < umax[j]) umin[j] = tmin;
690 if (tmax < umax[j] && tmax > umin[j]) umax[j] = tmax;
691 }
692 }
693 xmin = umin[0];
694 xmax = umax[0];
695 ymin = umin[1];
696 ymax = umax[1];
697 return true;
698}

◆ PlotLimitsFromUserBox()

bool Garfield::ViewBase::PlotLimitsFromUserBox ( double &  xmin,
double &  ymin,
double &  xmax,
double &  ymax 
) const
protected

Definition at line 663 of file ViewBase.cc.

664 {
665
666 std::array<double, 3> bbmin = {m_xMinBox, m_yMinBox, m_zMinBox};
667 std::array<double, 3> bbmax = {m_xMaxBox, m_yMaxBox, m_zMaxBox};
668 return PlotLimits(bbmin, bbmax, xmin, ymin, xmax, ymax);
669}

Referenced by Garfield::ViewGeometry::Plot2d().

◆ RangeSet()

bool Garfield::ViewBase::RangeSet ( TPad *  pad)
protected

Definition at line 83 of file ViewBase.cc.

83 {
84
85 if (pad->GetListOfPrimitives()->GetSize() == 0 &&
86 pad->GetX1() == 0 && pad->GetX2() == 1 &&
87 pad->GetY1() == 0 && pad->GetY2() == 1) {
88 return false;
89 }
90 return true;
91}

Referenced by Garfield::ViewFEMesh::Plot(), Garfield::ViewDrift::Plot2d(), and Garfield::ViewField::PlotFieldLines().

◆ Rotate()

void Garfield::ViewBase::Rotate ( const double  angle)

Rotate the viewing plane (angle in radian).

Definition at line 255 of file ViewBase.cc.

255 {
256 // Rotate the axes
257 double auxu[3], auxv[3];
258 const double ctheta = cos(theta);
259 const double stheta = sin(theta);
260 for (int i = 0; i < 3; ++i) {
261 auxu[i] = ctheta * m_proj[0][i] - stheta * m_proj[1][i];
262 auxv[i] = stheta * m_proj[0][i] + ctheta * m_proj[1][i];
263 }
264 for (int i = 0; i < 3; ++i) {
265 m_proj[0][i] = auxu[i];
266 m_proj[1][i] = auxv[i];
267 }
268
270}
void UpdateProjectionMatrix()
Definition: ViewBase.cc:320
DoubleAc cos(const DoubleAc &f)
Definition: DoubleAc.cpp:432
DoubleAc sin(const DoubleAc &f)
Definition: DoubleAc.cpp:384

◆ SetArea() [1/3]

void Garfield::ViewBase::SetArea ( )
inline

Use default x- and y-axis limits (based on the bounding box of the sensor/component, if applicable).

Definition at line 43 of file ViewBase.hh.

43 {
44 m_userBox = false;
45 m_userPlotLimits = false;
46 }

◆ SetArea() [2/3]

void Garfield::ViewBase::SetArea ( const double  xmin,
const double  ymin,
const double  xmax,
const double  ymax 
)

Set the x- and y-axis limits (in local coordinates of the current viewing plane, if applicable).

Definition at line 109 of file ViewBase.cc.

110 {
111 // Check range, assign if non-null.
112 if (xmin == xmax || ymin == ymax) {
113 std::cerr << m_className << "::SetArea: Null area is not permitted.\n"
114 << " " << xmin << " < x < " << xmax << "\n"
115 << " " << ymin << " < y < " << ymax << "\n";
116 return;
117 }
118 m_xMinPlot = std::min(xmin, xmax);
119 m_yMinPlot = std::min(ymin, ymax);
120 m_xMaxPlot = std::max(xmin, xmax);
121 m_yMaxPlot = std::max(ymin, ymax);
122 m_userPlotLimits = true;
123}

Referenced by main().

◆ SetArea() [3/3]

void Garfield::ViewBase::SetArea ( const double  xmin,
const double  ymin,
const double  zmin,
const double  xmax,
const double  ymax,
const double  zmax 
)
virtual

Set a bounding box (if applicable).

Definition at line 125 of file ViewBase.cc.

127 {
128 // Check range, assign if non-null
129 if (xmin == xmax || ymin == ymax || zmin == zmax) {
130 std::cerr << m_className << "::SetArea: Null area range not permitted.\n";
131 return;
132 }
133 m_xMinBox = std::min(xmin, xmax);
134 m_yMinBox = std::min(ymin, ymax);
135 m_zMinBox = std::min(zmin, zmax);
136 m_xMaxBox = std::max(xmin, xmax);
137 m_yMaxBox = std::max(ymin, ymax);
138 m_zMaxBox = std::max(zmin, zmax);
139 m_userBox = true;
140}

◆ SetCanvas() [1/2]

void Garfield::ViewBase::SetCanvas ( )
inline

Unset an external canvas.

Definition at line 30 of file ViewBase.hh.

30{ m_pad = nullptr; }

◆ SetCanvas() [2/2]

void Garfield::ViewBase::SetCanvas ( TPad *  pad)
inline

Set the canvas to be painted on.

Definition at line 28 of file ViewBase.hh.

28{ m_pad = pad; }

Referenced by main().

◆ SetPlane() [1/2]

void Garfield::ViewBase::SetPlane ( const double  fx,
const double  fy,
const double  fz,
const double  x0,
const double  y0,
const double  z0 
)
virtual

Set the projection (viewing plane), if applicable.

Parameters
fx,fy,fznormal vector
x0,y0,z0in-plane point

Reimplemented in Garfield::ViewFEMesh.

Definition at line 142 of file ViewBase.cc.

143 {
144 // Calculate two in-plane vectors for the normal vector
145 const double fnorm = sqrt(fx * fx + fy * fy + fz * fz);
146 if (fnorm > 0 && fx * fx + fz * fz > 0) {
147 const double fxz = sqrt(fx * fx + fz * fz);
148 m_proj[0][0] = fz / fxz;
149 m_proj[0][1] = 0;
150 m_proj[0][2] = -fx / fxz;
151 m_proj[1][0] = -fx * fy / (fxz * fnorm);
152 m_proj[1][1] = (fx * fx + fz * fz) / (fxz * fnorm);
153 m_proj[1][2] = -fy * fz / (fxz * fnorm);
154 m_proj[2][0] = x0;
155 m_proj[2][1] = y0;
156 m_proj[2][2] = z0;
157 } else if (fnorm > 0 && fy * fy + fz * fz > 0) {
158 const double fyz = sqrt(fy * fy + fz * fz);
159 m_proj[0][0] = (fy * fy + fz * fz) / (fyz * fnorm);
160 m_proj[0][1] = -fx * fz / (fyz * fnorm);
161 m_proj[0][2] = -fy * fz / (fyz * fnorm);
162 m_proj[1][0] = 0;
163 m_proj[1][1] = fz / fyz;
164 m_proj[1][2] = -fy / fyz;
165 m_proj[2][0] = x0;
166 m_proj[2][1] = y0;
167 m_proj[2][2] = z0;
168 } else {
169 std::cout << m_className << "::SetPlane:\n"
170 << " Normal vector has zero norm. No new projection set.\n";
171 }
172
173 // Store the plane description
174 m_plane[0] = fx;
175 m_plane[1] = fy;
176 m_plane[2] = fz;
177 m_plane[3] = fx * x0 + fy * y0 + fz * z0;
178
180}
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:314

Referenced by main(), and Garfield::ViewFEMesh::SetPlane().

◆ SetPlane() [2/2]

void Garfield::ViewBase::SetPlane ( const double  fx,
const double  fy,
const double  fz,
const double  x0,
const double  y0,
const double  z0,
const double  hx,
const double  hy,
const double  hz 
)
virtual

Set the projection plane specifying a hint for the in-plane x axis.

Reimplemented in Garfield::ViewFEMesh.

Definition at line 182 of file ViewBase.cc.

184 {
185
186 const double fnorm = sqrt(fx * fx + fy * fy + fz * fz);
187 if (fnorm < Small) {
188 std::cout << m_className << "::SetPlane:\n"
189 << " Normal vector has zero norm. No new projection set.\n";
190 return;
191 }
192 // Normalise the vector.
193 const double wx = fx / fnorm;
194 const double wy = fy / fnorm;
195 const double wz = fz / fnorm;
196 // Store the plane description.
197 m_plane[0] = wx;
198 m_plane[1] = wy;
199 m_plane[2] = wz;
200 m_plane[3] = wx * x0 + wy * y0 + wz * z0;
201
202 double d = hx * wx + hy * wy + hz * wz;
203 double ux = hx - d * wx;
204 double uy = hy - d * wy;
205 double uz = hz - d * wz;
206 double unorm = std::sqrt(ux * ux + uy * uy + uz * uz);
207 if (unorm < 1.e-10) {
208 // Wrong in-plane x hint (close to norm).
209 if (fy * fy + fz * fz > 0) {
210 // Taking global x as in-plane x hint.
211 ux = 1;
212 uy = 0;
213 uz = 0;
214 } else {
215 // Taking global y as in-plane x hint.
216 ux = 0;
217 uy = 1;
218 uz = 0;
219 }
220 d = ux * wx + uy * wy + uz * wz;
221 ux -= d * wx;
222 uy -= d * wy;
223 uz -= d * wz;
224 unorm = std::sqrt(ux * ux + uy * uy + uz * uz);
225 }
226 ux /= unorm;
227 uy /= unorm;
228 uz /= unorm;
229
230 m_prmat[0][0] = ux;
231 m_prmat[1][0] = uy;
232 m_prmat[2][0] = uz;
233 // In-plane y = cross product [z,x]
234 m_prmat[0][1] = wy * uz - wz * uy;
235 m_prmat[1][1] = wz * ux - wx * uz;
236 m_prmat[2][1] = wx * uy - wy * ux;
237 m_prmat[0][2] = wx;
238 m_prmat[1][2] = wy;
239 m_prmat[2][2] = wz;
240
241 for (unsigned int i = 0; i < 3; ++i) {
242 m_proj[0][i] = m_prmat[i][0];
243 m_proj[1][i] = m_prmat[i][1];
244 }
245 m_proj[2][0] = x0;
246 m_proj[2][1] = y0;
247 m_proj[2][2] = z0;
248 if (!Invert(m_prmat)) {
249 std::cerr << m_className << "::SetPlane:\n"
250 << " Inversion failed; reset to default.\n";
251 SetPlaneXY();
252 }
253}
std::array< std::array< double, 3 >, 3 > m_prmat
Definition: ViewBase.hh:101
void SetPlaneXY()
Set the viewing plane to x-y.
Definition: ViewBase.cc:272

◆ SetPlaneXY()

void Garfield::ViewBase::SetPlaneXY ( )

Set the viewing plane to x-y.

Definition at line 272 of file ViewBase.cc.

272 {
273 m_proj = {{{1, 0, 0}, {0, 1, 0}, {0, 0, 0}}};
274 m_plane = {0, 0, 1, 0};
275 m_prmat = {{{1, 0, 0}, {0, 1, 0}, {0, 0, 1}}};
276}

Referenced by SetPlane(), and UpdateProjectionMatrix().

◆ SetPlaneXZ()

void Garfield::ViewBase::SetPlaneXZ ( )

Set the viewing plane to x-z.

Definition at line 278 of file ViewBase.cc.

278 {
279 m_proj = {{{1, 0, 0}, {0, 0, 1}, {0, 0, 0}}};
280 m_plane = {0, 1, 0, 0};
282}

Referenced by main().

◆ SetPlaneYZ()

void Garfield::ViewBase::SetPlaneYZ ( )

Set the viewing plane to y-z.

Definition at line 284 of file ViewBase.cc.

284 {
285 m_proj = {{{0, 1, 0}, {0, 0, 1}, {0, 0, 0}}};
286 m_plane = {1, 0, 0, 0};
288}

◆ SetRange()

void Garfield::ViewBase::SetRange ( TPad *  pad,
const double  x0,
const double  y0,
const double  x1,
const double  y1 
)
protected

Definition at line 93 of file ViewBase.cc.

94 {
95 if (!pad) return;
96 const double bm = pad->GetBottomMargin();
97 const double lm = pad->GetLeftMargin();
98 const double rm = pad->GetRightMargin();
99 const double tm = pad->GetTopMargin();
100 const double dx = x1 - x0;
101 const double dy = y1 - y0;
102 pad->Range(x0 - dx * (lm / (1. - rm - lm)),
103 y0 - dy * (bm / (1. - tm - lm)),
104 x1 + dx * (rm / (1. - rm - lm)),
105 y1 + dy * (tm / (1. - tm - lm)));
106
107}

Referenced by Garfield::ViewFEMesh::Plot(), Garfield::ViewDrift::Plot2d(), and Garfield::ViewField::PlotFieldLines().

◆ ToPlane()

template<typename T >
void Garfield::ViewBase::ToPlane ( const T  x,
const T  y,
const T  z,
T &  xp,
T &  yp 
) const
inlineprotected

Definition at line 109 of file ViewBase.hh.

109 {
110 xp = m_prmat[0][0] * x + m_prmat[0][1] * y + m_prmat[0][2] * z;
111 yp = m_prmat[1][0] * x + m_prmat[1][1] * y + m_prmat[1][2] * z;
112 }

Referenced by DrawLine(), Garfield::ViewFEMesh::Plot(), Garfield::ViewGeometry::Plot2d(), and Garfield::ViewDrift::Plot2d().

◆ UpdateProjectionMatrix()

void Garfield::ViewBase::UpdateProjectionMatrix ( )
protected

Definition at line 320 of file ViewBase.cc.

320 {
321
322 m_prmat[0][0] = m_proj[0][0];
323 m_prmat[1][0] = m_proj[0][1];
324 m_prmat[2][0] = m_proj[0][2];
325 m_prmat[0][1] = m_proj[1][0];
326 m_prmat[1][1] = m_proj[1][1];
327 m_prmat[2][1] = m_proj[1][2];
328 const double vnorm = sqrt(m_plane[0] * m_plane[0] +
329 m_plane[1] * m_plane[1] +
330 m_plane[2] * m_plane[2]);
331 if (vnorm <= 0.) {
332 std::cerr << m_className << "::UpdateProjectionMatrix:\n"
333 << " Zero norm vector; reset to default.\n";
334 SetPlaneXY();
335 return;
336 }
337 m_prmat[0][2] = m_plane[0] / vnorm;
338 m_prmat[1][2] = m_plane[1] / vnorm;
339 m_prmat[2][2] = m_plane[2] / vnorm;
340 if (!Invert(m_prmat)) {
341 std::cerr << m_className << "::UpdateProjectionMatrix:\n"
342 << " Inversion failed; reset to default.\n";
343 SetPlaneXY();
344 }
345}

Referenced by Rotate(), SetPlane(), SetPlaneXZ(), and SetPlaneYZ().

Member Data Documentation

◆ m_className

std::string Garfield::ViewBase::m_className = "ViewBase"
protected

Definition at line 78 of file ViewBase.hh.

Referenced by Garfield::ViewDrift::AddDriftLinePoint(), Garfield::ViewDrift::AddTrackPoint(), Garfield::ViewFEMesh::CreateDefaultAxes(), Garfield::ViewField::EqualFluxIntervals(), Garfield::ViewField::FixedFluxIntervals(), GetCanvas(), Garfield::ViewFEMesh::Plot(), Garfield::ViewSignal::Plot(), Garfield::ViewCell::Plot2d(), Garfield::ViewGeometry::Plot2d(), Garfield::ViewDrift::Plot2d(), Garfield::ViewCell::Plot3d(), Garfield::ViewGeometry::Plot3d(), Garfield::ViewDrift::Plot3d(), Garfield::ViewMedium::PlotElectronCrossSections(), Garfield::ViewField::PlotFieldLines(), Garfield::ViewIsochrons::PlotIsochrons(), PlotLimits(), Garfield::ViewSignal::PlotSignal(), SetArea(), Garfield::ViewIsochrons::SetAspectRatioSwitch(), Garfield::ViewDrift::SetClusterMarkerSize(), Garfield::ViewDrift::SetCollisionMarkerSize(), Garfield::ViewField::SetComponent(), Garfield::ViewIsochrons::SetComponent(), Garfield::ViewCell::SetComponent(), Garfield::ViewFEMesh::SetComponent(), Garfield::ViewIsochrons::SetConnectionThreshold(), Garfield::ViewDrift::SetDriftLinePoint(), Garfield::ViewGeometry::SetGeometry(), Garfield::ViewIsochrons::SetLoopThreshold(), Garfield::ViewMedium::SetMedium(), SetPlane(), Garfield::ViewMedium::SetRangeA(), Garfield::ViewMedium::SetRangeB(), Garfield::ViewMedium::SetRangeE(), Garfield::ViewSignal::SetRangeX(), Garfield::ViewSignal::SetRangeY(), Garfield::ViewMedium::SetRangeY(), Garfield::ViewField::SetSensor(), Garfield::ViewIsochrons::SetSensor(), Garfield::ViewSignal::SetSensor(), Garfield::ViewDrift::SetTrackPoint(), and UpdateProjectionMatrix().

◆ m_debug

bool Garfield::ViewBase::m_debug = false
protected

◆ m_plane

◆ m_prmat

std::array<std::array<double, 3>, 3> Garfield::ViewBase::m_prmat
protected
Initial value:
{{
{{1, 0, 0}}, {{0, 1, 0}}, {{0, 0, 1}}
}}

Definition at line 101 of file ViewBase.hh.

Referenced by SetPlane(), SetPlaneXY(), ToPlane(), and UpdateProjectionMatrix().

◆ m_proj

std::array<std::array<double, 3>, 3> Garfield::ViewBase::m_proj
protected
Initial value:
{{
{{1, 0, 0}}, {{0, 1, 0}}, {{0, 0, 0}}
}}

Definition at line 96 of file ViewBase.hh.

Referenced by LabelX(), LabelY(), Garfield::ViewGeometry::Plot2d(), Garfield::ViewIsochrons::PlotIsochrons(), PlotLimits(), Rotate(), SetPlane(), SetPlaneXY(), SetPlaneXZ(), SetPlaneYZ(), and UpdateProjectionMatrix().

◆ m_userBox

bool Garfield::ViewBase::m_userBox = false
protected

◆ m_userPlotLimits

bool Garfield::ViewBase::m_userPlotLimits = false
protected

Definition at line 84 of file ViewBase.hh.

Referenced by Garfield::ViewGeometry::Plot2d(), and SetArea().

◆ m_xMaxBox

double Garfield::ViewBase::m_xMaxBox = 1.
protected

◆ m_xMaxPlot

◆ m_xMinBox

double Garfield::ViewBase::m_xMinBox = -1.
protected

◆ m_xMinPlot

◆ m_yMaxBox

double Garfield::ViewBase::m_yMaxBox = 1.
protected

◆ m_yMaxPlot

◆ m_yMinBox

double Garfield::ViewBase::m_yMinBox = -1.
protected

◆ m_yMinPlot

◆ m_zMaxBox

double Garfield::ViewBase::m_zMaxBox = 1.
protected

◆ m_zMinBox

double Garfield::ViewBase::m_zMinBox = -1.
protected

The documentation for this class was generated from the following files: