Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization measurement in gases and semiconductors
Loading...
Searching...
No Matches
ViewCell.hh
Go to the documentation of this file.
1#ifndef G_VIEW_CELL
2#define G_VIEW_CELL
3
4#include <TGeoManager.h>
5
6#include <memory>
7
9
10namespace Garfield {
11
14
16
17class ViewCell : public ViewBase {
18 public:
26 ~ViewCell() = default;
27
31
33 void Plot2d();
35 void Plot3d();
36
39 void EnableWireMarkers(const bool on = true) { m_useWireMarker = on; }
41
42 private:
43 bool m_useWireMarker = true;
44
47
48 // 3D geometry.
49 std::unique_ptr<TGeoManager> m_geo;
50
51 bool Plot(const bool twod);
52 // Draw a wire in 2D.
53 void PlotWire(const double x, const double y, const double d, const int type);
54 // Draw a wire in 3D.
55 void PlotWire(const double x, const double y, const double d, const int type,
56 const double lz);
57 // Draw a tube in 2D.
58 void PlotTube(const double x0, const double y0, const double r, const int n);
59 // Draw a tube in 3D.
60 void PlotTube(const double x0, const double y0, const double r1,
61 const double r2, const int n, const double lz);
62 // Draw a plane in 2D.
63 void PlotPlane(const double x0, const double y0, const double x1,
64 const double y1);
65 // Draw a plane in 3D.
66 void PlotPlane(const double dx, const double dy, const double dz,
67 const double x0, const double y0);
68 // Draw a neBEM 2D layout.
69 bool PlotNeBem(const bool twod);
70 // Setup the TGeoManager.
71 void SetupGeo(const double dx, const double dy, const double dz);
72};
73} // namespace Garfield
74#endif
Semi-analytic calculation of two-dimensional configurations consisting of wires, planes,...
Two-dimensional implementation of the nearly exact Boundary Element Method.
ViewBase()=delete
Default constructor.
~ViewCell()=default
Destructor.
ViewCell(ComponentNeBem2d *cmp)
Constructor from two-dimenstional neBEM component.
ViewCell(ComponentAnalyticField *cmp)
Constructor from analytic-field component.
void PlotWire(const double x, const double y, const double d, const int type)
void PlotPlane(const double x0, const double y0, const double x1, const double y1)
void PlotTube(const double x0, const double y0, const double r1, const double r2, const int n, const double lz)
ComponentAnalyticField * m_component
Definition ViewCell.hh:45
void Plot3d()
Make a three-dimensional drawing of the cell layout (using TGeo).
bool PlotNeBem(const bool twod)
void PlotWire(const double x, const double y, const double d, const int type, const double lz)
bool Plot(const bool twod)
ComponentNeBem2d * m_nebem
Definition ViewCell.hh:46
void PlotTube(const double x0, const double y0, const double r, const int n)
void SetComponent(ComponentAnalyticField *comp)
Set the component for which to draw the cell geometry.
std::unique_ptr< TGeoManager > m_geo
Definition ViewCell.hh:49
void SetComponent(ComponentNeBem2d *comp)
void SetupGeo(const double dx, const double dy, const double dz)
ViewCell()
Default constructor.
void Plot2d()
Make a two-dimensional drawing of the cell layout.
void EnableWireMarkers(const bool on=true)
Visualize wirers using markers (default setting) or as a circle with the actual wire radius.
Definition ViewCell.hh:39
void DisableWireMarkers()
Definition ViewCell.hh:40
void PlotPlane(const double dx, const double dy, const double dz, const double x0, const double y0)