]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TRD/AliTRDsim.cxx
Access function to local momenta renamed.
[u/mrichter/AliRoot.git] / TRD / AliTRDsim.cxx
CommitLineData
46d29e70 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16/*
17$Log$
2ab0c725 18Revision 1.7 2000/12/20 13:00:45 cblume
19Modifications for the HP-compiler
20
d1b06c24 21Revision 1.6 2000/12/12 10:20:10 cblume
22Initialize fSepctrum = 0 in ctors
23
fa5e892a 24Revision 1.5 2000/10/15 23:40:01 cblume
25Remove AliTRDconst
26
0e9c2ad5 27Revision 1.4 2000/10/06 16:49:46 cblume
28Made Getters const
29
46d29e70 30Revision 1.3.2.1 2000/09/18 13:45:30 cblume
31New class AliTRDsim that simulates TR photons
32
33Revision 1.2 1999/09/29 09:24:35 fca
34Introduction of the Copyright and cvs Log
35
36*/
37
38///////////////////////////////////////////////////////////////////////////////
39// //
40// TRD simulation - multimodule (regular rad.) //
41// after: M. CASTELLANO et al., COMP. PHYS. COMM. 51 (1988) 431 //
42// + COMP. PHYS. COMM. 61 (1990) 395 //
43// //
44// 17.07.1998 - A.Andronic //
45// 08.12.1998 - simplified version //
46// 11.07.2000 - Adapted code to aliroot environment (C.Blume) //
47// //
48///////////////////////////////////////////////////////////////////////////////
49
50#include <stdlib.h>
51
0e9c2ad5 52#include <TH1.h>
53#include <TRandom.h>
54#include <TMath.h>
55#include <TParticle.h>
46d29e70 56
46d29e70 57#include "AliModule.h"
58
0e9c2ad5 59#include "AliTRDsim.h"
60
46d29e70 61ClassImp(AliTRDsim)
62
63//_____________________________________________________________________________
64AliTRDsim::AliTRDsim():TObject()
65{
66 //
67 // AliTRDsim default constructor
68 //
69
fa5e892a 70 fSpectrum = 0;
d1b06c24 71 fSigma = 0;
fa5e892a 72
46d29e70 73 Init();
74
75}
76
77//_____________________________________________________________________________
78AliTRDsim::AliTRDsim(AliModule *mod, Int_t foil, Int_t gap)
79{
80 //
81 // AliTRDsim constructor. Takes the material properties of the radiator
82 // foils and the gas in the gaps from AliModule <mod>.
83 // The default number of foils is 100 with a thickness of 20 mu. The
84 // thickness of the gaps is 500 mu.
85 //
86
87 Float_t aFoil, zFoil, rhoFoil;
88 Float_t aGap, zGap, rhoGap;
89 Float_t rad, abs;
90 Char_t name[21];
91
fa5e892a 92 fSpectrum = 0;
d1b06c24 93 fSigma = 0;
fa5e892a 94
46d29e70 95 Init();
96
97 mod->AliGetMaterial(foil,name,aFoil,zFoil,rhoFoil,rad,abs);
98 mod->AliGetMaterial(gap ,name,aGap ,zGap ,rhoGap ,rad,abs);
99
100 fFoilDens = rhoFoil;
101 fFoilA = aFoil;
102 fFoilZ = zFoil;
103 fFoilOmega = Omega(fFoilDens,fFoilZ,fFoilA);
104
105 fGapDens = rhoGap;
106 fGapA = aGap;
107 fGapZ = zGap;
108 fGapOmega = Omega(fGapDens ,fGapZ ,fGapA );
109
110}
111
112//_____________________________________________________________________________
113AliTRDsim::AliTRDsim(const AliTRDsim &s)
114{
115 //
116 // AliTRDsim copy constructor
117 //
118
119 ((AliTRDsim &) s).Copy(*this);
120
121}
122
123//_____________________________________________________________________________
124AliTRDsim::~AliTRDsim()
125{
126 //
127 // AliTRDsim destructor
128 //
129
130 if (fSpectrum) delete fSpectrum;
131 if (fSigma) delete fSigma;
132
133}
134
135//_____________________________________________________________________________
136AliTRDsim &AliTRDsim::operator=(const AliTRDsim &s)
137{
138 //
139 // Assignment operator
140 //
141
142 if (this != &s) ((AliTRDsim &) s).Copy(*this);
143 return *this;
144
145}
146
147//_____________________________________________________________________________
148void AliTRDsim::Copy(TObject &s)
149{
150 //
151 // Copy function
152 //
153
154 ((AliTRDsim &) s).fNFoils = fNFoils;
155 ((AliTRDsim &) s).fFoilThick = fFoilThick;
156 ((AliTRDsim &) s).fFoilDens = fFoilDens;
157 ((AliTRDsim &) s).fFoilOmega = fFoilOmega;
158 ((AliTRDsim &) s).fFoilZ = fFoilZ;
159 ((AliTRDsim &) s).fFoilA = fFoilA;
160 ((AliTRDsim &) s).fGapThick = fGapThick;
161 ((AliTRDsim &) s).fGapDens = fGapDens;
162 ((AliTRDsim &) s).fGapOmega = fGapOmega;
163 ((AliTRDsim &) s).fGapZ = fGapZ;
164 ((AliTRDsim &) s).fGapA = fGapA;
165 ((AliTRDsim &) s).fTemp = fTemp;
166 ((AliTRDsim &) s).fSpNBins = fSpNBins;
167 ((AliTRDsim &) s).fSpRange = fSpRange;
168 ((AliTRDsim &) s).fSpBinWidth = fSpBinWidth;
169 ((AliTRDsim &) s).fSpLower = fSpLower;
170 ((AliTRDsim &) s).fSpUpper = fSpUpper;
171
172 if (((AliTRDsim &) s).fSigma) delete ((AliTRDsim &) s).fSigma;
173 ((AliTRDsim &) s).fSigma = new Double_t[fSpNBins];
174 for (Int_t iBin = 0; iBin < fSpNBins; iBin++) {
175 ((AliTRDsim &) s).fSigma[iBin] = fSigma[iBin];
176 }
177
178 fSpectrum->Copy(*((AliTRDsim &) s).fSpectrum);
179
180}
181
182//_____________________________________________________________________________
183void AliTRDsim::Init()
184{
185 //
186 // Initialization
187 // The default radiator are 100 prolypropilene foils of 20 mu thickness
188 // with gaps of 500 mu filled with CO2.
189 //
190 //
191
192 fNFoils = 100;
193
194 fFoilThick = 0.0020;
195 fFoilDens = 0.92;
196 fFoilZ = 5.28571;
197 fFoilA = 10.4286;
198 fFoilOmega = Omega(fFoilDens,fFoilZ,fFoilA);
199
200 fGapThick = 0.0500;
201 fGapDens = 0.001977;
202 fGapZ = 7.45455;
203 fGapA = 14.9091;
204 fGapOmega = Omega(fGapDens ,fGapZ ,fGapA );
205
206 fTemp = 293.16;
207
208 fSpNBins = 200;
209 fSpRange = 100;
210 fSpBinWidth = fSpRange / fSpNBins;
211 fSpLower = 1.0 - 0.5 * fSpBinWidth;
212 fSpUpper = fSpLower + fSpRange;
213
214 if (fSpectrum) delete fSpectrum;
215 fSpectrum = new TH1D("TRspectrum","TR spectrum",fSpNBins,fSpLower,fSpUpper);
216
217 // Set the sigma values
218 SetSigma();
219
220}
221
222//_____________________________________________________________________________
223Int_t AliTRDsim::CreatePhotons(Int_t pdg, Float_t p
224 , Int_t &nPhoton, Float_t *ePhoton)
225{
226 //
227 // Create TRD photons for a charged particle of type <pdg> with the total
228 // momentum <p>.
229 // Number of produced TR photons: <nPhoton>
230 // Energies of the produced TR photons: <ePhoton>
231 //
232
233 // PDG codes
234 const Int_t kPdgEle = 11;
235 const Int_t kPdgMuon = 13;
236 const Int_t kPdgPion = 211;
237 const Int_t kPdgKaon = 321;
238
239 Float_t mass = 0;
240 switch (TMath::Abs(pdg)) {
241 case kPdgEle:
242 mass = 5.11e-4;
243 break;
244 case kPdgMuon:
245 mass = 0.10566;
246 break;
247 case kPdgPion:
248 mass = 0.13957;
249 break;
250 case kPdgKaon:
251 mass = 0.4937;
252 break;
253 default:
254 return 0;
255 break;
256 };
257
258 // Calculate gamma
259 Double_t gamma = TMath::Sqrt(p*p + mass*mass) / mass;
260
261 // Calculate the TR photons
262 return TrPhotons(gamma, nPhoton, ePhoton);
263
264}
265
266//_____________________________________________________________________________
267Int_t AliTRDsim::TrPhotons(Double_t gamma, Int_t &nPhoton, Float_t *ePhoton)
268{
269 //
270 // Produces TR photons.
271 //
272
273 const Double_t kAlpha = 0.0072973;
274 const Int_t kSumMax = 10;
275
276 Double_t kappa = fGapThick / fFoilThick;
277
278 fSpectrum->Reset();
279
280 // The TR spectrum
281 Double_t stemp = 0;
282 for (Int_t iBin = 0; iBin < fSpNBins; iBin++) {
283
284 // keV -> eV
285 Double_t energyeV = (fSpBinWidth * iBin + 1.0) * 1e3;
286
287 Double_t csFoil = fFoilOmega / energyeV;
288 Double_t csGap = fGapOmega / energyeV;
289
290 Double_t rho1 = energyeV * fFoilThick * 1e4 * 2.5
291 * (1.0 / (gamma*gamma) + csFoil*csFoil);
292 Double_t rho2 = energyeV * fFoilThick * 1e4 * 2.5
293 * (1.0 / (gamma*gamma) + csGap *csGap);
294
295 // Calculate the sum
296 Double_t sum = 0;
297 for (Int_t iSum = 0; iSum < kSumMax; iSum++) {
298 Double_t tetan = (TMath::Pi() * 2.0 * (iSum+1) - (rho1 + kappa * rho2))
299 / (kappa + 1.0);
300 if (tetan < 0.0) tetan = 0.0;
301 Double_t aux = 1.0 / (rho1 + tetan) - 1.0 / (rho2 + tetan);
302 sum += tetan * (aux*aux) * (1.0 - TMath::Cos(rho1 + tetan));
303 }
304
305 // Absorbtion
306 Double_t conv = 1.0 - TMath::Exp(-fNFoils * fSigma[iBin]);
307
308 // eV -> keV
309 Float_t energykeV = energyeV * 0.001;
310
311 // dN / domega
312 Double_t wn = kAlpha * 4.0 / (fSigma[iBin] * (kappa + 1.0))
313 * conv * sum / energykeV;
314 fSpectrum->SetBinContent(iBin,wn);
315
316 stemp += wn;
317
318 }
319
320 // <nTR> (binsize corr.)
321 Float_t ntr = stemp * fSpBinWidth;
322 // Number of TR photons from Poisson distribution with mean <ntr>
323 nPhoton = gRandom->Poisson(ntr);
324 // Energy of the TR photons
325 for (Int_t iPhoton = 0; iPhoton < nPhoton; iPhoton++) {
326 ePhoton[iPhoton] = fSpectrum->GetRandom();
327 }
328
329 return 1;
330
331}
332
333//_____________________________________________________________________________
334void AliTRDsim::SetSigma()
335{
336 //
337 // Sets the absorbtion crosssection for the energies of the TR spectrum
338 //
339
340 if (fSigma) delete fSigma;
341 fSigma = new Double_t[fSpNBins];
342 for (Int_t iBin = 0; iBin < fSpNBins; iBin++) {
343 Double_t energykeV = iBin * fSpBinWidth + 1.0;
344 fSigma[iBin] = Sigma(energykeV);
345 //printf("SetSigma(): iBin = %d fSigma %g\n",iBin,fSigma[iBin]);
346 }
347
348}
349
350//_____________________________________________________________________________
351Double_t AliTRDsim::Sigma(Double_t energykeV)
352{
353 //
354 // Calculates the absorbtion crosssection for a one-foil-one-gap-radiator
355 //
356
357 // Gas at 0 C
358 const Double_t kTemp0 = 273.16;
359
360 // keV -> MeV
361 Double_t energyMeV = energykeV * 0.001;
362 if (energyMeV >= 0.001) {
363 return(GetMuPo(energyMeV) * fFoilDens * fFoilThick +
364 GetMuCO(energyMeV) * fGapDens * fGapThick * fTemp/kTemp0);
365 }
366 else {
367 return 1e6;
368 }
369
370}
371
372//_____________________________________________________________________________
373Double_t AliTRDsim::GetMuPo(Double_t energyMeV)
374{
375 //
376 // Returns the photon absorbtion cross section for polypropylene
377 //
378
379 const Int_t kN = 36;
380
381 Double_t mu[kN] = { 1.894E+03, 5.999E+02, 2.593E+02
382 , 7.743E+01, 3.242E+01, 1.643E+01
383 , 9.432E+00, 3.975E+00, 2.088E+00
384 , 7.452E-01, 4.315E-01, 2.706E-01
385 , 2.275E-01, 2.084E-01, 1.970E-01
386 , 1.823E-01, 1.719E-01, 1.534E-01
387 , 1.402E-01, 1.217E-01, 1.089E-01
388 , 9.947E-02, 9.198E-02, 8.078E-02
389 , 7.262E-02, 6.495E-02, 5.910E-02
390 , 5.064E-02, 4.045E-02, 3.444E-02
391 , 3.045E-02, 2.760E-02, 2.383E-02
392 , 2.145E-02, 1.819E-02, 1.658E-02 };
393
394 Double_t en[kN] = { 1.000E-03, 1.500E-03, 2.000E-03
395 , 3.000E-03, 4.000E-03, 5.000E-03
396 , 6.000E-03, 8.000E-03, 1.000E-02
397 , 1.500E-02, 2.000E-02, 3.000E-02
398 , 4.000E-02, 5.000E-02, 6.000E-02
399 , 8.000E-02, 1.000E-01, 1.500E-01
400 , 2.000E-01, 3.000E-01, 4.000E-01
401 , 5.000E-01, 6.000E-01, 8.000E-01
402 , 1.000E+00, 1.250E+00, 1.500E+00
403 , 2.000E+00, 3.000E+00, 4.000E+00
404 , 5.000E+00, 6.000E+00, 8.000E+00
405 , 1.000E+01, 1.500E+01, 2.000E+01 };
406
407 return Interpolate(energyMeV,en,mu,kN);
408
409}
410
411//_____________________________________________________________________________
412Double_t AliTRDsim::GetMuCO(Double_t energyMeV)
413{
414 //
415 // Returns the photon absorbtion cross section for CO2
416 //
417
418 const Int_t kN = 36;
419
420 Double_t mu[kN] = { 0.39383E+04, 0.13166E+04, 0.58750E+03
421 , 0.18240E+03, 0.77996E+02, 0.40024E+02
422 , 0.23116E+02, 0.96997E+01, 0.49726E+01
423 , 0.15543E+01, 0.74915E+00, 0.34442E+00
424 , 0.24440E+00, 0.20589E+00, 0.18632E+00
425 , 0.16578E+00, 0.15394E+00, 0.13558E+00
426 , 0.12336E+00, 0.10678E+00, 0.95510E-01
427 , 0.87165E-01, 0.80587E-01, 0.70769E-01
428 , 0.63626E-01, 0.56894E-01, 0.51782E-01
429 , 0.44499E-01, 0.35839E-01, 0.30825E-01
430 , 0.27555E-01, 0.25269E-01, 0.22311E-01
431 , 0.20516E-01, 0.18184E-01, 0.17152E-01 };
432
433 Double_t en[kN] = { 0.10000E-02, 0.15000E-02, 0.20000E-02
434 , 0.30000E-02, 0.40000E-02, 0.50000E-02
435 , 0.60000E-02, 0.80000E-02, 0.10000E-01
436 , 0.15000E-01, 0.20000E-01, 0.30000E-01
437 , 0.40000E-01, 0.50000E-01, 0.60000E-01
438 , 0.80000E-01, 0.10000E+00, 0.15000E+00
439 , 0.20000E+00, 0.30000E+00, 0.40000E+00
440 , 0.50000E+00, 0.60000E+00, 0.80000E+00
441 , 0.10000E+01, 0.12500E+01, 0.15000E+01
442 , 0.20000E+01, 0.30000E+01, 0.40000E+01
443 , 0.50000E+01, 0.60000E+01, 0.80000E+01
444 , 0.10000E+02, 0.15000E+02, 0.20000E+02 };
445
446 return Interpolate(energyMeV,en,mu,kN);
447
448}
449
450//_____________________________________________________________________________
451Double_t AliTRDsim::GetMuXe(Double_t energyMeV)
452{
453 //
454 // Returns the photon absorbtion cross section for xenon
455 //
456
457 const Int_t kN = 48;
458
459 Double_t mu[kN] = { 9.413E+03, 8.151E+03, 7.035E+03
460 , 7.338E+03, 4.085E+03, 2.088E+03
461 , 7.780E+02, 3.787E+02, 2.408E+02
462 , 6.941E+02, 6.392E+02, 6.044E+02
463 , 8.181E+02, 7.579E+02, 6.991E+02
464 , 8.064E+02, 6.376E+02, 3.032E+02
465 , 1.690E+02, 5.743E+01, 2.652E+01
466 , 8.930E+00, 6.129E+00, 3.316E+01
467 , 2.270E+01, 1.272E+01, 7.825E+00
468 , 3.633E+00, 2.011E+00, 7.202E-01
469 , 3.760E-01, 1.797E-01, 1.223E-01
470 , 9.699E-02, 8.281E-02, 6.696E-02
471 , 5.785E-02, 5.054E-02, 4.594E-02
472 , 4.078E-02, 3.681E-02, 3.577E-02
473 , 3.583E-02, 3.634E-02, 3.797E-02
474 , 3.987E-02, 4.445E-02, 4.815E-02 };
475
476 Double_t en[kN] = { 1.00000E-03, 1.07191E-03, 1.14900E-03
477 , 1.14900E-03, 1.50000E-03, 2.00000E-03
478 , 3.00000E-03, 4.00000E-03, 4.78220E-03
479 , 4.78220E-03, 5.00000E-03, 5.10370E-03
480 , 5.10370E-03, 5.27536E-03, 5.45280E-03
481 , 5.45280E-03, 6.00000E-03, 8.00000E-03
482 , 1.00000E-02, 1.50000E-02, 2.00000E-02
483 , 3.00000E-02, 3.45614E-02, 3.45614E-02
484 , 4.00000E-02, 5.00000E-02, 6.00000E-02
485 , 8.00000E-02, 1.00000E-01, 1.50000E-01
486 , 2.00000E-01, 3.00000E-01, 4.00000E-01
487 , 5.00000E-01, 6.00000E-01, 8.00000E-01
488 , 1.00000E+00, 1.25000E+00, 1.50000E+00
489 , 2.00000E+00, 3.00000E+00, 4.00000E+00
490 , 5.00000E+00, 6.00000E+00, 8.00000E+00
491 , 1.00000E+01, 1.50000E+01, 2.00000E+01 };
492
493 return Interpolate(energyMeV,en,mu,kN);
494
495}
496
497//_____________________________________________________________________________
498Double_t AliTRDsim::GetMuBu(Double_t energyMeV)
499{
500 //
501 // Returns the photon absorbtion cross section for isobutane
502 //
503
504 const Int_t kN = 36;
505
506 Double_t mu[kN] = { 0.38846E+03, 0.12291E+03, 0.53225E+02
507 , 0.16091E+02, 0.69114E+01, 0.36541E+01
508 , 0.22282E+01, 0.11149E+01, 0.72887E+00
509 , 0.45053E+00, 0.38167E+00, 0.33920E+00
510 , 0.32155E+00, 0.30949E+00, 0.29960E+00
511 , 0.28317E+00, 0.26937E+00, 0.24228E+00
512 , 0.22190E+00, 0.19289E+00, 0.17288E+00
513 , 0.15789E+00, 0.14602E+00, 0.12829E+00
514 , 0.11533E+00, 0.10310E+00, 0.93790E-01
515 , 0.80117E-01, 0.63330E-01, 0.53229E-01
516 , 0.46390E-01, 0.41425E-01, 0.34668E-01
517 , 0.30267E-01, 0.23910E-01, 0.20509E-01 };
518
519 Double_t en[kN] = { 0.10000E-02, 0.15000E-02, 0.20000E-02
520 , 0.30000E-02, 0.40000E-02, 0.50000E-02
521 , 0.60000E-02, 0.80000E-02, 0.10000E-01
522 , 0.15000E-01, 0.20000E-01, 0.30000E-01
523 , 0.40000E-01, 0.50000E-01, 0.60000E-01
524 , 0.80000E-01, 0.10000E+00, 0.15000E+00
525 , 0.20000E+00, 0.30000E+00, 0.40000E+00
526 , 0.50000E+00, 0.60000E+00, 0.80000E+00
527 , 0.10000E+01, 0.12500E+01, 0.15000E+01
528 , 0.20000E+01, 0.30000E+01, 0.40000E+01
529 , 0.50000E+01, 0.60000E+01, 0.80000E+01
530 , 0.10000E+02, 0.15000E+02, 0.20000E+02 };
531
532 return Interpolate(energyMeV,en,mu,kN);
533
534}
535
536//_____________________________________________________________________________
537Double_t AliTRDsim::GetMuMy(Double_t energyMeV)
538{
539 //
540 // Returns the photon absorbtion cross section for mylar
541 //
542
543 const Int_t kN = 36;
544
545 Double_t mu[kN] = { 2.911E+03, 9.536E+02, 4.206E+02
546 , 1.288E+02, 5.466E+01, 2.792E+01
547 , 1.608E+01, 6.750E+00, 3.481E+00
548 , 1.132E+00, 5.798E-01, 3.009E-01
549 , 2.304E-01, 2.020E-01, 1.868E-01
550 , 1.695E-01, 1.586E-01, 1.406E-01
551 , 1.282E-01, 1.111E-01, 9.947E-02
552 , 9.079E-02, 8.395E-02, 7.372E-02
553 , 6.628E-02, 5.927E-02, 5.395E-02
554 , 4.630E-02, 3.715E-02, 3.181E-02
555 , 2.829E-02, 2.582E-02, 2.257E-02
556 , 2.057E-02, 1.789E-02, 1.664E-02 };
557
558 Double_t en[kN] = { 1.00000E-03, 1.50000E-03, 2.00000E-03
559 , 3.00000E-03, 4.00000E-03, 5.00000E-03
560 , 6.00000E-03, 8.00000E-03, 1.00000E-02
561 , 1.50000E-02, 2.00000E-02, 3.00000E-02
562 , 4.00000E-02, 5.00000E-02, 6.00000E-02
563 , 8.00000E-02, 1.00000E-01, 1.50000E-01
564 , 2.00000E-01, 3.00000E-01, 4.00000E-01
565 , 5.00000E-01, 6.00000E-01, 8.00000E-01
566 , 1.00000E+00, 1.25000E+00, 1.50000E+00
567 , 2.00000E+00, 3.00000E+00, 4.00000E+00
568 , 5.00000E+00, 6.00000E+00, 8.00000E+00
569 , 1.00000E+01, 1.50000E+01, 2.00000E+01 };
570
571 return Interpolate(energyMeV,en,mu,kN);
572
573}
574
575//_____________________________________________________________________________
576Double_t AliTRDsim::GetMuN2(Double_t energyMeV)
577{
578 //
579 // Returns the photon absorbtion cross section for nitrogen
580 //
581
582 const Int_t kN = 36;
583
584 Double_t mu[kN] = { 3.311E+03, 1.083E+03, 4.769E+02
585 , 1.456E+02, 6.166E+01, 3.144E+01
586 , 1.809E+01, 7.562E+00, 3.879E+00
587 , 1.236E+00, 6.178E-01, 3.066E-01
588 , 2.288E-01, 1.980E-01, 1.817E-01
589 , 1.639E-01, 1.529E-01, 1.353E-01
590 , 1.233E-01, 1.068E-01, 9.557E-02
591 , 8.719E-02, 8.063E-02, 7.081E-02
592 , 6.364E-02, 5.693E-02, 5.180E-02
593 , 4.450E-02, 3.579E-02, 3.073E-02
594 , 2.742E-02, 2.511E-02, 2.209E-02
595 , 2.024E-02, 1.782E-02, 1.673E-02 };
596
597 Double_t en[kN] = { 1.00000E-03, 1.50000E-03, 2.00000E-03
598 , 3.00000E-03, 4.00000E-03, 5.00000E-03
599 , 6.00000E-03, 8.00000E-03, 1.00000E-02
600 , 1.50000E-02, 2.00000E-02, 3.00000E-02
601 , 4.00000E-02, 5.00000E-02, 6.00000E-02
602 , 8.00000E-02, 1.00000E-01, 1.50000E-01
603 , 2.00000E-01, 3.00000E-01, 4.00000E-01
604 , 5.00000E-01, 6.00000E-01, 8.00000E-01
605 , 1.00000E+00, 1.25000E+00, 1.50000E+00
606 , 2.00000E+00, 3.00000E+00, 4.00000E+00
607 , 5.00000E+00, 6.00000E+00, 8.00000E+00
608 , 1.00000E+01, 1.50000E+01, 2.00000E+01 };
609
610 return Interpolate(energyMeV,en,mu,kN);
611
612}
613
614//_____________________________________________________________________________
615Double_t AliTRDsim::GetMuO2(Double_t energyMeV)
616{
617 //
618 // Returns the photon absorbtion cross section for oxygen
619 //
620
621 const Int_t kN = 36;
622
623 Double_t mu[kN] = { 4.590E+03, 1.549E+03, 6.949E+02
624 , 2.171E+02, 9.315E+01, 4.790E+01
625 , 2.770E+01, 1.163E+01, 5.952E+00
626 , 1.836E+00, 8.651E-01, 3.779E-01
627 , 2.585E-01, 2.132E-01, 1.907E-01
628 , 1.678E-01, 1.551E-01, 1.361E-01
629 , 1.237E-01, 1.070E-01, 9.566E-02
630 , 8.729E-02, 8.070E-02, 7.087E-02
631 , 6.372E-02, 5.697E-02, 5.185E-02
632 , 4.459E-02, 3.597E-02, 3.100E-02
633 , 2.777E-02, 2.552E-02, 2.263E-02
634 , 2.089E-02, 1.866E-02, 1.770E-02 };
635
636 Double_t en[kN] = { 1.00000E-03, 1.50000E-03, 2.00000E-03
637 , 3.00000E-03, 4.00000E-03, 5.00000E-03
638 , 6.00000E-03, 8.00000E-03, 1.00000E-02
639 , 1.50000E-02, 2.00000E-02, 3.00000E-02
640 , 4.00000E-02, 5.00000E-02, 6.00000E-02
641 , 8.00000E-02, 1.00000E-01, 1.50000E-01
642 , 2.00000E-01, 3.00000E-01, 4.00000E-01
643 , 5.00000E-01, 6.00000E-01, 8.00000E-01
644 , 1.00000E+00, 1.25000E+00, 1.50000E+00
645 , 2.00000E+00, 3.00000E+00, 4.00000E+00
646 , 5.00000E+00, 6.00000E+00, 8.00000E+00
647 , 1.00000E+01, 1.50000E+01, 2.00000E+01 };
648
649 return Interpolate(energyMeV,en,mu,kN);
650
651}
652
653//_____________________________________________________________________________
654Double_t AliTRDsim::GetMuHe(Double_t energyMeV)
655{
656 //
657 // Returns the photon absorbtion cross section for helium
658 //
659
660 const Int_t kN = 36;
661
662 Double_t mu[kN] = { 6.084E+01, 1.676E+01, 6.863E+00
663 , 2.007E+00, 9.329E-01, 5.766E-01
664 , 4.195E-01, 2.933E-01, 2.476E-01
665 , 2.092E-01, 1.960E-01, 1.838E-01
666 , 1.763E-01, 1.703E-01, 1.651E-01
667 , 1.562E-01, 1.486E-01, 1.336E-01
668 , 1.224E-01, 1.064E-01, 9.535E-02
669 , 8.707E-02, 8.054E-02, 7.076E-02
670 , 6.362E-02, 5.688E-02, 5.173E-02
671 , 4.422E-02, 3.503E-02, 2.949E-02
672 , 2.577E-02, 2.307E-02, 1.940E-02
673 , 1.703E-02, 1.363E-02, 1.183E-02 };
674
675 Double_t en[kN] = { 1.00000E-03, 1.50000E-03, 2.00000E-03
676 , 3.00000E-03, 4.00000E-03, 5.00000E-03
677 , 6.00000E-03, 8.00000E-03, 1.00000E-02
678 , 1.50000E-02, 2.00000E-02, 3.00000E-02
679 , 4.00000E-02, 5.00000E-02, 6.00000E-02
680 , 8.00000E-02, 1.00000E-01, 1.50000E-01
681 , 2.00000E-01, 3.00000E-01, 4.00000E-01
682 , 5.00000E-01, 6.00000E-01, 8.00000E-01
683 , 1.00000E+00, 1.25000E+00, 1.50000E+00
684 , 2.00000E+00, 3.00000E+00, 4.00000E+00
685 , 5.00000E+00, 6.00000E+00, 8.00000E+00
686 , 1.00000E+01, 1.50000E+01, 2.00000E+01 };
687
688 return Interpolate(energyMeV,en,mu,kN);
689
690}
691
692//_____________________________________________________________________________
693Double_t AliTRDsim::Interpolate(Double_t energyMeV
694 , Double_t *en, Double_t *mu, Int_t n)
695{
696 //
697 // Interpolates the photon absorbtion cross section
698 // for a given energy <energyMeV>.
699 //
700
701 Double_t de = 0;
702 Int_t index = 0;
703 Int_t istat = Locate(en,n,energyMeV,index,de);
704 if (istat == 0) {
705 return (mu[index] - de * (mu[index] - mu[index+1])
706 / (en[index+1] - en[index] ));
707 }
708 else {
709 return 0.0;
710 }
711
712}
713
714//_____________________________________________________________________________
715Int_t AliTRDsim::Locate(Double_t *xv, Int_t n, Double_t xval
716 , Int_t &kl, Double_t &dx)
717{
718 //
719 // Locates a point (xval) in a 1-dim grid (xv(n))
720 //
721
722 if (xval >= xv[n-1]) return 1;
723 if (xval < xv[0]) return -1;
724
725 Int_t km;
726 Int_t kh = n - 1;
727
728 kl = 0;
729 while (kh - kl > 1) {
730 if (xval < xv[km = (kl+kh)/2]) kh = km;
731 else kl = km;
732 }
733 if (xval < xv[kl] || xval > xv[kl+1] || kl >= n-1) {
734 printf("Locate failed xv[%d] %f xval %f xv[%d] %f!!!\n"
735 ,kl,xv[kl],xval,kl+1,xv[kl+1]);
736 exit(1);
737 }
738
739 dx = xval - xv[kl];
740
741 return 0;
742
743}