Garfield++ 3.0
A toolkit for the detailed simulation of particle detectors based on ionisation 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 <memory>
5#include <string>
6
7#include <TGeoManager.h>
8
9#include "ViewBase.hh"
10
11namespace Garfield {
12
13class ComponentAnalyticField;
14class ComponentNeBem2d;
15
16/// Visualize the "cell" defined in an analytic-field component.
17
18class ViewCell : public ViewBase {
19 public:
20 /// Constructor
21 ViewCell();
22 /// Destructor
23 ~ViewCell() = default;
24
25 /// Set the component for which to draw the cell geometry.
28
29 /// Set the plot range explicitly.
30 void SetArea(const double xmin, const double ymin, const double zmin,
31 const double xmax, const double ymax, const double zmax);
32 /// Take the plot range from the bounding box of the component class.
33 void SetArea() { m_hasUserArea = false; }
34
35 /// Make a two-dimensional drawing of the cell geometry.
36 void Plot2d();
37 /// Make a three-dimensional drawing of the cell geometry (using TGeo).
38 void Plot3d();
39
40 /// Visualize wirers using markers or as a circle with the actual wire radius.
41 /// The default is markers.
42 void EnableWireMarkers(const bool on = true) { m_useWireMarker = on; }
44
45 private:
46 // Options
47 bool m_useWireMarker = true;
48
49 std::string m_label = "Cell Layout";
50
51 // Box dimensions
52 bool m_hasUserArea = false;
53 double m_xMin = -1., m_yMin = -1., m_zMin = -1.;
54 double m_xMax = 1., m_yMax = 1., m_zMax = 1.;
55
56 ComponentAnalyticField* m_component = nullptr;
57 ComponentNeBem2d* m_nebem = nullptr;
58
59 // 3D geometry.
60 std::unique_ptr<TGeoManager> m_geo;
61
62 bool Plot(const bool use3d);
63 // Draw a wire in 2D.
64 void PlotWire(const double x, const double y, const double d, const int type);
65 // Draw a wire in 3D.
66 void PlotWire(const double x, const double y, const double d, const int type,
67 const double lz);
68 // Draw a tube in 2D.
69 void PlotTube(const double x0, const double y0, const double r, const int n);
70 // Draw a tube in 3D.
71 void PlotTube(const double x0, const double y0,
72 const double r1, const double r2, const int n,
73 const double lz);
74 // Draw a plane in 2D.
75 void PlotPlane(const double x0, const double y0,
76 const double x1, const double y1);
77 // Draw a plane in 3D.
78 void PlotPlane(const double dx, const double dy, const double dz,
79 const double x0, const double y0);
80 // Draw a neBEM 2D layout.
81 bool PlotNeBem(const bool use3d);
82 // Setup the TGeoManager.
83 void SetupGeo(const double dx, const double dy, const double dz);
84};
85}
86#endif
Two-dimensional implementation of the nearly exact Boundary Element Method.
Base class for visualization classes.
Definition: ViewBase.hh:10
Visualize the "cell" defined in an analytic-field component.
Definition: ViewCell.hh:18
~ViewCell()=default
Destructor.
void Plot3d()
Make a three-dimensional drawing of the cell geometry (using TGeo).
Definition: ViewCell.cc:60
void SetArea()
Take the plot range from the bounding box of the component class.
Definition: ViewCell.hh:33
void SetComponent(ComponentAnalyticField *comp)
Set the component for which to draw the cell geometry.
Definition: ViewCell.cc:20
ViewCell()
Constructor.
Definition: ViewCell.cc:18
void Plot2d()
Make a two-dimensional drawing of the cell geometry.
Definition: ViewCell.cc:54
void EnableWireMarkers(const bool on=true)
Definition: ViewCell.hh:42
void DisableWireMarkers()
Definition: ViewCell.hh:43