Garfield++ v2r0
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

Plot transport coefficients as function of electric and magnetic field. More...

#include <ViewMedium.hh>

Public Types

enum  Property {
  ElectronVelocityE , ElectronTransverseDiffusion , ElectronLongitudinalDiffusion , ElectronTownsend ,
  ElectronAttachment , ElectronLorentzAngle , HoleVelocityE = 10 , HoleTransverseDiffusion ,
  HoleLongitudinalDiffusion , HoleTownsend , HoleAttachment , IonVelocity = 20 ,
  IonTransverseDiffusion , IonLongitudinalDiffusion , ElectronVelocityB , ElectronVelocityExB ,
  HoleVelocityB , HoleVelocityExB
}
 

Public Member Functions

 ViewMedium ()
 Constructor.
 
 ~ViewMedium ()
 Destructor.
 
void SetCanvas (TCanvas *c)
 Set the canvas to be painted on.
 
void SetMedium (Medium *m)
 
void SetElectricFieldRange (const double emin, const double emax)
 Set the limits of the electric field.
 
void SetMagneticFieldRange (const double bmin, const double bmax)
 Set the limits of the magnetic field.
 
void SetBAngleRange (const double amin, const double amax)
 Set the limits of the angle between electric and magnetic field.
 
void SetFunctionRange (const double vmin, const double vmax)
 Set the range of the function (velocity etc.) to be plotted.
 
void SetFunctionRange ()
 Use the default function range.
 
void PlotElectronVelocity (const char xaxis, const double e, const double b, const double a=HalfPi)
 
void PlotHoleVelocity (const char xaxis, const double e, const double b, const double a=HalfPi)
 
void PlotIonVelocity (const char xaxis, const double e, const double b, const double a=HalfPi)
 
void PlotElectronDiffusion (const char xaxis, const double e, const double b, const double a=HalfPi)
 
void PlotHoleDiffusion (const char xaxis, const double e, const double b, const double a=HalfPi)
 
void PlotIonDiffusion (const char xaxis, const double e, const double b, const double a=HalfPi)
 
void PlotElectronTownsend (const char xaxis, const double e, const double b, const double a=HalfPi)
 
void PlotHoleTownsend (const char xaxis, const double e, const double b, const double a=HalfPi)
 
void PlotElectronAttachment (const char xaxis, const double e, const double b, const double a=HalfPi)
 
void PlotHoleAttachment (const char xaxis, const double e, const double b, const double a=HalfPi)
 
void PlotElectronLorentzAngle (const char xaxis, const double e, const double b, const double a=HalfPi)
 
void PlotElectronCrossSections ()
 
double EvaluateFunction (double *pos, double *par)
 

Detailed Description

Plot transport coefficients as function of electric and magnetic field.

Definition at line 19 of file ViewMedium.hh.

Member Enumeration Documentation

◆ Property

Enumerator
ElectronVelocityE 
ElectronTransverseDiffusion 
ElectronLongitudinalDiffusion 
ElectronTownsend 
ElectronAttachment 
ElectronLorentzAngle 
HoleVelocityE 
HoleTransverseDiffusion 
HoleLongitudinalDiffusion 
HoleTownsend 
HoleAttachment 
IonVelocity 
IonTransverseDiffusion 
IonLongitudinalDiffusion 
ElectronVelocityB 
ElectronVelocityExB 
HoleVelocityB 
HoleVelocityExB 

Definition at line 68 of file ViewMedium.hh.

68 {
75 HoleVelocityE = 10,
80 IonVelocity = 20,
87 };

Constructor & Destructor Documentation

◆ ViewMedium()

Garfield::ViewMedium::ViewMedium ( )

Constructor.

Definition at line 15 of file ViewMedium.cc.

16 : m_debug(false),
17 m_canvas(NULL),
18 m_hasExternalCanvas(false),
19 m_medium(NULL),
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_labele("electric field [V/cm]"),
35 m_labelb("magnetic field [T]"),
36 m_labela("magnetic field angle [rad]"),
37 m_labelv("drift velocity [cm/ns]"),
38 m_labeld("diffusion coefficient [#sqrt{cm}]") {
39
40 m_className = "ViewMedium";
42}
PlottingEngineRoot plottingEngine

◆ ~ViewMedium()

Garfield::ViewMedium::~ViewMedium ( )

Destructor.

Definition at line 44 of file ViewMedium.cc.

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

Member Function Documentation

◆ EvaluateFunction()

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

Definition at line 619 of file ViewMedium.cc.

619 {
620 // to be modified to include B and angle
621
622 if (m_medium == 0) return 0.;
623 int type = int(par[0]);
624 char xaxis = char(par[1]);
625 const double x = pos[0];
626 double y = 0.;
627
628 // Auxiliary variables
629 double value = 0., a = 0., b = 0., c = 0.;
630 double alongx = 0., alongy = 0., alongz = 0.;
631
632 switch (type) {
634 if (xaxis == 'e') {
635 // plot with respect to E field
636 if (!m_medium->ElectronVelocity(x, 0, 0, m_bfield * cos(m_angle),
637 m_bfield * sin(m_angle), 0, a, b, c))
638 return 0.;
639 } else if (xaxis == 'b') {
640 // plot wrt B field
641 if (!m_medium->ElectronVelocity(m_efield, 0, 0, x * cos(m_angle),
642 x * sin(m_angle), 0, a, b, c))
643 return 0.;
644 } else if (xaxis == 'a') {
645 // plot wrt angle
646 if (!m_medium->ElectronVelocity(m_efield, 0, 0, m_bfield * cos(x),
647 m_bfield * sin(x), 0, a, b, c))
648 return 0.;
649 }
650 y = fabs(a);
651 break;
653 if (!m_medium->ElectronDiffusion(x, 0, 0, 0, 0, 0, a, b)) return 0.;
654 y = b;
655 break;
657 if (!m_medium->ElectronDiffusion(x, 0, 0, 0, 0, 0, a, b)) return 0.;
658 y = a;
659 break;
660 case ElectronTownsend:
661 if (!m_medium->ElectronTownsend(x, 0, 0, 0, 0, 0, a)) return 0.;
662 y = a;
663 break;
665 if (!m_medium->ElectronAttachment(x, 0, 0, 0, 0, 0, a)) return 0.;
666 y = a;
667 break;
669 if (!m_medium->ElectronLorentzAngle(x, 0, 0, 0, 0, 0, a)) return 0.;
670 y = a;
671 break;
672 case HoleVelocityE:
673 if (!m_medium->HoleVelocity(x, 0, 0, 0, 0, 0, a, b, c)) return 0.;
674 y = a;
675 break;
677 if (!m_medium->HoleDiffusion(x, 0, 0, 0, 0, 0, a, b)) return 0.;
678 y = b;
679 break;
681 if (!m_medium->HoleDiffusion(x, 0, 0, 0, 0, 0, a, b)) return 0.;
682 y = a;
683 break;
684 case HoleTownsend:
685 if (!m_medium->HoleTownsend(x, 0, 0, 0, 0, 0, a)) return 0.;
686 y = a;
687 break;
688 case HoleAttachment:
689 if (!m_medium->HoleAttachment(x, 0, 0, 0, 0, 0, a)) return 0.;
690 y = a;
691 break;
692 case IonVelocity:
693 if (!m_medium->IonVelocity(x, 0, 0, 0, 0, 0, a, b, c)) return 0.;
694 y = fabs(a);
695 break;
697 if (!m_medium->IonDiffusion(x, 0, 0, 0, 0, 0, a, b)) return 0.;
698 y = b;
699 break;
701 if (!m_medium->IonDiffusion(x, 0, 0, 0, 0, 0, a, b)) return 0.;
702 y = a;
703 break;
705 if (xaxis == 'e') {
706 // plot with respect to E field
707 if (!m_medium->ElectronVelocity(x * cos(m_angle), x * sin(m_angle), 0,
708 m_bfield, 0, 0, value, alongy, alongz))
709 return 0.;
710 } else if (xaxis == 'b') {
711 // plot wrt B field
712 if (!m_medium->ElectronVelocity(m_efield * cos(m_angle),
713 m_efield * sin(m_angle), 0, x, 0, 0,
714 value, alongy, alongz))
715 return 0.;
716 } else if (xaxis == 'a') {
717 // plot wrt angle
718 if (!m_medium->ElectronVelocity(m_efield * cos(x), m_efield * sin(x), 0,
719 m_bfield, 0, 0, value, alongy, alongz))
720 return 0.;
721 }
722 y = fabs(value);
723 break;
725 if (xaxis == 'e') {
726 // plot with respect to E field
727 if (!m_medium->ElectronVelocity(x * cos(m_angle), x * sin(m_angle), 0,
728 m_bfield, 0, 0, alongx, alongy, value))
729 return 0.;
730 } else if (xaxis == 'b') {
731 // plot wrt B field
732 if (!m_medium->ElectronVelocity(m_efield * cos(m_angle),
733 m_efield * sin(m_angle), 0, x, 0, 0,
734 alongx, alongy, value))
735 return 0.;
736 } else if (xaxis == 'a') {
737 // plot wrt angle
738 if (!m_medium->ElectronVelocity(m_efield * cos(x), m_efield * sin(x), 0,
739 m_bfield, 0, 0, alongx, alongy, value))
740 return 0.;
741 }
742 y = fabs(value);
743 break;
744 case HoleVelocityB:
745 if (xaxis == 'e') {
746 // plot with respect to E field
747 if (!m_medium->HoleVelocity(x * cos(m_angle), x * sin(m_angle), 0,
748 m_bfield, 0, 0, value, alongy, alongz))
749 return 0.;
750 } else if (xaxis == 'b') {
751 // plot wrt B field
752 if (!m_medium->HoleVelocity(m_efield * cos(m_angle),
753 m_efield * sin(m_angle), 0, x, 0, 0, value,
754 alongy, alongz))
755 return 0.;
756 } else if (xaxis == 'a') {
757 // plot wrt angle
758 if (!m_medium->HoleVelocity(m_efield * cos(x), m_efield * sin(x), 0,
759 m_bfield, 0, 0, value, alongy, alongz))
760 return 0.;
761 }
762 y = fabs(value);
763 break;
764 case HoleVelocityExB:
765 if (xaxis == 'e') {
766 // plot with respect to E field
767 if (!m_medium->HoleVelocity(x * cos(m_angle), x * sin(m_angle), 0,
768 m_bfield, 0, 0, alongx, alongy, value))
769 return 0.;
770 } else if (xaxis == 'b') {
771 // plot wrt B field
772 if (!m_medium->HoleVelocity(m_efield * cos(m_angle),
773 m_efield * sin(m_angle), 0, x, 0, 0, alongx,
774 alongy, value))
775 return 0.;
776 } else if (xaxis == 'a') {
777 // plot wrt angle
778 if (!m_medium->HoleVelocity(m_efield * cos(x), m_efield * sin(x), 0,
779 m_bfield, 0, 0, alongx, alongy, value))
780 return 0.;
781 }
782 y = fabs(value);
783 break;
784 default:
785 std::cerr << m_className << "::EvaluateFunction:\n";
786 std::cerr << " Unknown type of transport coefficient requested.\n";
787 std::cerr << " Program bug!\n";
788 return 0.;
789 }
790
791 return y;
792}
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:971
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:704
virtual bool ElectronLorentzAngle(const double ex, const double ey, const double ez, const double bx, const double by, const double bz, double &lor)
Definition: Medium.cc:599
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:204
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:1077
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:362
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:490
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:544
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:846
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:1023
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:1122
DoubleAc cos(const DoubleAc &f)
Definition: DoubleAc.cpp:432
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:615
DoubleAc sin(const DoubleAc &f)
Definition: DoubleAc.cpp:384

◆ PlotElectronAttachment()

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

Definition at line 251 of file ViewMedium.cc.

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

◆ PlotElectronCrossSections()

void Garfield::ViewMedium::PlotElectronCrossSections ( )

Definition at line 283 of file ViewMedium.cc.

283 {
284
285 std::cerr << m_className << "::PlotElectronCrossSections:\n";
286 std::cerr << " Function not yet implemented.\n";
287 SetupCanvas();
288}

◆ PlotElectronDiffusion()

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

Definition at line 191 of file ViewMedium.cc.

192 {
193
194 SetupCanvas();
195 bool keep = false;
196 AddFunction(m_eMin, m_eMax, m_vMin, m_vMax, keep, m_labele, m_labeld,
197 ElectronTransverseDiffusion, xaxis, e, b, a);
198 keep = true;
199 AddFunction(m_eMin, m_eMax, m_vMin, m_vMax, keep, m_labele, m_labeld,
200 ElectronLongitudinalDiffusion, xaxis, e, b, a);
201
202 m_canvas->Update();
203}

◆ PlotElectronLorentzAngle()

void Garfield::ViewMedium::PlotElectronLorentzAngle ( const char  xaxis,
const double  e,
const double  b,
const double  a = HalfPi 
)

Definition at line 272 of file ViewMedium.cc.

273 {
274
275 bool keep = false;
276 SetupCanvas();
277 AddFunction(m_eMin, m_eMax, 0., 0., keep, m_labele,
278 "Lorentz angle",
279 ElectronLorentzAngle, xaxis, e, b, a);
280 m_canvas->Update();
281}

◆ PlotElectronTownsend()

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

Definition at line 231 of file ViewMedium.cc.

232 {
233
234 bool keep = false;
235 SetupCanvas();
236 AddFunction(m_eMin, m_eMax, 0., 0., keep, m_labele,
237 "Townsend coefficient [1/cm]", ElectronTownsend, xaxis, e, b, a);
238 m_canvas->Update();
239}

◆ PlotElectronVelocity()

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

Definition at line 119 of file ViewMedium.cc.

120 {
121
122 SetupCanvas();
123 double min = 0., max = 0.;
124 std::string title = "";
125 if (xaxis == 'e') {
126 title = m_labele;
127 min = m_eMin;
128 max = m_eMax;
129 } else if (xaxis == 'b') {
130 title = m_labelb;
131 min = m_bMin;
132 max = m_bMax;
133 } else if (xaxis == 'a') {
134 title = m_labela;
135 min = m_aMin;
136 max = m_aMax;
137 }
138 bool keep = false;
139 AddFunction(min, max, m_vMin, m_vMax, keep, title, m_labelv,
140 ElectronVelocityE, xaxis, e, b, a);
141 keep = true;
142 AddFunction(min, max, m_vMin, m_vMax, keep, title, m_labelv,
143 ElectronVelocityB, xaxis, e, b, a);
144 keep = true;
145 AddFunction(min, max, m_vMin, m_vMax, keep, title, m_labelv,
146 ElectronVelocityExB, xaxis, e, b, a);
147 m_canvas->Update();
148}

◆ PlotHoleAttachment()

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

Definition at line 262 of file ViewMedium.cc.

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

◆ PlotHoleDiffusion()

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

Definition at line 205 of file ViewMedium.cc.

206 {
207
208 SetupCanvas();
209 bool keep = false;
210 AddFunction(m_eMin, m_eMax, m_vMin, m_vMax, keep, m_labele, m_labeld,
211 HoleTransverseDiffusion, xaxis, e, b, a);
212 keep = true;
213 AddFunction(m_eMin, m_eMax, m_vMin, m_vMax, keep, m_labele, m_labeld,
214 HoleLongitudinalDiffusion, xaxis, e, b, a);
215 m_canvas->Update();
216}

◆ PlotHoleTownsend()

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

Definition at line 241 of file ViewMedium.cc.

242 {
243
244 bool keep = false;
245 SetupCanvas();
246 AddFunction(m_eMin, m_eMax, 0., 0., keep, m_labele,
247 "Townsend coefficient [1/cm]", HoleTownsend, xaxis, e, b, a);
248 m_canvas->Update();
249}

◆ PlotHoleVelocity()

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

Definition at line 150 of file ViewMedium.cc.

151 {
152
153 SetupCanvas();
154 double min = 0., max = 0.;
155 std::string title = "";
156 if (xaxis == 'e') {
157 title = m_labele;
158 min = m_eMin;
159 max = m_eMax;
160 } else if (xaxis == 'b') {
161 title = m_labelb;
162 min = m_bMin;
163 max = m_bMax;
164 } else if (xaxis == 'a') {
165 title = m_labela;
166 min = m_aMin;
167 max = m_aMax;
168 }
169 bool keep = false;
170 AddFunction(min, max, m_vMin, m_vMax, keep, title, m_labelv,
171 HoleVelocityE, xaxis, e, b, a);
172 keep = true;
173 AddFunction(min, max, m_vMin, m_vMax, keep, title, m_labelv,
174 HoleVelocityB, xaxis, e, b, a);
175 keep = true;
176 AddFunction(min, max, m_vMin, m_vMax, keep, title, m_labelv,
177 HoleVelocityExB, xaxis, e, b, a);
178 m_canvas->Update();
179}

◆ PlotIonDiffusion()

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

Definition at line 218 of file ViewMedium.cc.

219 {
220
221 SetupCanvas();
222 bool keep = false;
223 AddFunction(m_eMin, m_eMax, m_vMin, m_vMax, keep, m_labele, m_labeld,
224 IonTransverseDiffusion, xaxis, e, b, a);
225 keep = true;
226 AddFunction(m_eMin, m_eMax, m_vMin, m_vMax, keep, m_labele, m_labeld,
227 IonLongitudinalDiffusion, xaxis, e, b, a);
228 m_canvas->Update();
229}

◆ PlotIonVelocity()

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

Definition at line 181 of file ViewMedium.cc.

182 {
183
184 SetupCanvas();
185 bool keep = false;
186 AddFunction(m_eMin, m_eMax, m_vMin, m_vMax, keep, m_labele, m_labelv,
187 IonVelocity, xaxis, e, b, a);
188 m_canvas->Update();
189}

◆ SetBAngleRange()

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

Set the limits of the angle between electric and magnetic field.

Definition at line 95 of file ViewMedium.cc.

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

◆ SetCanvas()

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

Set the canvas to be painted on.

Definition at line 49 of file ViewMedium.cc.

49 {
50
51 if (!c) {
52 std::cerr << m_className << "::SetCanvas: Null pointer.\n";
53 return;
54 }
55 if (!m_hasExternalCanvas && m_canvas) {
56 delete m_canvas;
57 m_canvas = NULL;
58 }
59 m_canvas = c;
60 m_hasExternalCanvas = true;
61}

◆ SetElectricFieldRange()

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

Set the limits of the electric field.

Definition at line 73 of file ViewMedium.cc.

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

◆ SetFunctionRange() [1/2]

void Garfield::ViewMedium::SetFunctionRange ( )

Use the default function range.

Definition at line 117 of file ViewMedium.cc.

117{ m_vMin = m_vMax = 0.; }

◆ SetFunctionRange() [2/2]

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

Set the range of the function (velocity etc.) to be plotted.

Definition at line 106 of file ViewMedium.cc.

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

◆ SetMagneticFieldRange()

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

Set the limits of the magnetic field.

Definition at line 84 of file ViewMedium.cc.

84 {
85
86 if (bmin >= bmax || bmin < 0.) {
87 std::cerr << m_className << "::SetMagneticFieldRange: Incorrect 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 63 of file ViewMedium.cc.

63 {
64
65 if (!m) {
66 std::cerr << m_className << "::SetMedium: Null pointer.\n";
67 return;
68 }
69
70 m_medium = m;
71}

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