Garfield++ v1r0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
Garfield::ViewMedium Class Reference

#include <ViewMedium.hh>

Public Member Functions

 ViewMedium ()
 
 ~ViewMedium ()
 
void SetCanvas (TCanvas *c)
 
void SetMedium (Medium *m)
 
void SetElectricFieldRange (const double emin, const double emax)
 
void SetMagneticFieldRange (const double bmin, const double bmax)
 
void SetBAngleRange (const double amin, const double amax)
 
void SetFunctionRange (const double vmin, const double vmax)
 
void SetFunctionRange ()
 
void PlotElectronVelocity (const char xaxis, const double e, const double b, const double a)
 
void PlotHoleVelocity (const char xaxis, const double e, const double b, const double a)
 
void PlotIonVelocity (const char xaxis, const double e, const double b, const double a)
 
void PlotElectronDiffusion (const char xaxis, const double e, const double b, const double a)
 
void PlotHoleDiffusion (const char xaxis, const double e, const double b, const double a)
 
void PlotIonDiffusion (const char xaxis, const double e, const double b, const double a)
 
void PlotElectronTownsend (const char xaxis, const double e, const double b, const double a)
 
void PlotHoleTownsend (const char xaxis, const double e, const double b, const double a)
 
void PlotElectronAttachment (const char xaxis, const double e, const double b, const double a)
 
void PlotHoleAttachment (const char xaxis, const double e, const double b, const double a)
 
void PlotElectronCrossSections ()
 
double EvaluateFunction (double *pos, double *par)
 

Detailed Description

Definition at line 15 of file ViewMedium.hh.

Constructor & Destructor Documentation

◆ ViewMedium()

Garfield::ViewMedium::ViewMedium ( )

Definition at line 15 of file ViewMedium.cc.

16 : m_debug(false),
17 m_canvas(0),
18 m_hasExternalCanvas(false),
19 m_medium(0),
20 m_eMin(0.),
21 m_eMax(1000.),
22 m_bMin(0.),
23 m_bMax(5.),
24 m_aMin(0.),
25 m_aMax(3.14),
26 m_vMin(0.),
27 m_vMax(0.),
28 m_efield(500.),
29 m_bfield(1.e2),
30 m_angle(0.),
31 m_etolerance(1.),
32 m_btolerance(0.01),
33 m_atolerance(0.05),
34 m_nFunctions(0),
35 m_nGraphs(0) {
36
37 m_className = "ViewMedium";
38
39 m_functions.clear();
40 m_graphs.clear();
42}
PlottingEngineRoot plottingEngine

◆ ~ViewMedium()

Garfield::ViewMedium::~ViewMedium ( )

Definition at line 44 of file ViewMedium.cc.

44 {
45
46 if (!m_hasExternalCanvas && m_canvas != 0) delete m_canvas;
47}

Member Function Documentation

◆ EvaluateFunction()

double Garfield::ViewMedium::EvaluateFunction ( double *  pos,
double *  par 
)

Definition at line 632 of file ViewMedium.cc.

632 {
633 // to be modified to include B and angle
634
635 if (m_medium == 0) return 0.;
636 int type = int(par[0]);
637 char xaxis = char(par[1]);
638 const double x = pos[0];
639 double y = 0.;
640
641 // Auxiliary variables
642 double value = 0., a = 0., b = 0., c = 0., alongx = 0., alongy = 0.,
643 alongz = 0.;
644
645 switch (type) {
646 case 0:
647 // Electron drift velocity
648 if (xaxis == 'e') { // plot with respect to E field
649 if (!m_medium->ElectronVelocity(x, 0, 0, m_bfield * cos(m_angle),
650 m_bfield * sin(m_angle), 0, a, b, c))
651 return 0.;
652 } else if (xaxis == 'b') { // plot wrt B field
653 if (!m_medium->ElectronVelocity(m_efield, 0, 0, x * cos(m_angle),
654 x * sin(m_angle), 0, a, b, c))
655 return 0.;
656 } else if (xaxis == 'a') { // plot wrt angle
657 if (!m_medium->ElectronVelocity(m_efield, 0, 0, m_bfield * cos(x),
658 m_bfield * sin(x), 0, a, b, c))
659 return 0.;
660 }
661 y = fabs(a);
662 break;
663 case 1:
664 // Electron transverse diffusion
665 if (!m_medium->ElectronDiffusion(x, 0, 0, 0, 0, 0, a, b)) return 0.;
666 y = b;
667 break;
668 case 2:
669 // Electron longitudinal diffusion
670 if (!m_medium->ElectronDiffusion(x, 0, 0, 0, 0, 0, a, b)) return 0.;
671 y = a;
672 break;
673 case 3:
674 // Electron Townsend coefficient
675 if (!m_medium->ElectronTownsend(x, 0, 0, 0, 0, 0, a)) return 0.;
676 y = a;
677 break;
678 case 4:
679 // Electron attachment coefficient
680 if (!m_medium->ElectronAttachment(x, 0, 0, 0, 0, 0, a)) return 0.;
681 y = a;
682 break;
683 case 10:
684 // Hole drift velocity
685 if (!m_medium->HoleVelocity(x, 0, 0, 0, 0, 0, a, b, c)) return 0.;
686 y = a;
687 break;
688 case 11:
689 // Hole transverse diffusion
690 if (!m_medium->HoleDiffusion(x, 0, 0, 0, 0, 0, a, b)) return 0.;
691 y = b;
692 break;
693 case 12:
694 // Hole longitudinal diffusion
695 if (!m_medium->HoleDiffusion(x, 0, 0, 0, 0, 0, a, b)) return 0.;
696 y = a;
697 break;
698 case 13:
699 // Hole Townsend coefficient
700 if (!m_medium->HoleTownsend(x, 0, 0, 0, 0, 0, a)) return 0.;
701 y = a;
702 break;
703 case 14:
704 // Hole attachment coefficient
705 if (!m_medium->HoleAttachment(x, 0, 0, 0, 0, 0, a)) return 0.;
706 y = a;
707 break;
708 case 20:
709 // Ion drift velocity
710 if (!m_medium->IonVelocity(x, 0, 0, 0, 0, 0, a, b, c)) return 0.;
711 y = fabs(a);
712 break;
713 case 21:
714 // Ion transverse diffusion
715 if (!m_medium->IonDiffusion(x, 0, 0, 0, 0, 0, a, b)) return 0.;
716 y = b;
717 break;
718 case 22:
719 // Ion longitudinal diffusion
720 if (!m_medium->IonDiffusion(x, 0, 0, 0, 0, 0, a, b)) return 0.;
721 y = a;
722 break;
723 case 23:
724 // Electron drift velocity along B
725 if (xaxis == 'e') { // plot with respect to E field
726 if (!m_medium->ElectronVelocity(x * cos(m_angle), x * sin(m_angle), 0,
727 m_bfield, 0, 0, value, alongy, alongz))
728 return 0.;
729 } else if (xaxis == 'b') { // plot wrt B field
730 if (!m_medium->ElectronVelocity(m_efield * cos(m_angle),
731 m_efield * sin(m_angle), 0, x, 0, 0,
732 value, alongy, alongz))
733 return 0.;
734 } else if (xaxis == 'a') { // plot wrt angle
735 if (!m_medium->ElectronVelocity(m_efield * cos(x), m_efield * sin(x), 0,
736 m_bfield, 0, 0, value, alongy, alongz))
737 return 0.;
738 }
739 y = fabs(value);
740 break;
741 case 24:
742 // Electron drift velocity along ExB
743 if (xaxis == 'e') { // plot with respect to E field
744 if (!m_medium->ElectronVelocity(x * cos(m_angle), x * sin(m_angle), 0,
745 m_bfield, 0, 0, alongx, alongy, value))
746 return 0.;
747 } else if (xaxis == 'b') { // plot wrt B field
748 if (!m_medium->ElectronVelocity(m_efield * cos(m_angle),
749 m_efield * sin(m_angle), 0, x, 0, 0,
750 alongx, alongy, value))
751 return 0.;
752 } else if (xaxis == 'a') { // plot wrt angle
753 if (!m_medium->ElectronVelocity(m_efield * cos(x), m_efield * sin(x), 0,
754 m_bfield, 0, 0, alongx, alongy, value))
755 return 0.;
756 }
757 y = fabs(value);
758 break;
759 case 25:
760 // Hole drift velocity along B
761 if (xaxis == 'e') { // plot with respect to E field
762 if (!m_medium->HoleVelocity(x * cos(m_angle), x * sin(m_angle), 0,
763 m_bfield, 0, 0, value, alongy, alongz))
764 return 0.;
765 } else if (xaxis == 'b') { // plot wrt B field
766 if (!m_medium->HoleVelocity(m_efield * cos(m_angle),
767 m_efield * sin(m_angle), 0, x, 0, 0, value,
768 alongy, alongz))
769 return 0.;
770 } else if (xaxis == 'a') { // plot wrt angle
771 if (!m_medium->HoleVelocity(m_efield * cos(x), m_efield * sin(x), 0,
772 m_bfield, 0, 0, value, alongy, alongz))
773 return 0.;
774 }
775 y = fabs(value);
776 break;
777 case 26:
778 // Hole drift velocity along ExB
779 if (xaxis == 'e') { // plot with respect to E field
780 if (!m_medium->HoleVelocity(x * cos(m_angle), x * sin(m_angle), 0,
781 m_bfield, 0, 0, alongx, alongy, value))
782 return 0.;
783 } else if (xaxis == 'b') { // plot wrt B field
784 if (!m_medium->HoleVelocity(m_efield * cos(m_angle),
785 m_efield * sin(m_angle), 0, x, 0, 0, alongx,
786 alongy, value))
787 return 0.;
788 } else if (xaxis == 'a') { // plot wrt angle
789 if (!m_medium->HoleVelocity(m_efield * cos(x), m_efield * sin(x), 0,
790 m_bfield, 0, 0, alongx, alongy, value))
791 return 0.;
792 }
793 y = fabs(value);
794 break;
795 default:
796 std::cerr << m_className << "::EvaluateFunction:\n";
797 std::cerr << " Unknown type of transport coefficient requested.\n";
798 std::cerr << " Program bug!\n";
799 return 0.;
800 }
801
802 return y;
803}
DoubleAc cos(const DoubleAc &f)
Definition: DoubleAc.cpp:431
DoubleAc sin(const DoubleAc &f)
Definition: DoubleAc.cpp:383
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:616
virtual bool HoleTownsend(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
Definition: Medium.cc:1079
virtual bool HoleVelocity(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &vx, double &vy, double &vz)
Definition: Medium.cc:773
virtual bool ElectronVelocity(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &vx, double &vy, double &vz)
Definition: Medium.cc:217
virtual bool IonVelocity(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &vx, double &vy, double &vz)
Definition: Medium.cc:1211
virtual bool ElectronDiffusion(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &dl, double &dt)
Definition: Medium.cc:388
virtual bool ElectronTownsend(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &alpha)
Definition: Medium.cc:542
virtual bool ElectronAttachment(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
Definition: Medium.cc:609
virtual bool HoleDiffusion(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &dl, double &dt)
Definition: Medium.cc:928
virtual bool HoleAttachment(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &eta)
Definition: Medium.cc:1144
virtual bool IonDiffusion(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &dl, double &dt)
Definition: Medium.cc:1268

◆ PlotElectronAttachment()

void Garfield::ViewMedium::PlotElectronAttachment ( const char  xaxis,
const double  e,
const double  b,
const double  a 
)

Definition at line 253 of file ViewMedium.cc.

254 {
255
256 bool keep = false;
257 SetupCanvas();
258 AddFunction(m_eMin, m_eMax, 0., 0., keep, "electric field [V/cm]",
259 "Attachment coefficient [1/cm]", 4, xaxis, e, b, a);
260 m_canvas->Update();
261}

◆ PlotElectronCrossSections()

void Garfield::ViewMedium::PlotElectronCrossSections ( )

Definition at line 273 of file ViewMedium.cc.

273 {
274
275 std::cerr << m_className << "::PlotElectronCrossSections:\n";
276 std::cerr << " Function not yet implemented.\n";
277 SetupCanvas();
278}

◆ PlotElectronDiffusion()

void Garfield::ViewMedium::PlotElectronDiffusion ( const char  xaxis,
const double  e,
const double  b,
const double  a 
)

Definition at line 193 of file ViewMedium.cc.

194 {
195
196 bool keep = false;
197 SetupCanvas();
198 AddFunction(m_eMin, m_eMax, m_vMin, m_vMax, keep, "electric field [V/cm]",
199 "diffusion coefficient [#sqrt{cm}]", 1, xaxis, e, b, a);
200 keep = true;
201 AddFunction(m_eMin, m_eMax, m_vMin, m_vMax, keep, "electric field [V/cm]",
202 "diffusion coefficient [#sqrt{cm}]", 2, xaxis, e, b, a);
203
204 m_canvas->Update();
205}

◆ PlotElectronTownsend()

void Garfield::ViewMedium::PlotElectronTownsend ( const char  xaxis,
const double  e,
const double  b,
const double  a 
)

Definition at line 233 of file ViewMedium.cc.

234 {
235
236 bool keep = false;
237 SetupCanvas();
238 AddFunction(m_eMin, m_eMax, 0., 0., keep, "electric field [V/cm]",
239 "Townsend coefficient [1/cm]", 3, xaxis, e, b, a);
240 m_canvas->Update();
241}

◆ PlotElectronVelocity()

void Garfield::ViewMedium::PlotElectronVelocity ( const char  xaxis,
const double  e,
const double  b,
const double  a 
)

Definition at line 121 of file ViewMedium.cc.

122 {
123
124 bool keep = false;
125 SetupCanvas();
126 double min = 0., max = 0.;
127 std::string title = "";
128 if (xaxis == 'e') {
129 title = "electric field [V/cm]";
130 min = m_eMin;
131 max = m_eMax;
132 } else if (xaxis == 'b') {
133 title = "magnetic field [T]";
134 min = m_bMin;
135 max = m_bMax;
136 } else if (xaxis == 'a') {
137 title = "magnetic field angle [rad]";
138 min = m_aMin;
139 max = m_aMax;
140 }
141 AddFunction(min, max, m_vMin, m_vMax, keep, title, "drift velocity [cm/ns]",
142 0, xaxis, e, b, a);
143 keep = true;
144 AddFunction(min, max, m_vMin, m_vMax, keep, title, "drift velocity [cm/ns]",
145 23, xaxis, e, b, a);
146 keep = true;
147 AddFunction(min, max, m_vMin, m_vMax, keep, title, "drift velocity [cm/ns]",
148 24, xaxis, e, b, a);
149 m_canvas->Update();
150}

◆ PlotHoleAttachment()

void Garfield::ViewMedium::PlotHoleAttachment ( const char  xaxis,
const double  e,
const double  b,
const double  a 
)

Definition at line 263 of file ViewMedium.cc.

264 {
265
266 bool keep = false;
267 SetupCanvas();
268 AddFunction(m_eMin, m_eMax, 0., 0., keep, "electric field [V/cm]",
269 "Attachment coefficient [1/cm]", 14, xaxis, e, b, a);
270 m_canvas->Update();
271}

◆ PlotHoleDiffusion()

void Garfield::ViewMedium::PlotHoleDiffusion ( const char  xaxis,
const double  e,
const double  b,
const double  a 
)

Definition at line 207 of file ViewMedium.cc.

208 {
209
210 bool keep = false;
211 SetupCanvas();
212 AddFunction(m_eMin, m_eMax, m_vMin, m_vMax, keep, "electric field [V/cm]",
213 "diffusion coefficient [#sqrt{cm}]", 11, xaxis, e, b, a);
214 keep = true;
215 AddFunction(m_eMin, m_eMax, m_vMin, m_vMax, keep, "electric field [V/cm]",
216 "diffusion coefficient [#sqrt{cm}]", 12, xaxis, e, b, a);
217 m_canvas->Update();
218}

◆ PlotHoleTownsend()

void Garfield::ViewMedium::PlotHoleTownsend ( const char  xaxis,
const double  e,
const double  b,
const double  a 
)

Definition at line 243 of file ViewMedium.cc.

244 {
245
246 bool keep = false;
247 SetupCanvas();
248 AddFunction(m_eMin, m_eMax, 0., 0., keep, "electric field [V/cm]",
249 "Townsend coefficient [1/cm]", 13, xaxis, e, b, a);
250 m_canvas->Update();
251}

◆ PlotHoleVelocity()

void Garfield::ViewMedium::PlotHoleVelocity ( const char  xaxis,
const double  e,
const double  b,
const double  a 
)

Definition at line 152 of file ViewMedium.cc.

153 {
154
155 bool keep = false;
156 SetupCanvas();
157 double min = 0., max = 0.;
158 std::string title = "";
159 if (xaxis == 'e') {
160 title = "electric field [V/cm]";
161 min = m_eMin;
162 max = m_eMax;
163 } else if (xaxis == 'b') {
164 title = "magnetic field [T]";
165 min = m_bMin;
166 max = m_bMax;
167 } else if (xaxis == 'a') {
168 title = "magnetic field angle [rad]";
169 min = m_aMin;
170 max = m_aMax;
171 }
172 AddFunction(min, max, m_vMin, m_vMax, keep, title,
173 "drift velocity [cm/ns]", 10, xaxis, e, b, a);
174 keep = true;
175 AddFunction(min, max, m_vMin, m_vMax, keep, title,
176 "drift velocity [cm/ns]", 25, xaxis, e, b, a);
177 keep = true;
178 AddFunction(min, max, m_vMin, m_vMax, keep, title,
179 "drift velocity [cm/ns]", 26, xaxis, e, b, a);
180 m_canvas->Update();
181}

◆ PlotIonDiffusion()

void Garfield::ViewMedium::PlotIonDiffusion ( const char  xaxis,
const double  e,
const double  b,
const double  a 
)

Definition at line 220 of file ViewMedium.cc.

221 {
222
223 bool keep = false;
224 SetupCanvas();
225 AddFunction(m_eMin, m_eMax, m_vMin, m_vMax, keep, "electric field [V/cm]",
226 "diffusion coefficient [#sqrt{cm}]", 21, xaxis, e, b, a);
227 keep = true;
228 AddFunction(m_eMin, m_eMax, m_vMin, m_vMax, keep, "electric field [V/cm]",
229 "diffusion coefficient [#sqrt{cm}]", 22, xaxis, e, b, a);
230 m_canvas->Update();
231}

◆ PlotIonVelocity()

void Garfield::ViewMedium::PlotIonVelocity ( const char  xaxis,
const double  e,
const double  b,
const double  a 
)

Definition at line 183 of file ViewMedium.cc.

184 {
185
186 bool keep = false;
187 SetupCanvas();
188 AddFunction(m_eMin, m_eMax, m_vMin, m_vMax, keep, "electric field [V/cm]",
189 "drift velocity [cm/ns]", 20, xaxis, e, b, a);
190 m_canvas->Update();
191}

◆ SetBAngleRange()

void Garfield::ViewMedium::SetBAngleRange ( const double  amin,
const double  amax 
)

Definition at line 95 of file ViewMedium.cc.

95 {
96
97 if (amin >= amax || amin < 0.) {
98 std::cerr << m_className << "::SetBAngleRange:\n";
99 std::cerr << " Incorrect field range.\n";
100 return;
101 }
102
103 m_aMin = amin;
104 m_aMax = amax;
105}

◆ SetCanvas()

void Garfield::ViewMedium::SetCanvas ( TCanvas *  c)

Definition at line 49 of file ViewMedium.cc.

49 {
50
51 if (c == 0) return;
52 if (!m_hasExternalCanvas && m_canvas != 0) {
53 delete m_canvas;
54 m_canvas = 0;
55 }
56 m_canvas = c;
57 m_hasExternalCanvas = true;
58}

◆ SetElectricFieldRange()

void Garfield::ViewMedium::SetElectricFieldRange ( const double  emin,
const double  emax 
)

Definition at line 71 of file ViewMedium.cc.

71 {
72
73 if (emin >= emax || emin < 0.) {
74 std::cerr << m_className << "::SetElectricFieldRange:\n";
75 std::cerr << " Incorrect field range.\n";
76 return;
77 }
78
79 m_eMin = emin;
80 m_eMax = emax;
81}

◆ SetFunctionRange() [1/2]

void Garfield::ViewMedium::SetFunctionRange ( )

Definition at line 119 of file ViewMedium.cc.

119{ m_vMin = m_vMax = 0.; }

◆ SetFunctionRange() [2/2]

void Garfield::ViewMedium::SetFunctionRange ( const double  vmin,
const double  vmax 
)

Definition at line 107 of file ViewMedium.cc.

107 {
108
109 if (vmin >= vmax || vmin < 0.) {
110 std::cerr << m_className << "::SetFunctionRange:\n";
111 std::cerr << " Incorrect range.\n";
112 return;
113 }
114
115 m_vMin = vmin;
116 m_vMax = vmax;
117}

◆ SetMagneticFieldRange()

void Garfield::ViewMedium::SetMagneticFieldRange ( const double  bmin,
const double  bmax 
)

Definition at line 83 of file ViewMedium.cc.

83 {
84
85 if (bmin >= bmax || bmin < 0.) {
86 std::cerr << m_className << "::SetMagneticFieldRange:\n";
87 std::cerr << " Incorrect field range.\n";
88 return;
89 }
90
91 m_bMin = bmin;
92 m_bMax = bmax;
93}

◆ SetMedium()

void Garfield::ViewMedium::SetMedium ( Medium m)

Definition at line 60 of file ViewMedium.cc.

60 {
61
62 if (m == 0) {
63 std::cerr << m_className << "::SetMedium:\n";
64 std::cerr << " Medium pointer is null.\n";
65 return;
66 }
67
68 m_medium = m;
69}

The documentation for this class was generated from the following files: