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::ViewFEMesh Class Reference

Draw the mesh of a field-map component. More...

#include <ViewFEMesh.hh>

+ Inheritance diagram for Garfield::ViewFEMesh:

Public Member Functions

 ViewFEMesh ()
 Constructor.
 
 ~ViewFEMesh ()=default
 Destructor.
 
void SetComponent (ComponentFieldMap *cmp)
 Set the component from which to retrieve the mesh and field.
 
void SetPlane (const double fx, const double fy, const double fz, const double x0, const double y0, const double z0) override
 
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) override
 Set the projection plane specifying a hint for the in-plane x axis.
 
void SetXaxis (TGaxis *ax)
 
void SetYaxis (TGaxis *ay)
 
void SetXaxisTitle (const std::string &xtitle)
 
void SetYaxisTitle (const std::string &ytitle)
 
void EnableAxes ()
 
void DisableAxes ()
 
bool Plot ()
 Plot method to be called by user.
 
void SetFillMesh (const bool f)
 Element fill switch; 2D only, set false for wireframe mesh.
 
void SetDrawViewRegion (bool do_draw)
 Display intersection of projection plane with viewing area.
 
bool GetDrawViewRegion (void) const
 
void SetColor (int matID, int colorID)
 
void SetFillColor (int matID, int colorID)
 
void SetViewDrift (ViewDrift *vd)
 Set the optional associated ViewDrift.
 
void SetFillMeshWithBorders ()
 Show filled mesh elements.
 
void CreateDefaultAxes ()
 Create a default set of custom-made axes.
 
void DisableMaterial (int materialID)
 Disable a material so that its mesh cells are not drawn.
 
- Public Member Functions inherited from Garfield::ViewBase
 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.
 

Additional Inherited Members

- Static Public Member Functions inherited from Garfield::ViewBase
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 inherited from Garfield::ViewBase
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 inherited from Garfield::ViewBase
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

Draw the mesh of a field-map component.

Definition at line 21 of file ViewFEMesh.hh.

Constructor & Destructor Documentation

◆ ViewFEMesh()

Garfield::ViewFEMesh::ViewFEMesh ( )

Constructor.

Definition at line 19 of file ViewFEMesh.cc.

19: ViewBase("ViewFEMesh") {}
ViewBase()=delete
Default constructor.

◆ ~ViewFEMesh()

Garfield::ViewFEMesh::~ViewFEMesh ( )
default

Destructor.

Member Function Documentation

◆ CreateDefaultAxes()

void Garfield::ViewFEMesh::CreateDefaultAxes ( )

Create a default set of custom-made axes.

Definition at line 224 of file ViewFEMesh.cc.

224 {
225 // Create a new x and y axis.
226 if (!GetPlotLimits()) {
227 std::cerr << m_className << "::CreateDefaultAxes:\n"
228 << " Cannot determine the axis limits.\n";
229 return;
230 }
231 const double dx = std::abs(m_xMaxPlot - m_xMinPlot) * 0.1;
232 const double dy = std::abs(m_yMaxPlot - m_yMinPlot) * 0.1;
233 const double x0 = m_xMinPlot + dx;
234 const double y0 = m_yMinPlot + dy;
235 const double x1 = m_xMaxPlot - dx;
236 const double y1 = m_yMaxPlot - dy;
237 m_xaxis = new TGaxis(x0, y0, x1, y0, x0, x1, 2405, "x");
238 m_yaxis = new TGaxis(x0, y0, x0, y1, y0, y1, 2405, "y");
239
240 // Label sizes
241 m_xaxis->SetLabelSize(0.025);
242 m_yaxis->SetLabelSize(0.025);
243
244 // Titles
245 m_xaxis->SetTitleSize(0.03);
246 m_xaxis->SetTitle(LabelX().c_str());
247 m_yaxis->SetTitleSize(0.03);
248 m_yaxis->SetTitle(LabelY().c_str());
249}
std::string LabelY()
Definition: ViewBase.cc:483
std::string LabelX()
Definition: ViewBase.cc:420
std::string m_className
Definition: ViewBase.hh:78

◆ DisableAxes()

void Garfield::ViewFEMesh::DisableAxes ( )
inline

Definition at line 43 of file ViewFEMesh.hh.

43{ m_drawAxes = false; }

◆ DisableMaterial()

void Garfield::ViewFEMesh::DisableMaterial ( int  materialID)
inline

Disable a material so that its mesh cells are not drawn.

Definition at line 75 of file ViewFEMesh.hh.

75 {
76 m_disabledMaterial[materialID] = true;
77 }

◆ EnableAxes()

void Garfield::ViewFEMesh::EnableAxes ( )
inline

Definition at line 42 of file ViewFEMesh.hh.

42{ m_drawAxes = true; }

Referenced by main().

◆ GetDrawViewRegion()

bool Garfield::ViewFEMesh::GetDrawViewRegion ( void  ) const
inline

Definition at line 53 of file ViewFEMesh.hh.

53{ return m_drawViewRegion; }

◆ Plot()

bool Garfield::ViewFEMesh::Plot ( )

Plot method to be called by user.

Definition at line 32 of file ViewFEMesh.cc.

32 {
33 if (!m_component) {
34 std::cerr << m_className << "::Plot: Component is not defined.\n";
35 return false;
36 }
37
38 double pmin = 0., pmax = 0.;
39 if (!m_component->GetVoltageRange(pmin, pmax)) {
40 std::cerr << m_className << "::Plot: Component is not ready.\n";
41 return false;
42 }
43
44 if (!GetPlotLimits()) return false;
45
46 auto pad = GetCanvas();
47 pad->cd();
48 if (!RangeSet(pad)) {
50 }
51
52 if (m_drawAxes) {
53 if (!m_xaxis && !m_yaxis) {
54 // Draw default axes.
55 auto frame = pad->DrawFrame(m_xMinPlot, m_yMinPlot,
57 if (m_xaxisTitle.empty()) {
58 frame->GetXaxis()->SetTitle(LabelX().c_str());
59 } else {
60 frame->GetXaxis()->SetTitle(m_xaxisTitle.c_str());
61 }
62 if (m_yaxisTitle.empty()) {
63 frame->GetYaxis()->SetTitle(LabelY().c_str());
64 } else {
65 frame->GetYaxis()->SetTitle(m_yaxisTitle.c_str());
66 }
67 } else {
68 // Draw custom axes.
69 if (m_xaxis) m_xaxis->Draw();
70 if (m_yaxis) m_yaxis->Draw();
71 }
72 }
73
74 // Plot the mesh elements.
75 ComponentCST* componentCST = dynamic_cast<ComponentCST*>(m_component);
76 if (componentCST) {
77 std::cout << m_className << "::Plot: CST component. Calling DrawCST.\n";
78 DrawCST(componentCST);
79 } else {
80 DrawElements();
81 }
82
83 // If we have an associated ViewDrift object, plot the drift lines.
84 if (m_viewDrift) {
85 // Plot a 2D projection of the drift line.
86 for (const auto& dline : m_viewDrift->m_driftLines) {
87 TGraph gr;
88 if (dline.second == ViewDrift::Particle::Electron) {
89 gr.SetLineColor(kOrange - 3);
90 } else {
91 gr.SetLineColor(kRed + 1);
92 }
93 std::vector<float> xgr;
94 std::vector<float> ygr;
95 // Loop over the points.
96 for (const auto& point : dline.first) {
97 // Project this point onto the plane.
98 float xp = 0., yp = 0.;
99 ToPlane(point[0], point[1], point[2], xp, yp);
100 // Add this point if it is within the view.
101 if (InView(xp, yp)) {
102 xgr.push_back(xp);
103 ygr.push_back(yp);
104 }
105 }
106 if (!xgr.empty()) {
107 gr.DrawGraph(xgr.size(), xgr.data(), ygr.data(), "lsame");
108 }
109 }
110 }
111
112 if (m_drawViewRegion && !m_viewRegionX.empty()) {
113 TPolyLine poly;
114 poly.SetLineColor(kSpring + 4);
115 poly.SetLineWidth(3);
116 std::vector<double> xv = m_viewRegionX;
117 std::vector<double> yv = m_viewRegionY;
118 // Close the polygon.
119 xv.push_back(m_viewRegionX[0]);
120 yv.push_back(m_viewRegionY[0]);
121 poly.DrawPolyLine(xv.size(), xv.data(), yv.data(), "same");
122 }
123 gPad->Update();
124 // Draw axes again so they are on top
125 gPad->RedrawAxis("g");
126 return true;
127}
bool GetVoltageRange(double &vmin, double &vmax) override
Calculate the voltage range [V].
bool RangeSet(TPad *)
Definition: ViewBase.cc:83
TPad * GetCanvas()
Retrieve the canvas.
Definition: ViewBase.cc:74
void SetRange(TPad *pad, const double x0, const double y0, const double x1, const double y1)
Definition: ViewBase.cc:93
void ToPlane(const T x, const T y, const T z, T &xp, T &yp) const
Definition: ViewBase.hh:109

Referenced by main().

◆ SetColor()

void Garfield::ViewFEMesh::SetColor ( int  matID,
int  colorID 
)
inline

Associate a color with each element material map ID; Uses ROOT color numberings

Definition at line 57 of file ViewFEMesh.hh.

57{ m_colorMap[matID] = colorID; }

Referenced by main().

◆ SetComponent()

void Garfield::ViewFEMesh::SetComponent ( ComponentFieldMap cmp)

Set the component from which to retrieve the mesh and field.

Definition at line 21 of file ViewFEMesh.cc.

21 {
22 if (!cmp) {
23 std::cerr << m_className << "::SetComponent: Null pointer.\n";
24 return;
25 }
26
27 m_component = cmp;
28}

Referenced by main().

◆ SetDrawViewRegion()

void Garfield::ViewFEMesh::SetDrawViewRegion ( bool  do_draw)
inline

Display intersection of projection plane with viewing area.

Definition at line 52 of file ViewFEMesh.hh.

52{ m_drawViewRegion = do_draw; }

◆ SetFillColor()

void Garfield::ViewFEMesh::SetFillColor ( int  matID,
int  colorID 
)
inline

Definition at line 58 of file ViewFEMesh.hh.

58 {
59 m_colorMap_fill[matID] = colorID;
60 }

◆ SetFillMesh()

void Garfield::ViewFEMesh::SetFillMesh ( const bool  f)
inline

Element fill switch; 2D only, set false for wireframe mesh.

Definition at line 49 of file ViewFEMesh.hh.

49{ m_fillMesh = f; }

Referenced by main().

◆ SetFillMeshWithBorders()

void Garfield::ViewFEMesh::SetFillMeshWithBorders ( )
inline

Show filled mesh elements.

Definition at line 66 of file ViewFEMesh.hh.

66 {
67 m_plotMeshBorders = true;
68 m_fillMesh = true;
69 }

◆ SetPlane() [1/2]

void Garfield::ViewFEMesh::SetPlane ( const double  fx,
const double  fy,
const double  fz,
const double  x0,
const double  y0,
const double  z0 
)
overridevirtual

Set the projection (viewing plane), if applicable.

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

Reimplemented from Garfield::ViewBase.

Definition at line 202 of file ViewFEMesh.cc.

203 {
204 if (fy * fy + fz * fz > 0) {
205 SetPlane(fx, fy, fz, x0, y0, z0, 1, 0, 0);
206 } else {
207 SetPlane(fx, fy, fz, x0, y0, z0, 0, 1, 0);
208 }
209}
void SetPlane(const double fx, const double fy, const double fz, const double x0, const double y0, const double z0) override
Definition: ViewFEMesh.cc:202

Referenced by main(), and SetPlane().

◆ SetPlane() [2/2]

void Garfield::ViewFEMesh::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 
)
overridevirtual

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

Reimplemented from Garfield::ViewBase.

Definition at line 211 of file ViewFEMesh.cc.

213 {
214 ViewBase::SetPlane(fx, fy, fz, x0, y0, z0, hx, hy, hz);
215}
virtual void SetPlane(const double fx, const double fy, const double fz, const double x0, const double y0, const double z0)
Definition: ViewBase.cc:142

◆ SetViewDrift()

void Garfield::ViewFEMesh::SetViewDrift ( ViewDrift vd)
inline

Set the optional associated ViewDrift.

Definition at line 63 of file ViewFEMesh.hh.

63{ m_viewDrift = vd; }

Referenced by main().

◆ SetXaxis()

void Garfield::ViewFEMesh::SetXaxis ( TGaxis *  ax)

Definition at line 218 of file ViewFEMesh.cc.

218{ m_xaxis = ax; }

◆ SetXaxisTitle()

void Garfield::ViewFEMesh::SetXaxisTitle ( const std::string &  xtitle)
inline

Definition at line 40 of file ViewFEMesh.hh.

40{ m_xaxisTitle = xtitle; }

Referenced by main().

◆ SetYaxis()

void Garfield::ViewFEMesh::SetYaxis ( TGaxis *  ay)

Definition at line 221 of file ViewFEMesh.cc.

221{ m_yaxis = ay; }

◆ SetYaxisTitle()

void Garfield::ViewFEMesh::SetYaxisTitle ( const std::string &  ytitle)
inline

Definition at line 41 of file ViewFEMesh.hh.

41{ m_yaxisTitle = ytitle; }

Referenced by main().


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