]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - ITS/AliITSresponse.h
Remove warnings with gcc4.0 (F.Carminati)
[u/mrichter/AliRoot.git] / ITS / AliITSresponse.h
... / ...
CommitLineData
1#ifndef ALIITSRESPONSE_H
2#define ALIITSRESPONSE_H
3/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * See cxx source for full Copyright notice */
5
6/* $Id$ */
7
8#include <TObject.h>
9#include <TString.h>
10
11class AliITSsegmentation;
12class TF1;
13class AliITSgeom;
14
15//----------------------------------------------
16//
17// ITS response virtual base class
18//
19class AliITSresponse : public TObject {
20 public:
21 // Default Constructor
22 AliITSresponse();
23 // Standard Constructor
24 AliITSresponse(Double_t Thickness);
25 // Destructor.
26 virtual ~AliITSresponse() {}
27 //
28 // Configuration methods
29 //
30 // fGeVcharge is set by default 3.6e-9 GeV See for ex. PDG 2004.
31 virtual void SetGeVToCharge(Double_t gc=3.6e-9){fGeVcharge = gc;}
32 // Returns the value fGeVcharge
33 virtual Double_t GetGeVToCharge() const {return fGeVcharge;}
34 // Converts deposited energy to number of electrons liberated
35 virtual Double_t GeVToCharge(Double_t gev) const {return gev/fGeVcharge;}
36 // Diffusion coefficient
37 virtual void SetDiffCoeff(Double_t, Double_t) = 0;
38 // Get diffusion coefficients
39 virtual void DiffCoeff(Double_t &,Double_t &) const = 0;
40
41 // Temperature in [degree K]
42 virtual void SetTemperature(Double_t t=300.0) {fT = t;}
43 // Get temperature [degree K]
44 virtual Double_t Temperature() const {return fT;}
45 // Set the impurity concentrations in [#/cm^3]
46 virtual void SetImpurity(Double_t n=0.0){fN = n;}
47 // Returns the impurity consentration in [#/cm^3]
48 virtual Double_t Impurity() const {return fN;}
49 // Sets the applied ratio distance/voltage [cm/volt]
50 virtual void SetDistanceOverVoltage(Double_t d,Double_t v){fdv = d/v;}
51 // Sets the applied ration distance/voltage [cm/volt]. Default value
52 // is 300E-4cm/80 volts = 0.000375 cm/volts
53 virtual void SetDistanceOverVoltage(Double_t dv=0.000375){fdv = dv;}
54 // Returns the ration distance/voltage
55 virtual Double_t DistanceOverVoltage() const {return fdv;}
56
57 // Get data type
58 virtual const char *DataType() const {return fDataType.Data();}
59 // Type of data - real or simulated
60 virtual void SetDataType(const char *data="simulated") {fDataType=data;}
61 // Set parameters options: "same" or read from "file" or "SetInvalid" or...
62 virtual void SetParamOptions(const char*,const char*) = 0;
63 // Set noise parameters
64 virtual void SetNoiseParam(Double_t, Double_t) = 0;
65 // Number of parameters to be set
66 virtual void SetNDetParam(Int_t) = 0;
67 // Set detector parameters: gain, coupling ...
68 virtual void SetDetParam(Double_t *) = 0;
69
70 // Parameters options
71 virtual void ParamOptions(char *,char*) const = 0;
72 virtual Int_t NDetParam() const = 0;
73 virtual void GetDetParam(Double_t *) const = 0;
74 virtual void GetNoiseParam(Double_t&, Double_t&) const = 0;
75
76 // Zero-suppression option - could be 1D, 2D or non-ZeroSuppressed
77 virtual void SetZeroSupp(const char*) = 0;
78 // Get zero-suppression option
79 virtual const char *ZeroSuppOption() const = 0;
80 // Set thresholds
81 virtual void SetThresholds(Double_t, Double_t) = 0;
82 virtual void Thresholds(Double_t &, Double_t &) const = 0;
83
84 // Set filenames
85 virtual void SetFilenames(const char *f1="",const char *f2="",
86 const char *f3=""){
87 // Set filenames - input, output, parameters ....
88 fFileName1=f1; fFileName2=f2; fFileName3=f3;}
89 // Filenames
90 virtual void Filenames(char* input,char* baseline,char* param) {
91 strcpy(input,fFileName1.Data()); strcpy(baseline,fFileName2.Data());
92 strcpy(param,fFileName3.Data());}
93
94 virtual Double_t DriftSpeed() const {return SpeedElectron();};
95 // set output option
96 virtual void SetOutputOption(Bool_t write=kFALSE) {fWrite = write;}
97
98 virtual Bool_t OutputOption() const {return fWrite;}
99 virtual Bool_t Do10to8() const {return kTRUE;}
100 virtual void GiveCompressParam(Int_t *) const =0;
101 //
102 // Detector type response methods
103 // Set number of sigmas over which cluster disintegration is performed
104 virtual void SetNSigmaIntegration(Double_t) = 0;
105 // Get number of sigmas over which cluster disintegration is performed
106 virtual Double_t NSigmaIntegration() const = 0;
107 // Set number of bins for the gaussian lookup table
108 virtual void SetNLookUp(Int_t) = 0;
109 // Get number of bins for the gaussian lookup table
110 virtual Int_t GausNLookUp() const {return 0;}
111 // Get scaling factor for bin i-th from the gaussian lookup table
112 virtual Double_t GausLookUp(Int_t) const {return 0.;}
113 // Set sigmas of the charge spread function
114 virtual void SetSigmaSpread(Double_t, Double_t) = 0;
115 // Get sigmas for the charge spread
116 virtual void SigmaSpread(Double_t &,Double_t &) const = 0;
117 // Pulse height from scored quantity (eloss)
118 virtual Double_t IntPH(Double_t) const {return 0.;}
119 // Charge disintegration
120 virtual Double_t IntXZ(AliITSsegmentation *) const {return 0.;}
121 // Electron mobility in Si. [cm^2/(Volt Sec)]. T in degree K, N in #/cm^3
122 virtual Double_t MobilityElectronSiEmp() const ;
123 // Hole mobility in Si. [cm^2/(Volt Sec)] T in degree K, N in #/cm^3
124 virtual Double_t MobilityHoleSiEmp() const ;
125 // Einstein relation for Diffusion Coefficient of Electrons. [cm^2/sec]
126 // T in degree K, N in #/cm^3
127 virtual Double_t DiffusionCoefficientElectron() const ;
128 // Einstein relation for Diffusion Coefficient of Holes. [cm^2/sec]
129 // T in [degree K], N in [#/cm^3]
130 virtual Double_t DiffusionCoefficientHole() const ;
131 // Electron <speed> under an applied electric field E=Volts/cm. [cm/sec]
132 // d distance-thickness in [cm], v in [volts], T in [degree K],
133 // N in [#/cm^3]
134 virtual Double_t SpeedElectron() const ;
135 // Holes <speed> under an applied electric field E=Volts/cm. [cm/sec]
136 // d distance-thickness in [cm], v in [volts], T in [degree K],
137 // N in [#/cm^3]
138 virtual Double_t SpeedHole() const ;
139 // Returns the Gaussian sigma == <x^2+z^2> [cm^2] due to the defusion of
140 // electrons or holes through a distance l [cm] caused by an applied
141 // voltage v [volt] through a distance d [cm] in any material at a
142 // temperature T [degree K].
143 virtual Double_t SigmaDiffusion3D(Double_t l) const;
144 // Returns the Gaussian sigma == <x^2 +y^2+z^2> [cm^2] due to the
145 // defusion of electrons or holes through a distance l [cm] caused by an
146 // applied voltage v [volt] through a distance d [cm] in any material at a
147 // temperature T [degree K].
148 virtual Double_t SigmaDiffusion2D(Double_t l) const;
149 // Returns the Gaussian sigma == <x^2+z^2> [cm^2] due to the defusion of
150 // electrons or holes through a distance l [cm] caused by an applied
151 // voltage v [volt] through a distance d [cm] in any material at a
152 // temperature T [degree K].
153 virtual Double_t SigmaDiffusion1D(Double_t l) const;
154 // Compute the thickness of the depleted region in a Si detector, version A
155 virtual Double_t DepletedRegionThicknessA(Double_t dopCons,
156 Double_t voltage,
157 Double_t elecCharge,
158 Double_t voltBuiltIn=0.5)const;
159 // Compute the thickness of the depleted region in a Si detector, version B
160 virtual Double_t DepletedRegionThicknessB(Double_t resist,Double_t voltage,
161 Double_t mobility,
162 Double_t voltBuiltIn=0.5,
163 Double_t dielConst=1.E-12)const;
164 // Computes the temperature dependance of the reverse bias current
165 virtual Double_t ReverseBiasCurrent(Double_t temp,Double_t revBiasCurT1,
166 Double_t tempT1,Double_t energy=1.2)const;
167 // Prints out the content of this class in ASCII format.
168 virtual void Print(ostream *os) const;
169 // Reads in the content of this class in the format of Print
170 virtual void Read(istream *is);
171 virtual void Print(Option_t *option="") const {TObject::Print(option);}
172 virtual Int_t Read(const char *name) {return TObject::Read(name);}
173 protected:
174 void NotImplemented(const char *method) const {if(gDebug>0)
175 Warning(method,"This method is not implemented for this sub-class");}
176
177 TString fDataType; // data type - real or simulated
178 private:
179 Double_t fdv; // The parameter d/v where d is the disance over which the
180 // the potential v is applied d/v [cm/volts]
181 Double_t fN; // the impurity consentration of the material in #/cm^3
182 Double_t fT; // The temperature of the Si in Degree K.
183 Double_t fGeVcharge; // Energy to ionize (free an electron) in GeV
184 TString fFileName1; // input keys : run, module #
185 TString fFileName2; // baseline & noise val or output code
186 // signal or monitored bgr.
187 TString fFileName3; // param values or output coded signal
188 Bool_t fWrite; // Write option for the compression algorithms
189
190 ClassDef(AliITSresponse,2) // Detector type response virtual base class
191};
192// Input and output function for standard C++ input/output.
193ostream& operator<<(ostream &os,AliITSresponse &source);
194istream& operator>>(istream &os,AliITSresponse &source);
195#endif