Garfield++ 5.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#include <map>
7
8#include <TArrayD.h>
9#include <TGaxis.h>
10#include <TGeoManager.h>
11#include <TMatrixD.h>
12
13#include "ComponentCST.hh"
14#include "ComponentFieldMap.hh"
15#include "ViewBase.hh"
16#include "ViewDrift.hh"
17
18namespace Garfield {
19
20/// Draw the mesh of a field-map component.
21
22class ViewFEMesh : public ViewBase {
23 public:
24 /// Constructor.
25 ViewFEMesh();
26 /// Destructor.
28
29 /// Set the component from which to retrieve the mesh and field.
31
32 void SetPlane(const double fx, const double fy, const double fz,
33 const double x0, const double y0, const double z0) override;
34 void SetPlane(const double fx, const double fy, const double fz,
35 const double x0, const double y0, const double z0,
36 const double hx, const double hy, const double hz) override;
37
38 // Axes
39 void SetXaxis(TGaxis* ax);
40 void SetYaxis(TGaxis* ay);
41 void SetXaxisTitle(const std::string& xtitle) { m_xaxisTitle = xtitle; }
42 void SetYaxisTitle(const std::string& ytitle) { m_yaxisTitle = ytitle; }
43 void EnableAxes() { m_drawAxes = true; }
44 void DisableAxes() { m_drawAxes = false; }
45
46 /// Plot method to be called by user
47 bool Plot(const bool twod = true);
48
49 /// Element fill switch; 2D only, set false for wireframe mesh
50 void SetFillMesh(const bool f) { m_fillMesh = f; }
51
52 /// Display intersection of projection plane with viewing area
53 void SetDrawViewRegion(bool do_draw) { m_drawViewRegion = do_draw; }
54 bool GetDrawViewRegion(void) const { return m_drawViewRegion; }
55
56 /// Associate a color with each element material map ID;
57 /// Uses ROOT color numberings
58 void SetColor(int matID, int colorID) { m_colorMap[matID] = colorID; }
59 void SetFillColor(int matID, int colorID) {
60 m_colorMap_fill[matID] = colorID;
61 }
62
63 /// Set the optional associated ViewDrift
64 void SetViewDrift(ViewDrift* vd) { m_viewDrift = vd; }
65
66 /// Show filled mesh elements
68 m_plotMeshBorders = true;
69 m_fillMesh = true;
70 }
71
72 /// Create a default set of custom-made axes.
73 void CreateDefaultAxes();
74
75 /// Disable a material so that its mesh cells are not drawn
76 void DisableMaterial(int materialID) {
77 m_disabledMaterial[materialID] = true;
78 }
79
80 private:
81 // Options
82 bool m_fillMesh = false;
83
84 // Intersection of viewing plane with plotted area in planar coordinates
85 bool m_drawViewRegion = false;
86 std::vector<double> m_viewRegionX;
87 std::vector<double> m_viewRegionY;
88
89 // The field map object
90 ComponentFieldMap* m_cmp = nullptr;
91
92 // Optional associated ViewDrift object
93 ViewDrift* m_viewDrift = nullptr;
94 bool m_plotMeshBorders = false;
95
96 // Axes
97 TGaxis* m_xaxis = nullptr;
98 TGaxis* m_yaxis = nullptr;
99 std::string m_xaxisTitle = "";
100 std::string m_yaxisTitle = "";
101 bool m_drawAxes = false;
102
103 // The color map
104 std::map<int, int> m_colorMap;
105 std::map<int, int> m_colorMap_fill;
106
107 // Disabled materials -> not shown in the mesh view
108 std::map<int, bool> m_disabledMaterial;
109
110 std::vector<TGeoVolume*> m_volumes;
111 std::vector<TGeoMedium*> m_media;
112 std::unique_ptr<TGeoManager> m_geoManager;
113
114 // Element plotting methods
115 void DrawElements2d();
116 void DrawElements3d();
117 void DrawCST(ComponentCST* componentCST);
118
119 void DrawDriftLines2d();
120 void DrawDriftLines3d();
121
122 bool GetPlotLimits();
123
124 /// Return true if the specified point is in the view region.
125 bool InView(const double x, const double y) const;
126
127 bool LinesCrossed(double x1, double y1, double x2, double y2, double u1,
128 double v1, double u2, double v2, double& xc,
129 double& yc) const;
130 bool IntersectPlaneArea(double& xmin, double& ymin,
131 double& xmax, double& ymax);
132 bool OnLine(double x1, double y1, double x2, double y2, double u,
133 double v) const;
134 void RemoveCrossings(std::vector<double>& x, std::vector<double>& y);
135 bool PlaneCut(double x1, double y1, double z1, double x2, double y2,
136 double z2, TMatrixD& xMat);
137 void ClipToView(std::vector<double>& px, std::vector<double>& py,
138 std::vector<double>& cx, std::vector<double>& cy);
139 bool IsInPolygon(double x, double y, const std::vector<double>& px,
140 const std::vector<double>& py, bool& edge) const;
141
142 void Reset();
143};
144} // namespace Garfield
145#endif
Base class for components based on finite-element field maps.
ViewBase()=delete
Default constructor.
Visualize drift lines and tracks.
Definition ViewDrift.hh:19
bool GetDrawViewRegion(void) const
Definition ViewFEMesh.hh:54
bool Plot(const bool twod=true)
Plot method to be called by user.
Definition ViewFEMesh.cc:60
void CreateDefaultAxes()
Create a default set of custom-made axes.
void SetYaxisTitle(const std::string &ytitle)
Definition ViewFEMesh.hh:42
ViewFEMesh()
Constructor.
Definition ViewFEMesh.cc:22
void SetFillMeshWithBorders()
Show filled mesh elements.
Definition ViewFEMesh.hh:67
void SetViewDrift(ViewDrift *vd)
Set the optional associated ViewDrift.
Definition ViewFEMesh.hh:64
void SetFillColor(int matID, int colorID)
Definition ViewFEMesh.hh:59
void SetXaxisTitle(const std::string &xtitle)
Definition ViewFEMesh.hh:41
void DisableMaterial(int materialID)
Disable a material so that its mesh cells are not drawn.
Definition ViewFEMesh.hh:76
void SetPlane(const double fx, const double fy, const double fz, const double x0, const double y0, const double z0) override
void SetDrawViewRegion(bool do_draw)
Display intersection of projection plane with viewing area.
Definition ViewFEMesh.hh:53
void SetComponent(ComponentFieldMap *cmp)
Set the component from which to retrieve the mesh and field.
Definition ViewFEMesh.cc:49
void SetColor(int matID, int colorID)
Definition ViewFEMesh.hh:58
~ViewFEMesh()
Destructor.
Definition ViewFEMesh.cc:24
void SetXaxis(TGaxis *ax)
void SetYaxis(TGaxis *ay)
void SetFillMesh(const bool f)
Element fill switch; 2D only, set false for wireframe mesh.
Definition ViewFEMesh.hh:50