PandaRoot
PndSttSingleStraw.h
Go to the documentation of this file.
1 #ifndef __PNDSTTSINGLESTRAW__
2 #define __PNDSTTSINGLESTRAW__
3 
4 #include "TNamed.h"
5 #include "TMatrix.h"
6 #include "TVector3.h"
7 #include "TH2.h"
8 #include <vector>
9 
10 #include "TMatrixD.h"
11 
12 // ClassImp(TMatrixDBase)
13 // class TMatrixD;;
14 
15 class TH2F;
16 
17 class PndSttSingleStraw : public TNamed {
18 
19  // --------------------------------------------------------------------
20  public:
22  virtual ~PndSttSingleStraw(){};
24 
25  // Get methods
26 
27  Double_t GetCDist(Int_t k) { return CDist[k]; };
28  Double_t GetCNele(Int_t k) { return CNele[k]; };
29  Double_t GetCNeleT(Int_t k) { return CNeleT[k]; };
30  Double_t GetTeleTime(Int_t k) { return TeleTime[k]; };
31  Double_t GetPulse(Int_t k) { return Pulse[k]; };
32  Double_t GetPulseT(Int_t k) { return PulseT[k]; };
33  Double_t GetWi() { return Wi * 1.e-09; }; // in GeV
34  Double_t GetGasGain() { return GasGain; };
35  Double_t GetXmax() { return Xmax; };
36  Double_t GetDx() { return Dx; };
37  Double_t GetEmed() { return Emed; };
38  Double_t GetEmin() { return Emin; };
39  Double_t GetNcl() { return Ncl; };
40  Double_t GetCsi() { return Csi; };
41  Double_t GetDelta() { return Delta; };
42  Double_t GetEmax() { return Emax; };
43  Double_t GetIAr() { return IAr; };
44  Int_t GetNNClus() { return NNClus; };
45  Int_t GetNchann() { return Nchann; };
46  Double_t GetPulseMax() { return PulseMax; };
47  Double_t GetGamma() { return gamma; };
48  Double_t GetBeta() { return beta; };
49  Double_t GetpSTP() { return pSTP; };
50  Double_t GetPulseTime() { return PulseTime; };
51  Double_t GetbPolya() { return bPolya; };
52  Double_t GetSigUrb() { return SigUrb; };
53  Double_t GetNUrban() { return NUrban; };
54  Double_t GetEup() { return Eup; };
55  Double_t GetAvUrb() { return AvUrb; };
56 
57  // coordinates
58  Double_t GetXin() { return Xin; };
59  Double_t GetYin() { return Yin; };
60  Double_t GetZin() { return Zin; };
61  Double_t GetXout() { return Xout; };
62  Double_t GetYout() { return Yout; };
63  Double_t GetZout() { return Zout; };
64  Double_t GetRpath() { return Rpath; };
65 
66  TVector3 GetWDist(Int_t k) { return WDist[k]; };
67 
68  // put methods
69  void PutTrackXYZ(Double_t v1, Double_t v2, Double_t v3, Double_t v4, Double_t v5, Double_t v6)
70  {
71  Xin = v1;
72  Yin = v2;
73  Zin = v3;
74  Xout = v4;
75  Yout = v5;
76  Zout = v6;
77  };
78  void PutRpath(Double_t value) { Rpath = value; };
79 
80  void PutPolya(Double_t par) { bPolya = par; };
81 
82  void PutRadius(Double_t value) { Radius = value; };
83  void PutPress(Double_t value) { pSTP = value; };
84 
85  // ----------------------------------------------------------------------------
86  // calls at each track for MC applications
87 
88  // define the wire
89  void PutWireXYZ(Double_t w1, Double_t w2, Double_t w3, Double_t w4, Double_t w5, Double_t w6);
90  // straw time
91  Double_t PartToTime(Double_t Mass, Double_t Momentum, Double_t InOut[]);
92  // ADC signal corresponding to the energy loss of StrawCharge
93  Double_t PartToADC();
94 
95  // fast reconstructed distance in cm
96  Double_t FastRec(Double_t TrueDcm, Int_t Flag);
97 
98  // ADC signal charge fast simulation
99  Double_t FastPartToADC();
100 
101  // ----------------------------------------------------------------------------
102  // standard calls
103 
104  // once to set constants
105 
106  void TConst(Double_t Radius, Double_t pSTP, Double_t ArP, Double_t CO2P);
107 
108  // to define the track
109  void TInit(Double_t Mass, Double_t Momentum, Double_t InOut[]); // for each track
110 
111  // other possible calls
112 
113  Double_t StrawCharge(); // total discharge electron calculation
114  // energy loss in GeV
115  Int_t StrawSignal(Int_t nsteps); // oscilloscope signal generation (ns)
116  // Pulse[PulseT] in ns
117  // return number of primary electrons
118  Int_t StrawTime(); // output time of the straw (ns) &
119  // PulseTime[Int_t] in ns
120  Int_t StrawTot(); // time window over threshold of the signal
121  Double_t TimnsToDiscm(Double_t time); // from time (ns) to radius in cm
122 
123  //-----------------------------------------------------------------------------
124  // utility methods
125 
126  void TDirCos(); // track director cosines
127  Double_t TrueDist(Double_t Point[]); // true distance wire-track
128  Double_t DiffLong(Double_t distcm); // longituinal diffusion in microns
129  Double_t DiffTran(Double_t distcm); // transverse diffusion in microns
130  void Polya(Double_t bpar); // Polya cumulative calculation
131  Double_t PolyaSamp(); // sampling from Polya distribution
132  Double_t RRise(Double_t gamma); // relativistic rise calculation
133  Int_t Cluster(); // cluster generation
134  Int_t Eject(); // primary electron generation
135  TVector3 WDistCalc(Double_t d); // distance electron-wire
136  Double_t Signal(Double_t t, Double_t t0); // signal functional form
137  Double_t STEloss(); // Landau energy loss
138  Double_t STUrban(); // Urban energy loss
139  Int_t TimeEle(); // arrivals times of all the electrons
140  Double_t DistEle(Double_t tns); // distance of all the electrons
141  //------------------------------------------------------------------------------
142 
143  private:
144  //-----------------------------------------------------------------
145 
146  // distance, number of electrons (with delta rays),
147  // distance from wire of the cluster
148  std::vector<Double_t> CDist, CDistC, CNele, CNeleT, TeleTime, AmplSig;
149  std::vector<Double_t> Pulse, PulseT;
150  std::vector<TVector3> WDist;
151 
152  // set constants
153  // masses in GeV, energies in GeV, cgs system
154 
155  // Input for the medium
156 
157  // Double_t PClus[20], CO2Clus[20], CumClus[21], CH4Clus[20];
158  Double_t CumClus[21], CH4Clus[20]; // shadows deleted
159 
160  Double_t Wi;
161  Double_t ArPerc, CO2Perc, CH4Perc; // volume percentages
162  Double_t ArWPerc, CO2WPerc, CH4WPerc; // weight percentages
163  Double_t pSTP; // pressure (STP reference)
164  Double_t Radius; // straw radius
165  Int_t Field; // flag for magnetic field --> if =1, magnetic field on
166  Double_t AAr; // Argon
167  Double_t ZAr; // Argon
168  Double_t RhoAr; // g/cm3 (1.78 mg/cm3)
169  Double_t NclAr; // clusters/cm
170 
171  Double_t EmedAr;
172  Double_t EminAr;
173  Double_t EmpAr;
174  Double_t CsiAr;
175  Double_t IAr; // ionization potential (188 eV)
176  Double_t WiAr; // energy to reate an ion pair in Argon
177  Double_t Ncl; // mean number of cluster in the mixture
178  Double_t Ecl; // electron per cluster
179  Double_t Lcl; // mean free path between clusters
180  Double_t Ntote; // mean total number of eletrons
181  Double_t GasGain; // gain of the gas
182  Double_t Cutoff; // limit the number of primry electrons
183  // CO2
184  Double_t EmedCO2;
185  Double_t EminCO2;
186  Double_t EmpCO2;
187  Double_t CsiCO2;
188  Double_t ACO2; // CO2 (39.948)
189  Double_t ZCO2; // CO2 (18)
190  Double_t RhoCO2; // g/cm3 CO2 (1.98 mg/cm3)
191  Double_t ICO2; // ionization potential (GeV) (188 eV)
192  Double_t WiCO2; // energy to create an ion pair
193  Double_t NclCO2; // clusters/cm
194 
195  // Methane CH4 ------------------------------------------------
196  Double_t EmedCH4;
197  Double_t EmpCH4;
198  Double_t CsiCH4;
199  Double_t EminCH4;
200  Double_t ACH4; // CO2 (39.948)
201  Double_t ZCH4; // CO2 (18)
202  Double_t RhoCH4; // g/cm3 CO2 (0.71 mg/cm3)
203  Double_t ICH4; // ionization potential (GeV) (188 eV)
204  Double_t WiCH4; // energy to create an ion pair
205  Double_t NclCH4; // clusters/cm
206 
207  Double_t RhoMixCO2;
208  Double_t RhoMixCH4;
209  //----------------------------------------------------------------------
210  // Input for the particle (Gev, energy loses in Kev
211 
212  Double_t PZeta; // charge
213  Double_t piMass; // particle mass (GeV)
214  Double_t PMass; // incident particle mass
215  Double_t PMom; // particle momentum (GeV)
216  Double_t Dx; // distance travelled in gas (cm)
217  Double_t eMass; // electron mass (GeV) (0.511 MeV)
218  Double_t prMass; // proton mass
219  Double_t Delta; // polarization Sternheimer parameter
220  Int_t CNumb; // current number of clusters
221 
222  // ------------------------------------------------------
223  // quantities for each track
224  // calculated in TInit
225 
226  Double_t PEn; // particle energy GeV
227  Double_t beta;
228  Double_t gamma;
229  Double_t Emed; // GeV
230  Double_t Emin; // GeV/cm (2.7 keV)
231  Double_t Csi;
232  Double_t Emax;
233  Double_t Emp; // most probable energy
234  Int_t NNClus; // number of clusters
235 
236  // ---------------------------------------------------------------------
237  // mathematical and statistical parameters
238 
239  Double_t PolyaCum[100], Xs[100];
240  Double_t Xin, Yin, Zin, Xout, Yout, Zout, Rpath;
241  Int_t NPolya;
242  Double_t Xmax;
243  Double_t bPolya;
244  Double_t Calpha, Cbeta, Cgamma;
245  Double_t NUrban; // total number of collision in the Urban model
246  Double_t SigUrb; // std dev of the Urban distribution (routine STUrban)
247  Double_t Eup; // upper limit of the Urban distribution (routine STUrban)
248  Double_t AvUrb; // average of the Urban distribution (routine STUrban)
249 
250  // --------------------------------------------------------------
251  // for the straws
252 
253  Double_t Wx1, Wy1, Wz1, Wx2, Wy2, Wz2; // wire coordinates
254  Double_t Wp, Wq, Wr; // director cosine of the wire
255  Double_t PulseMax;
256  Double_t PulseTime, PulseTime1;
257  Double_t Thresh1; // first threshold
258  Double_t Thresh2;
259  Double_t OffT; // time offset
260  Int_t Nchann; // number of channels for the straw signal
261 
262  // ----------------------------------------------------------------------
263  // dummy
264 
265  Double_t Out1, Out4;
266  Int_t Out2, Out3;
267 };
268 
269 #endif
ClassDef(PndSttSingleStraw, 1)
virtual ~PndSttSingleStraw()
Double_t GetTeleTime(Int_t k)
void PutWireXYZ(Double_t w1, Double_t w2, Double_t w3, Double_t w4, Double_t w5, Double_t w6)
Double_t GetPulse(Int_t k)
void PutPress(Double_t value)
Double_t GetCDist(Int_t k)
TVector3 WDistCalc(Double_t d)
Double_t TrueDist(Double_t Point[])
Double_t GetCNele(Int_t k)
Double_t FastPartToADC()
TVector3 GetWDist(Int_t k)
Double_t PolyaSamp()
Int_t StrawSignal(Int_t nsteps)
Double_t Signal(Double_t t, Double_t t0)
Double_t DiffLong(Double_t distcm)
void PutRadius(Double_t value)
Double_t DistEle(Double_t tns)
void PutPolya(Double_t par)
Double_t DiffTran(Double_t distcm)
Double_t TimnsToDiscm(Double_t time)
void TInit(Double_t Mass, Double_t Momentum, Double_t InOut[])
Double_t GetCNeleT(Int_t k)
Double_t GetPulseT(Int_t k)
void Polya(Double_t bpar)
Double_t PartToADC()
Double_t StrawCharge()
Double_t FastRec(Double_t TrueDcm, Int_t Flag)
void TConst(Double_t Radius, Double_t pSTP, Double_t ArP, Double_t CO2P)
Double_t STEloss()
Double_t RRise(Double_t gamma)
Double_t PartToTime(Double_t Mass, Double_t Momentum, Double_t InOut[])
void PutTrackXYZ(Double_t v1, Double_t v2, Double_t v3, Double_t v4, Double_t v5, Double_t v6)
void PutRpath(Double_t value)
Double_t STUrban()