Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization 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 <TGeoManager.h>
5#include <TMatrixDfwd.h>
6
7#include <map>
8#include <string>
9#include <vector>
10
11#include "Garfield/ViewBase.hh"
12
13class TGaxis;
14class TGeoVolume;
15class TGeoMedium;
16
17namespace Garfield {
18
19class Component;
20class ComponentCST;
21class ViewDrift;
22
24
25class ViewFEMesh : public ViewBase {
26 public:
28 ViewFEMesh() : ViewFEMesh(nullptr) {}
30 ViewFEMesh(Component* cmp);
33
35 void SetComponent(Component* cmp);
36
37 void SetPlane(const double fx, const double fy, const double fz,
38 const double x0, const double y0, const double z0) override;
39 void SetPlane(const double fx, const double fy, const double fz,
40 const double x0, const double y0, const double z0,
41 const double hx, const double hy, const double hz) override;
42
43 // Axes
44 void SetXaxis(TGaxis* ax);
45 void SetYaxis(TGaxis* ay);
46 void SetXaxisTitle(const std::string& xtitle) { m_xaxisTitle = xtitle; }
47 void SetYaxisTitle(const std::string& ytitle) { m_yaxisTitle = ytitle; }
48 void EnableAxes() { m_drawAxes = true; }
49 void DisableAxes() { m_drawAxes = false; }
50
52 bool Plot(const bool twod = true, const bool outline = false);
53
55 void SetFillMesh(const bool f) { m_fillMesh = f; }
56
58 void SetDrawViewRegion(bool do_draw) { m_drawViewRegion = do_draw; }
59 bool GetDrawViewRegion(void) const { return m_drawViewRegion; }
60
63 void SetColor(int matID, int colorID) { m_colorMap[matID] = colorID; }
64 void SetColors(const std::map<int, int>& colors) {
65 for (const auto& c : colors) SetColor(c.first, c.second);
66 }
67 void SetFillColor(int matID, int colorID) {
68 m_colorMap_fill[matID] = colorID;
69 }
70 void SetFillColors(const std::map<int, int>& colors) {
71 for (const auto& c : colors) SetFillColor(c.first, c.second);
72 }
73
75
78 m_plotMeshBorders = true;
79 m_fillMesh = true;
80 }
81
84
86 void DisableMaterial(int materialID) {
87 m_disabledMaterial[materialID] = true;
88 }
89
90 private:
91 // Options
92 bool m_fillMesh = false;
93
94 // Intersection of viewing plane with plotted area in planar coordinates
95 bool m_drawViewRegion = false;
96 std::vector<double> m_viewRegionX;
97 std::vector<double> m_viewRegionY;
98
99 // Field map object
100 Component* m_cmp = nullptr;
101
102 // Optional associated ViewDrift object
104 bool m_plotMeshBorders = false;
105
106 // Axes
107 TGaxis* m_xaxis = nullptr;
108 TGaxis* m_yaxis = nullptr;
109 std::string m_xaxisTitle = "";
110 std::string m_yaxisTitle = "";
111 bool m_drawAxes = false;
112
113 // The color map
114 std::map<int, int> m_colorMap;
115 std::map<int, int> m_colorMap_fill;
116
117 // Disabled materials -> not shown in the mesh view
118 std::map<int, bool> m_disabledMaterial;
119
120 std::vector<TGeoVolume*> m_volumes;
121 std::vector<TGeoMedium*> m_media;
122 std::unique_ptr<TGeoManager> m_geoManager;
123
124 // Element plotting methods
127 void DrawCST(ComponentCST* componentCST);
130
131 typedef std::vector<size_t> Facet;
132 void AddFacets(const size_t i,
133 const std::vector<std::vector<Facet> >& elementFacets,
134 const std::map<Facet, std::vector<size_t> >& facetElements,
135 std::vector<Facet>& facets, std::vector<bool>& done) const;
136 bool FacetSign(const Facet& f, const size_t element) const;
139
141
143 bool InView(const double x, const double y) const;
144
145 bool LinesCrossed(double x1, double y1, double x2, double y2, double u1,
146 double v1, double u2, double v2, double& xc,
147 double& yc) const;
148 bool IntersectPlaneArea(double& xmin, double& ymin, double& xmax,
149 double& ymax);
150 bool OnLine(double x1, double y1, double x2, double y2, double u,
151 double v) const;
152 void RemoveCrossings(std::vector<double>& x, std::vector<double>& y);
153 bool PlaneCut(double x1, double y1, double z1, double x2, double y2,
154 double z2, TMatrixD& xMat);
155 void ClipToView(std::vector<double>& px, std::vector<double>& py,
156 std::vector<double>& cx, std::vector<double>& cy);
157 bool IsInPolygon(double x, double y, const std::vector<double>& px,
158 const std::vector<double>& py, bool& edge) const;
159
160 void Reset();
161};
162
164
165} // namespace Garfield
166#endif
Component for importing and interpolating field maps from CST.
ViewBase()=delete
Default constructor.
Visualize drift lines and tracks.
Definition ViewDrift.hh:18
Draw the mesh of a field-map component.
Definition ViewFEMesh.hh:25
bool GetDrawViewRegion(void) const
Definition ViewFEMesh.hh:59
void CreateDefaultAxes()
Create a default set of custom-made axes.
bool PlaneCut(double x1, double y1, double z1, double x2, double y2, double z2, TMatrixD &xMat)
void SetYaxisTitle(const std::string &ytitle)
Definition ViewFEMesh.hh:47
ViewFEMesh()
Default constructor.
Definition ViewFEMesh.hh:28
std::vector< double > m_viewRegionY
Definition ViewFEMesh.hh:97
void SetFillMeshWithBorders()
Show filled mesh elements.
Definition ViewFEMesh.hh:77
bool Plot(const bool twod=true, const bool outline=false)
Plot method to be called by user.
void SetViewDrift(ViewDrift *vd)
Set the optional associated ViewDrift.
Definition ViewFEMesh.hh:74
std::vector< size_t > Facet
bool IntersectPlaneArea(double &xmin, double &ymin, double &xmax, double &ymax)
void SetFillColor(int matID, int colorID)
Definition ViewFEMesh.hh:67
void SetXaxisTitle(const std::string &xtitle)
Definition ViewFEMesh.hh:46
void DisableMaterial(int materialID)
Disable a material so that its mesh cells are not drawn.
Definition ViewFEMesh.hh:86
void RemoveCrossings(std::vector< double > &x, std::vector< double > &y)
ViewDrift * m_viewDrift
void SetPlane(const double fx, const double fy, const double fz, const double x0, const double y0, const double z0) override
Set the projection (viewing plane), if applicable.
std::string m_xaxisTitle
std::vector< double > m_viewRegionX
Definition ViewFEMesh.hh:96
std::vector< TGeoVolume * > m_volumes
std::unique_ptr< TGeoManager > m_geoManager
void SetDrawViewRegion(bool do_draw)
Display intersection of projection plane with viewing area.
Definition ViewFEMesh.hh:58
bool LinesCrossed(double x1, double y1, double x2, double y2, double u1, double v1, double u2, double v2, double &xc, double &yc) const
void DrawCST(ComponentCST *componentCST)
std::vector< TGeoMedium * > m_media
bool IsInPolygon(double x, double y, const std::vector< double > &px, const std::vector< double > &py, bool &edge) const
bool InView(const double x, const double y) const
Return true if the specified point is in the view region.
bool OnLine(double x1, double y1, double x2, double y2, double u, double v) const
std::map< int, bool > m_disabledMaterial
void AddFacets(const size_t i, const std::vector< std::vector< Facet > > &elementFacets, const std::map< Facet, std::vector< size_t > > &facetElements, std::vector< Facet > &facets, std::vector< bool > &done) const
void SetFillColors(const std::map< int, int > &colors)
Definition ViewFEMesh.hh:70
void SetColor(int matID, int colorID)
Associate a color with each element material map ID.
Definition ViewFEMesh.hh:63
void SetComponent(Component *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, const double hx, const double hy, const double hz) override
Set the projection plane specifying a hint for the in-plane x axis.
std::map< int, int > m_colorMap_fill
std::string m_yaxisTitle
~ViewFEMesh()
Destructor.
void SetXaxis(TGaxis *ax)
void SetYaxis(TGaxis *ay)
void SetColors(const std::map< int, int > &colors)
Definition ViewFEMesh.hh:64
bool FacetSign(const Facet &f, const size_t element) const
void SetFillMesh(const bool f)
Element fill switch; 2D only, set false for wireframe mesh.
Definition ViewFEMesh.hh:55
void ClipToView(std::vector< double > &px, std::vector< double > &py, std::vector< double > &cx, std::vector< double > &cy)
std::map< int, int > m_colorMap
ViewFEMesh(Component *cmp)
Constructor.
ViewFEMesh ViewMesh