Garfield++ 3.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
ViewFEMesh.hh
Go to the documentation of this file.
1#ifndef G_VIEW_FE_MESH
2#define G_VIEW_FE_MESH
3
4#include <memory>
5#include <string>
6#ifndef __CINT__
7#include <map>
8#endif
9
10#include <TArrayD.h>
11#include <TGaxis.h>
12#include <TH2D.h>
13#include <TMatrixD.h>
14#include <TPolyLine.h>
15#include <TPolyLine3D.h>
16
17#include "ComponentCST.hh"
18#include "ComponentFieldMap.hh"
19#include "ViewBase.hh"
20#include "ViewDrift.hh"
21
22namespace Garfield {
23
24/// Draw the mesh of a field-map component.
25
26class ViewFEMesh : public ViewBase {
27 public:
28 /// Constructor.
29 ViewFEMesh();
30 /// Destructor.
31 ~ViewFEMesh() = default;
32
33 /// Set the component from which to retrieve the mesh and field.
35
36 /// Set area to be plotted to the default.
37 void SetArea();
38 /// Set area to be plotted explicitly.
39 void SetArea(const double xmin, const double ymin, const double zmin,
40 const double xmax, const double ymax, const double zmax);
41 /// Reset the projection plane.
43 /// Set the projection plane.
44 void SetPlane(double fx, double fy, double fz, double x0, double y0,
45 double z0);
46 /// Set the projection plane specifying hint for in-plane x axis.
47 void SetPlane(double fx, double fy, double fz, double x0, double y0,
48 double z0, double hx, double hy, double hz);
49
50 // Axes
51 void SetXaxis(TGaxis* ax);
52 void SetYaxis(TGaxis* ay);
53 void SetXaxisTitle(const char* xtitle);
54 void SetYaxisTitle(const char* ytitle);
55 void EnableAxes() { m_drawAxes = true; }
56 void DisableAxes() { m_drawAxes = false; }
57
58 /// Plot method to be called by user
59 bool Plot();
60
61 /// Element fill switch; 2D only, set false for wireframe mesh
62 void SetFillMesh(const bool f) { m_fillMesh = f; }
63
64 /// Display intersection of projection plane with viewing area
65 void SetDrawViewRegion(bool do_draw) { m_drawViewRegion = do_draw; }
66 bool GetDrawViewRegion(void) const { return m_drawViewRegion; }
67
68 /// Associate a color with each element material map ID;
69 /// Uses ROOT color numberings
70 void SetColor(int matID, int colorID) { m_colorMap[matID] = colorID; }
71 void SetFillColor(int matID, int colorID) {
72 m_colorMap_fill[matID] = colorID;
73 }
74
75 /// Set the optional associated ViewDrift
76 void SetViewDrift(ViewDrift* vd) { m_viewDrift = vd; }
77
78 /// Show filled mesh elements
80 m_plotMeshBorders = true;
81 m_fillMesh = true;
82 }
83
84 /// Create a default set of custom-made axes.
85 void CreateDefaultAxes();
86
87 /// Disable a material so that its mesh cells are not drawn
88 void DisableMaterial(int materialID) {
89 m_disabledMaterial[materialID] = true;
90 }
91
92 private:
93 std::string m_label = "Mesh";
94
95 // Options
96 bool m_fillMesh = false;
97
98 // Viewing plane (plane normal is stored in m_proj[2])
99 double m_proj[3][3];
100 double m_dist;
101
102 // Box dimensions
103 bool m_hasUserArea = false;
104 double m_xMin = -1., m_yMin = -1., m_zMin = -1.;
105 double m_xMax = 1., m_yMax = 1., m_zMax = 1.;
106
107 // Intersection of viewing plane with plotted area in planar coordinates
108 bool m_drawViewRegion = false;
109 std::vector<TPolyLine> m_viewRegionLines;
110
111 // Viewing plane dimensions
112 double m_xPlaneMin = -1., m_xPlaneMax = 1.;
113 double m_yPlaneMin = -1., m_yPlaneMax = 1.;
114
115 // The field map object
116 ComponentFieldMap* m_component = nullptr;
117
118 // Optional associated ViewDrift object
119 ViewDrift* m_viewDrift = nullptr;
120 bool m_plotMeshBorders = false;
121
122 // Axes
123 TGaxis* m_xaxis = nullptr;
124 TGaxis* m_yaxis = nullptr;
125 std::unique_ptr<TH2D> m_axes;
126 bool m_drawAxes = false;
127
128 // The mesh, stored as a vector of TPolyLine(3D) objects
129 std::vector<TPolyLine> m_mesh;
130 std::vector<TPolyLine> m_driftLines;
131
132// The color map
133#ifndef __CINT__
134 std::map<int, int> m_colorMap;
135 std::map<int, int> m_colorMap_fill;
136
137 // Disabled materials -> not shown in the mesh view
138 std::map<int, bool> m_disabledMaterial;
139#endif
140 // Element plotting methods
141 void DrawElements();
142 void DrawCST(ComponentCST* componentCST);
143
144 /// Return true if the specified point is in the view region.
145 bool InView(const double x, const double y) const;
146
147 bool LinesCrossed(double x1, double y1, double x2, double y2, double u1,
148 double v1, double u2, double v2, double& xc,
149 double& yc) const;
150 std::string CreateAxisTitle(const double* norm) const;
151 bool IntersectPlaneArea(void);
152 bool PlaneVector(double& x, double& y, double& z) const;
153 bool OnLine(double x1, double y1, double x2, double y2, double u,
154 double v) const;
155 void RemoveCrossings(std::vector<double>& x, std::vector<double>& y);
156 bool PlaneCut(double x1, double y1, double z1, double x2, double y2,
157 double z2, TMatrixD& xMat);
158 bool PlaneCoords(double x, double y, double z, const TMatrixD& projMat,
159 TMatrixD& xMat);
160 void ClipToView(std::vector<double>& px, std::vector<double>& py,
161 std::vector<double>& cx, std::vector<double>& cy);
162 bool IsInPolygon(double x, double y, std::vector<double>& px,
163 std::vector<double>& py, bool& edge) const;
164
165 // Plot method to be called by Plot() for CST cubic elements
166 // available are "xy", "yz" and "xz"
167};
168} // namespace Garfield
169#endif
Base class for components based on finite-element field maps.
Base class for visualization classes.
Definition: ViewBase.hh:10
Visualize drift lines and tracks.
Definition: ViewDrift.hh:19
Draw the mesh of a field-map component.
Definition: ViewFEMesh.hh:26
void SetXaxisTitle(const char *xtitle)
Definition: ViewFEMesh.cc:184
bool GetDrawViewRegion(void) const
Definition: ViewFEMesh.hh:66
void CreateDefaultAxes()
Create a default set of custom-made axes.
Definition: ViewFEMesh.cc:194
ViewFEMesh()
Constructor.
Definition: ViewFEMesh.cc:15
void SetComponent(ComponentFieldMap *comp)
Set the component from which to retrieve the mesh and field.
Definition: ViewFEMesh.cc:40
void SetFillMeshWithBorders()
Show filled mesh elements.
Definition: ViewFEMesh.hh:79
void SetViewDrift(ViewDrift *vd)
Set the optional associated ViewDrift.
Definition: ViewFEMesh.hh:76
void SetFillColor(int matID, int colorID)
Definition: ViewFEMesh.hh:71
void DisableMaterial(int materialID)
Disable a material so that its mesh cells are not drawn.
Definition: ViewFEMesh.hh:88
void SetYaxisTitle(const char *ytitle)
Definition: ViewFEMesh.cc:189
bool Plot()
Plot method to be called by user.
Definition: ViewFEMesh.cc:71
void SetDrawViewRegion(bool do_draw)
Display intersection of projection plane with viewing area.
Definition: ViewFEMesh.hh:65
void SetPlane(double fx, double fy, double fz, double x0, double y0, double z0)
Set the projection plane.
Definition: ViewFEMesh.cc:120
void SetArea()
Set area to be plotted to the default.
Definition: ViewFEMesh.cc:67
~ViewFEMesh()=default
Destructor.
void SetColor(int matID, int colorID)
Definition: ViewFEMesh.hh:70
void SetXaxis(TGaxis *ax)
Definition: ViewFEMesh.cc:178
void SetYaxis(TGaxis *ay)
Definition: ViewFEMesh.cc:181
void SetFillMesh(const bool f)
Element fill switch; 2D only, set false for wireframe mesh.
Definition: ViewFEMesh.hh:62
void SetDefaultProjection()
Reset the projection plane.
Definition: ViewFEMesh.cc:25