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