]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TRD/AliTRDv1.cxx
update for non null magnetic field
[u/mrichter/AliRoot.git] / TRD / AliTRDv1.cxx
CommitLineData
4c039060 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
88cb7938 16/* $Id$ */
4c039060 17
030b4415 18////////////////////////////////////////////////////////////////////////////
19// //
20// Transition Radiation Detector version 1 -- slow simulator //
21// //
22////////////////////////////////////////////////////////////////////////////
fe4da5cc 23
769257f4 24#include <stdlib.h>
25
793ff80c 26#include <TF1.h>
1819f4bb 27#include <TLorentzVector.h>
88cb7938 28#include <TMath.h>
29#include <TRandom.h>
30#include <TVector.h>
31#include <TVirtualMC.h>
f57bb418 32#include <TGeoManager.h>
268f57b1 33#include <TGeoMatrix.h>
170c35f1 34#include <TGeoPhysicalNode.h>
fe4da5cc 35
d3f347ff 36#include "AliConst.h"
45160b1f 37#include "AliLog.h"
e6add757 38#include "AliTrackReference.h"
45160b1f 39#include "AliMC.h"
88cb7938 40#include "AliRun.h"
02cb65d6 41#include "AliGeomManager.h"
030b4415 42
88cb7938 43#include "AliTRDgeometry.h"
a076fc2f 44#include "AliTRDCommonParam.h"
793ff80c 45#include "AliTRDhit.h"
cb2f9e9b 46#include "AliTRDsimTR.h"
88cb7938 47#include "AliTRDv1.h"
851d3db9 48
fe4da5cc 49ClassImp(AliTRDv1)
8230f242 50
51//_____________________________________________________________________________
030b4415 52AliTRDv1::AliTRDv1()
53 :AliTRD()
54 ,fTRon(kFALSE)
55 ,fTR(NULL)
56 ,fTypeOfStepManager(0)
57 ,fStepSize(0)
f2979d08 58 ,fWion(0)
030b4415 59 ,fDeltaE(NULL)
60 ,fDeltaG(NULL)
61 ,fTrackLength0(0)
62 ,fPrimaryTrackPid(0)
8230f242 63{
64 //
65 // Default constructor
66 //
67
8230f242 68}
69
fe4da5cc 70//_____________________________________________________________________________
71AliTRDv1::AliTRDv1(const char *name, const char *title)
030b4415 72 :AliTRD(name,title)
73 ,fTRon(kTRUE)
74 ,fTR(NULL)
67c47633 75 ,fTypeOfStepManager(2)
030b4415 76 ,fStepSize(0.1)
f2979d08 77 ,fWion(0)
030b4415 78 ,fDeltaE(NULL)
79 ,fDeltaG(NULL)
80 ,fTrackLength0(0)
81 ,fPrimaryTrackPid(0)
fe4da5cc 82{
83 //
851d3db9 84 // Standard constructor for Transition Radiation Detector version 1
fe4da5cc 85 //
82bbf98a 86
5c7f4665 87 SetBufferSize(128000);
88
a076fc2f 89 if (AliTRDCommonParam::Instance()->IsXenon()) {
f2979d08 90 fWion = 23.53; // Ionization energy XeCO2 (85/15)
91 }
a076fc2f 92 else if (AliTRDCommonParam::Instance()->IsArgon()) {
f2979d08 93 fWion = 27.21; // Ionization energy ArCO2 (82/18)
94 }
95 else {
96 AliFatal("Wrong gas mixture");
97 exit(1);
98 }
99
5c7f4665 100}
101
102//_____________________________________________________________________________
103AliTRDv1::~AliTRDv1()
104{
dd9a6ee3 105 //
106 // AliTRDv1 destructor
107 //
82bbf98a 108
030b4415 109 if (fDeltaE) {
110 delete fDeltaE;
111 fDeltaE = 0;
112 }
113
114 if (fDeltaG) {
115 delete fDeltaG;
116 fDeltaG = 0;
117 }
118
119 if (fTR) {
120 delete fTR;
121 fTR = 0;
122 }
82bbf98a 123
fe4da5cc 124}
125
f57bb418 126//_____________________________________________________________________________
127void AliTRDv1::AddAlignableVolumes() const
128{
129 //
130 // Create entries for alignable volumes associating the symbolic volume
131 // name with the corresponding volume path. Needs to be syncronized with
132 // eventual changes in the geometry.
133 //
134
135 TString volPath;
136 TString symName;
137
92cd97ad 138 TString vpStr = "ALIC_1/B077_1/BSEGMO";
139 TString vpApp1 = "_1/BTRD";
140 TString vpApp2 = "_1";
141 TString vpApp3a = "/UTR1_1/UTS1_1/UTI1_1/UT";
142 TString vpApp3b = "/UTR2_1/UTS2_1/UTI2_1/UT";
143 TString vpApp3c = "/UTR3_1/UTS3_1/UTI3_1/UT";
f57bb418 144
92cd97ad 145 TString snStr = "TRD/sm";
146 TString snApp1 = "/st";
147 TString snApp2 = "/pl";
f57bb418 148
149 //
150 // The super modules
151 // The symbolic names are: TRD/sm00
152 // ...
153 // TRD/sm17
154 //
053767a4 155 for (Int_t isector = 0; isector < AliTRDgeometry::Nsector(); isector++) {
f57bb418 156
157 volPath = vpStr;
053767a4 158 volPath += isector;
f57bb418 159 volPath += vpApp1;
053767a4 160 volPath += isector;
f57bb418 161 volPath += vpApp2;
162
163 symName = snStr;
053767a4 164 symName += Form("%02d",isector);
f57bb418 165
166 gGeoManager->SetAlignableEntry(symName.Data(),volPath.Data());
167
168 }
169
170 //
171 // The readout chambers
172 // The symbolic names are: TRD/sm00/st0/pl0
173 // ...
174 // TRD/sm17/st4/pl5
175 //
02cb65d6 176 AliGeomManager::ELayerID idTRD1 = AliGeomManager::kTRD1;
177 Int_t layer, modUID;
178
053767a4 179 for (Int_t isector = 0; isector < AliTRDgeometry::Nsector(); isector++) {
8bf0cd64 180
053767a4 181 if (fGeometry->GetSMstatus(isector) == 0) continue;
8bf0cd64 182
053767a4 183 for (Int_t istack = 0; istack < AliTRDgeometry::Nstack(); istack++) {
184 for (Int_t ilayer = 0; ilayer < AliTRDgeometry::Nlayer(); ilayer++) {
f57bb418 185
053767a4 186 layer = idTRD1 + ilayer;
187 modUID = AliGeomManager::LayerToVolUIDSafe(layer,isector*5+istack);
188
189 Int_t idet = AliTRDgeometry::GetDetectorSec(ilayer,istack);
f57bb418 190
191 volPath = vpStr;
053767a4 192 volPath += isector;
f57bb418 193 volPath += vpApp1;
053767a4 194 volPath += isector;
f57bb418 195 volPath += vpApp2;
053767a4 196 switch (isector) {
92cd97ad 197 case 13:
198 case 14:
199 case 15:
053767a4 200 if (istack == 2) {
92cd97ad 201 continue;
202 }
203 volPath += vpApp3c;
204 break;
205 case 11:
206 case 12:
207 volPath += vpApp3b;
208 break;
209 default:
210 volPath += vpApp3a;
211 };
f57bb418 212 volPath += Form("%02d",idet);
213 volPath += vpApp2;
214
215 symName = snStr;
053767a4 216 symName += Form("%02d",isector);
f57bb418 217 symName += snApp1;
053767a4 218 symName += istack;
f57bb418 219 symName += snApp2;
053767a4 220 symName += ilayer;
f57bb418 221
51a5f1d0 222 TGeoPNEntry *alignableEntry =
02cb65d6 223 gGeoManager->SetAlignableEntry(symName.Data(),volPath.Data(),modUID);
f57bb418 224
170c35f1 225 // Add the tracking to local matrix following the TPC example
51a5f1d0 226 if (alignableEntry) {
02cb65d6 227 TGeoHMatrix *globMatrix = alignableEntry->GetGlobalOrig();
053767a4 228 Double_t sectorAngle = 20.0 * (isector % 18) + 10.0;
51a5f1d0 229 TGeoHMatrix *t2lMatrix = new TGeoHMatrix();
230 t2lMatrix->RotateZ(sectorAngle);
231 t2lMatrix->MultiplyLeft(&(globMatrix->Inverse()));
232 alignableEntry->SetMatrix(t2lMatrix);
233 }
234 else {
235 AliError(Form("Alignable entry %s is not valid!",symName.Data()));
236 }
8bf0cd64 237
f57bb418 238 }
239 }
240 }
241
242}
243
fe4da5cc 244//_____________________________________________________________________________
245void AliTRDv1::CreateGeometry()
246{
247 //
851d3db9 248 // Create the GEANT geometry for the Transition Radiation Detector - Version 1
5c7f4665 249 // This version covers the full azimuth.
d3f347ff 250 //
251
82bbf98a 252 // Check that FRAME is there otherwise we have no place where to put the TRD
8230f242 253 AliModule* frame = gAlice->GetModule("FRAME");
030b4415 254 if (!frame) {
255 AliError("TRD needs FRAME to be present\n");
256 return;
257 }
d3f347ff 258
82bbf98a 259 // Define the chambers
260 AliTRD::CreateGeometry();
d3f347ff 261
fe4da5cc 262}
263
264//_____________________________________________________________________________
265void AliTRDv1::CreateMaterials()
266{
267 //
851d3db9 268 // Create materials for the Transition Radiation Detector version 1
fe4da5cc 269 //
82bbf98a 270
d3f347ff 271 AliTRD::CreateMaterials();
82bbf98a 272
fe4da5cc 273}
274
793ff80c 275//_____________________________________________________________________________
276void AliTRDv1::CreateTRhit(Int_t det)
277{
278 //
279 // Creates an electron cluster from a TR photon.
280 // The photon is assumed to be created a the end of the radiator. The
281 // distance after which it deposits its energy takes into account the
282 // absorbtion of the entrance window and of the gas mixture in drift
283 // volume.
284 //
285
793ff80c 286 // Maximum number of TR photons per track
287 const Int_t kNTR = 50;
288
030b4415 289 TLorentzVector mom;
290 TLorentzVector pos;
793ff80c 291
ce0d6231 292 Float_t eTR[kNTR];
293 Int_t nTR;
793ff80c 294
ce0d6231 295 // Create TR photons
296 gMC->TrackMomentum(mom);
297 Float_t pTot = mom.Rho();
298 fTR->CreatePhotons(11,pTot,nTR,eTR);
299 if (nTR > kNTR) {
300 AliFatal(Form("Boundary error: nTR = %d, kNTR = %d",nTR,kNTR));
301 }
f73816f5 302
ce0d6231 303 // Loop through the TR photons
304 for (Int_t iTR = 0; iTR < nTR; iTR++) {
793ff80c 305
ce0d6231 306 Float_t energyMeV = eTR[iTR] * 0.001;
307 Float_t energyeV = eTR[iTR] * 1000.0;
308 Float_t absLength = 0.0;
309 Float_t sigma = 0.0;
793ff80c 310
ce0d6231 311 // Take the absorbtion in the entrance window into account
312 Double_t muMy = fTR->GetMuMy(energyMeV);
313 sigma = muMy * fFoilDensity;
314 if (sigma > 0.0) {
315 absLength = gRandom->Exp(1.0/sigma);
316 if (absLength < AliTRDgeometry::MyThick()) {
842287f2 317 continue;
318 }
ce0d6231 319 }
320 else {
321 continue;
322 }
793ff80c 323
ce0d6231 324 // The absorbtion cross sections in the drift gas
325 // Gas-mixture (Xe/CO2)
f2979d08 326 Double_t muNo = 0.0;
a076fc2f 327 if (AliTRDCommonParam::Instance()->IsXenon()) {
f2979d08 328 muNo = fTR->GetMuXe(energyMeV);
329 }
a076fc2f 330 else if (AliTRDCommonParam::Instance()->IsArgon()) {
f2979d08 331 muNo = fTR->GetMuAr(energyMeV);
332 }
ce0d6231 333 Double_t muCO = fTR->GetMuCO(energyMeV);
f2979d08 334 sigma = (fGasNobleFraction * muNo + (1.0 - fGasNobleFraction) * muCO)
335 * fGasDensity
336 * fTR->GetTemp();
ce0d6231 337
338 // The distance after which the energy of the TR photon
339 // is deposited.
340 if (sigma > 0.0) {
341 absLength = gRandom->Exp(1.0/sigma);
342 if (absLength > (AliTRDgeometry::DrThick()
343 + AliTRDgeometry::AmThick())) {
842287f2 344 continue;
345 }
ce0d6231 346 }
347 else {
348 continue;
349 }
793ff80c 350
ce0d6231 351 // The position of the absorbtion
352 Float_t posHit[3];
353 gMC->TrackPosition(pos);
354 posHit[0] = pos[0] + mom[0] / pTot * absLength;
355 posHit[1] = pos[1] + mom[1] / pTot * absLength;
356 posHit[2] = pos[2] + mom[2] / pTot * absLength;
793ff80c 357
ce0d6231 358 // Create the charge
f2979d08 359 Int_t q = ((Int_t) (energyeV / fWion));
793ff80c 360
ce0d6231 361 // Add the hit to the array. TR photon hits are marked
362 // by negative charge
363 AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber()
364 ,det
365 ,posHit
366 ,-q
25ca55ce 367 ,gMC->TrackTime()*1.0e06
d4c6453d 368 ,kTRUE);
793ff80c 369
370 }
371
372}
373
5c7f4665 374//_____________________________________________________________________________
375void AliTRDv1::Init()
376{
377 //
378 // Initialise Transition Radiation Detector after geometry has been built.
5c7f4665 379 //
380
381 AliTRD::Init();
382
45160b1f 383 AliDebug(1,"Slow simulator\n");
bd0f8685 384
385 // Switch on TR simulation as default
386 if (!fTRon) {
45160b1f 387 AliInfo("TR simulation off");
bd0f8685 388 }
389 else {
cb2f9e9b 390 fTR = new AliTRDsimTR();
bd0f8685 391 }
5c7f4665 392
f2979d08 393 // First ionization potential (eV) for the gas mixture
5c7f4665 394 const Float_t kPoti = 12.1;
395 // Maximum energy (50 keV);
396 const Float_t kEend = 50000.0;
397 // Ermilova distribution for the delta-ray spectrum
030b4415 398 Float_t poti = TMath::Log(kPoti);
399 Float_t eEnd = TMath::Log(kEend);
a328fff9 400
401 // Ermilova distribution for the delta-ray spectrum
c4214bc0 402 fDeltaE = new TF1("deltae" ,Ermilova ,poti,eEnd,0);
a328fff9 403
404 // Geant3 distribution for the delta-ray spectrum
c4214bc0 405 fDeltaG = new TF1("deltag",IntSpecGeant,2.421257,28.536469,0);
5c7f4665 406
45160b1f 407 AliDebug(1,"+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++");
5c7f4665 408
fe4da5cc 409}
410
5c7f4665 411//_____________________________________________________________________________
412void AliTRDv1::StepManager()
a328fff9 413{
414 //
c4214bc0 415 // Slow simulator. Every charged track produces electron cluster as hits
a328fff9 416 // along its path across the drift volume.
417 //
418
419 switch (fTypeOfStepManager) {
a6dd11e9 420 case 0:
421 StepManagerErmilova();
422 break;
423 case 1:
424 StepManagerGeant();
425 break;
426 case 2:
427 StepManagerFixedStep();
428 break;
429 default:
430 AliWarning("Not a valid Step Manager.");
a328fff9 431 }
432
433}
434
435//_____________________________________________________________________________
436void AliTRDv1::SelectStepManager(Int_t t)
437{
438 //
439 // Selects a step manager type:
440 // 0 - Ermilova
441 // 1 - Geant3
442 // 2 - Fixed step size
443 //
444
a328fff9 445 fTypeOfStepManager = t;
45160b1f 446 AliInfo(Form("Step Manager type %d was selected",fTypeOfStepManager));
a328fff9 447
448}
449
450//_____________________________________________________________________________
451void AliTRDv1::StepManagerGeant()
452{
453 //
c4214bc0 454 // Slow simulator. Every charged track produces electron cluster as hits
a328fff9 455 // along its path across the drift volume. The step size is set acording
456 // to Bethe-Bloch. The energy distribution of the delta electrons follows
457 // a spectrum taken from Geant3.
458 //
f2979d08 459 // Works only for Xe/CO2!!
460 //
f2e3a0b5 461 // Version by A. Bercuci
462 //
463
053767a4 464 Int_t layer = 0;
465 Int_t stack = 0;
466 Int_t sector = 0;
467 Int_t det = 0;
c4214bc0 468 Int_t iPdg;
469 Int_t qTot;
470
471 Float_t hits[3];
472 Float_t charge;
473 Float_t aMass;
474
030b4415 475 Double_t pTot = 0;
c4214bc0 476 Double_t eDelta;
030b4415 477 Double_t betaGamma;
478 Double_t pp;
f2e3a0b5 479 Double_t stepSize = 0;
c4214bc0 480
481 Bool_t drRegion = kFALSE;
482 Bool_t amRegion = kFALSE;
483
2c8bf4aa 484 TString cIdPath;
485 Char_t cIdSector[3];
486 cIdSector[2] = 0;
487
c4214bc0 488 TString cIdCurrent;
489 TString cIdSensDr = "J";
490 TString cIdSensAm = "K";
491 Char_t cIdChamber[3];
492 cIdChamber[2] = 0;
493
030b4415 494 TLorentzVector pos;
495 TLorentzVector mom;
c4214bc0 496
030b4415 497 TArrayI processes;
f2e3a0b5 498
053767a4 499 const Int_t kNlayer = AliTRDgeometry::Nlayer();
500 const Int_t kNstack = AliTRDgeometry::Nstack();
501 const Int_t kNdetsec = kNlayer * kNstack;
c4214bc0 502
030b4415 503 const Double_t kBig = 1.0e+12; // Infinitely big
c4214bc0 504 const Float_t kPTotMaxEl = 0.002; // Maximum momentum for e+ e- g
505
506 // Minimum energy for the step size adjustment
507 const Float_t kEkinMinStep = 1.0e-5;
508 // energy threshold for production of delta electrons
f2e3a0b5 509 const Float_t kECut = 1.0e4;
510 // Parameters entering the parametrized range for delta electrons
030b4415 511 const Float_t kRa = 5.37e-4;
f2e3a0b5 512 const Float_t kRb = 0.9815;
030b4415 513 const Float_t kRc = 3.123e-3;
f2e3a0b5 514 // Gas density -> To be made user adjustable !
030b4415 515 // [0.85*0.00549+0.15*0.00186 (Xe-CO2 85-15)]
516 const Float_t kRho = 0.004945 ;
a328fff9 517
c4214bc0 518 // Plateau value of the energy-loss for electron in xenon
030b4415 519 // The averaged value (26/3/99)
c4214bc0 520 const Float_t kPlateau = 1.55;
030b4415 521 // dN1/dx|min for the gas mixture (90% Xe + 10% CO2)
522 const Float_t kPrim = 19.34;
c4214bc0 523 // First ionization potential (eV) for the gas mixture (90% Xe + 10% CO2)
524 const Float_t kPoti = 12.1;
030b4415 525 // PDG code electron
526 const Int_t kPdgElectron = 11;
c4214bc0 527
528 // Set the maximum step size to a very large number for all
529 // neutral particles and those outside the driftvolume
530 gMC->SetMaxStep(kBig);
531
532 // Use only charged tracks
533 if (( gMC->TrackCharge() ) &&
c4214bc0 534 (!gMC->IsTrackDisappeared())) {
535
536 // Inside a sensitive volume?
537 drRegion = kFALSE;
538 amRegion = kFALSE;
539 cIdCurrent = gMC->CurrentVolName();
540 if (cIdSensDr == cIdCurrent[1]) {
541 drRegion = kTRUE;
542 }
543 if (cIdSensAm == cIdCurrent[1]) {
544 amRegion = kTRUE;
545 }
546 if (drRegion || amRegion) {
a328fff9 547
c4214bc0 548 // The hit coordinates and charge
549 gMC->TrackPosition(pos);
550 hits[0] = pos[0];
551 hits[1] = pos[1];
552 hits[2] = pos[2];
553
2c8bf4aa 554 // The sector number (0 - 17), according to standard coordinate system
555 cIdPath = gGeoManager->GetPath();
556 cIdSector[0] = cIdPath[21];
557 cIdSector[1] = cIdPath[22];
053767a4 558 sector = atoi(cIdSector);
c4214bc0 559
053767a4 560 // The layer and stack number
030b4415 561 cIdChamber[0] = cIdCurrent[2];
562 cIdChamber[1] = cIdCurrent[3];
c4214bc0 563 Int_t idChamber = (atoi(cIdChamber) % kNdetsec);
053767a4 564 stack = ((Int_t) idChamber / kNlayer);
565 layer = ((Int_t) idChamber % kNlayer);
c4214bc0 566
ce0d6231 567 // The detector number
053767a4 568 det = fGeometry->GetDetector(layer,stack,sector);
c4214bc0 569
ce0d6231 570 // Special hits only in the drift region
571 if ((drRegion) &&
572 (gMC->IsTrackEntering())) {
c4214bc0 573
ce0d6231 574 // Create a track reference at the entrance of each
575 // chamber that contains the momentum components of the particle
576 gMC->TrackMomentum(mom);
e6add757 577 AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber(), AliTrackReference::kTRD);
c4214bc0 578
ce0d6231 579 // Create the hits from TR photons if electron/positron is
580 // entering the drift volume
581 if ((fTR) &&
582 (TMath::Abs(gMC->TrackPid()) == kPdgElectron)) {
f2979d08 583 CreateTRhit(det);
ce0d6231 584 }
f2e3a0b5 585
ce0d6231 586 }
587 else if ((amRegion) &&
588 (gMC->IsTrackExiting())) {
f2e3a0b5 589
ce0d6231 590 // Create a track reference at the exit of each
591 // chamber that contains the momentum components of the particle
592 gMC->TrackMomentum(mom);
e6add757 593 AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber(), AliTrackReference::kTRD);
c4214bc0 594
ce0d6231 595 }
c4214bc0 596
ce0d6231 597 // Calculate the energy of the delta-electrons
598 // modified by Alex Bercuci (A.Bercuci@gsi.de) on 26.01.06
599 // take into account correlation with the underlying GEANT tracking
600 // mechanism. see
601 // http://www-linux.gsi.de/~abercuci/Contributions/TRD/index.html
602 //
603 // determine the most significant process (last on the processes list)
604 // which caused this hit
605 gMC->StepProcesses(processes);
606 Int_t nofprocesses = processes.GetSize();
607 Int_t pid;
608 if (!nofprocesses) {
609 pid = 0;
610 }
611 else {
053767a4 612 pid = processes[nofprocesses-1];
ce0d6231 613 }
f2e3a0b5 614
ce0d6231 615 // Generate Edep according to GEANT parametrisation
616 eDelta = TMath::Exp(fDeltaG->GetRandom()) - kPoti;
617 eDelta = TMath::Max(eDelta,0.0);
618 Float_t prRange = 0.0;
619 Float_t range = gMC->TrackLength() - fTrackLength0;
620 // merge GEANT tracker information with locally cooked one
621 if (gAlice->GetMCApp()->GetCurrentTrackNumber() == fPrimaryTrackPid) {
622 if (pid == 27) {
623 if (eDelta >= kECut) {
624 prRange = kRa * eDelta * 0.001
625 * (1.0 - kRb / (1.0 + kRc * eDelta * 0.001)) / kRho;
626 if (prRange >= (3.7 - range)) {
627 eDelta *= 0.1;
628 }
629 }
630 }
631 else if (pid == 1) {
632 if (eDelta < kECut) {
633 eDelta *= 0.5;
634 }
635 else {
636 prRange = kRa * eDelta * 0.001
637 * (1.0 - kRb / (1.0 + kRc * eDelta * 0.001)) / kRho;
638 if (prRange >= ((AliTRDgeometry::DrThick()
639 + AliTRDgeometry::AmThick()) - range)) {
640 eDelta *= 0.05;
641 }
642 else {
f2e3a0b5 643 eDelta *= 0.5;
ce0d6231 644 }
645 }
646 }
f2e3a0b5 647 else {
648 eDelta = 0.0;
ce0d6231 649 }
650 }
651 else {
652 eDelta = 0.0;
653 }
c4214bc0 654
ce0d6231 655 // Generate the electron cluster size
656 if (eDelta > 0.0) {
f2e3a0b5 657
f2979d08 658 qTot = ((Int_t) (eDelta / fWion) + 1);
ce0d6231 659
660 // Create a new dEdx hit
030b4415 661 AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber()
662 ,det
663 ,hits
664 ,qTot
25ca55ce 665 ,gMC->TrackTime()*1.0e06
030b4415 666 ,drRegion);
c4214bc0 667
ce0d6231 668 }
669
670 // Calculate the maximum step size for the next tracking step
671 // Produce only one hit if Ekin is below cutoff
672 aMass = gMC->TrackMass();
673 if ((gMC->Etot() - aMass) > kEkinMinStep) {
674
675 // The energy loss according to Bethe Bloch
676 iPdg = TMath::Abs(gMC->TrackPid());
677 if ((iPdg != kPdgElectron) ||
678 ((iPdg == kPdgElectron) &&
679 (pTot < kPTotMaxEl))) {
680 gMC->TrackMomentum(mom);
681 pTot = mom.Rho();
682 betaGamma = pTot / aMass;
683 pp = BetheBlochGeant(betaGamma);
684 // Take charge > 1 into account
685 charge = gMC->TrackCharge();
686 if (TMath::Abs(charge) > 1) {
687 pp = pp * charge*charge;
688 }
689 }
690 else {
691 // Electrons above 20 Mev/c are at the plateau
692 pp = kPrim * kPlateau;
693 }
f2e3a0b5 694
ce0d6231 695 Int_t nsteps = 0;
696 do {
697 nsteps = gRandom->Poisson(pp);
698 } while(!nsteps);
699 stepSize = 1.0 / nsteps;
700 gMC->SetMaxStep(stepSize);
f2e3a0b5 701
c4214bc0 702 }
f2e3a0b5 703
c4214bc0 704 }
f2e3a0b5 705
c4214bc0 706 }
f2e3a0b5 707
a328fff9 708}
709
710//_____________________________________________________________________________
711void AliTRDv1::StepManagerErmilova()
5c7f4665 712{
713 //
5c7f4665 714 // Slow simulator. Every charged track produces electron cluster as hits
715 // along its path across the drift volume. The step size is set acording
716 // to Bethe-Bloch. The energy distribution of the delta electrons follows
717 // a spectrum taken from Ermilova et al.
718 //
f2979d08 719 // Works only for Xe/CO2!!
720 //
5c7f4665 721
053767a4 722 Int_t layer = 0;
723 Int_t stack = 0;
724 Int_t sector = 0;
725 Int_t det = 0;
851d3db9 726 Int_t iPdg;
793ff80c 727 Int_t qTot;
5c7f4665 728
793ff80c 729 Float_t hits[3];
a5cadd36 730 Double_t random[1];
5c7f4665 731 Float_t charge;
732 Float_t aMass;
733
030b4415 734 Double_t pTot = 0.0;
5c7f4665 735 Double_t eDelta;
030b4415 736 Double_t betaGamma;
737 Double_t pp;
f73816f5 738 Double_t stepSize;
5c7f4665 739
332e9569 740 Bool_t drRegion = kFALSE;
741 Bool_t amRegion = kFALSE;
742
2c8bf4aa 743 TString cIdPath;
744 Char_t cIdSector[3];
745 cIdSector[2] = 0;
746
332e9569 747 TString cIdCurrent;
748 TString cIdSensDr = "J";
749 TString cIdSensAm = "K";
593a9fc3 750 Char_t cIdChamber[3];
751 cIdChamber[2] = 0;
332e9569 752
030b4415 753 TLorentzVector pos;
754 TLorentzVector mom;
82bbf98a 755
053767a4 756 const Int_t kNlayer = AliTRDgeometry::Nlayer();
757 const Int_t kNstack = AliTRDgeometry::Nstack();
758 const Int_t kNdetsec = kNlayer * kNstack;
e644678a 759
030b4415 760 const Double_t kBig = 1.0e+12; // Infinitely big
a328fff9 761 const Float_t kPTotMaxEl = 0.002; // Maximum momentum for e+ e- g
5c7f4665 762
f73816f5 763 // Minimum energy for the step size adjustment
764 const Float_t kEkinMinStep = 1.0e-5;
a328fff9 765
5c7f4665 766 // Plateau value of the energy-loss for electron in xenon
030b4415 767 // The averaged value (26/3/99)
a3c76cdc 768 const Float_t kPlateau = 1.55;
030b4415 769 // dN1/dx|min for the gas mixture (90% Xe + 10% CO2)
770 const Float_t kPrim = 48.0;
5c7f4665 771 // First ionization potential (eV) for the gas mixture (90% Xe + 10% CO2)
a3c76cdc 772 const Float_t kPoti = 12.1;
030b4415 773 // PDG code electron
774 const Int_t kPdgElectron = 11;
5c7f4665 775
776 // Set the maximum step size to a very large number for all
777 // neutral particles and those outside the driftvolume
778 gMC->SetMaxStep(kBig);
779
780 // Use only charged tracks
781 if (( gMC->TrackCharge() ) &&
5c7f4665 782 (!gMC->IsTrackDisappeared())) {
fe4da5cc 783
5c7f4665 784 // Inside a sensitive volume?
332e9569 785 drRegion = kFALSE;
786 amRegion = kFALSE;
787 cIdCurrent = gMC->CurrentVolName();
e6674585 788 if (cIdSensDr == cIdCurrent[1]) {
332e9569 789 drRegion = kTRUE;
790 }
e6674585 791 if (cIdSensAm == cIdCurrent[1]) {
332e9569 792 amRegion = kTRUE;
793 }
794 if (drRegion || amRegion) {
fe4da5cc 795
5c7f4665 796 // The hit coordinates and charge
797 gMC->TrackPosition(pos);
798 hits[0] = pos[0];
799 hits[1] = pos[1];
800 hits[2] = pos[2];
5c7f4665 801
2c8bf4aa 802 // The sector number (0 - 17), according to standard coordinate system
803 cIdPath = gGeoManager->GetPath();
804 cIdSector[0] = cIdPath[21];
805 cIdSector[1] = cIdPath[22];
053767a4 806 sector = atoi(cIdSector);
82bbf98a 807
332e9569 808 // The plane and chamber number
809 cIdChamber[0] = cIdCurrent[2];
810 cIdChamber[1] = cIdCurrent[3];
e644678a 811 Int_t idChamber = (atoi(cIdChamber) % kNdetsec);
053767a4 812 stack = ((Int_t) idChamber / kNlayer);
813 layer = ((Int_t) idChamber % kNlayer);
82bbf98a 814
ce0d6231 815 // The detector number
053767a4 816 det = fGeometry->GetDetector(layer,stack,sector);
5c7f4665 817
ce0d6231 818 // Special hits only in the drift region
819 if ((drRegion) &&
820 (gMC->IsTrackEntering())) {
5c7f4665 821
ce0d6231 822 // Create a track reference at the entrance of each
823 // chamber that contains the momentum components of the particle
824 gMC->TrackMomentum(mom);
e6add757 825 AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber(), AliTrackReference::kTRD);
793ff80c 826
ce0d6231 827 // Create the hits from TR photons if electron/positron is
828 // entering the drift volume
829 if ((fTR) &&
830 (TMath::Abs(gMC->TrackPid()) == kPdgElectron)) {
f2979d08 831 CreateTRhit(det);
ce0d6231 832 }
f73816f5 833
ce0d6231 834 }
835 else if ((amRegion) &&
836 (gMC->IsTrackExiting())) {
f73816f5 837
ce0d6231 838 // Create a track reference at the exit of each
839 // chamber that contains the momentum components of the particle
840 gMC->TrackMomentum(mom);
e6add757 841 AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber(), AliTrackReference::kTRD);
f73816f5 842
ce0d6231 843 }
844
845 // Calculate the energy of the delta-electrons
846 eDelta = TMath::Exp(fDeltaE->GetRandom()) - kPoti;
847 eDelta = TMath::Max(eDelta,0.0);
848
849 // Generate the electron cluster size
850 if (eDelta > 0.0) {
851
f2979d08 852 qTot = ((Int_t) (eDelta / fWion) + 1);
f73816f5 853
030b4415 854 // Create a new dEdx hit
332e9569 855 if (drRegion) {
a328fff9 856 AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber()
030b4415 857 ,det
858 ,hits
859 ,qTot
25ca55ce 860 ,gMC->TrackTime()*1.0e06
030b4415 861 ,kTRUE);
862 }
5c7f4665 863 else {
a328fff9 864 AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber()
030b4415 865 ,det
866 ,hits
867 ,qTot
25ca55ce 868 ,gMC->TrackTime()*1.0e06
030b4415 869 ,kFALSE);
870 }
f73816f5 871
ce0d6231 872 }
873
874 // Calculate the maximum step size for the next tracking step
875 // Produce only one hit if Ekin is below cutoff
876 aMass = gMC->TrackMass();
877 if ((gMC->Etot() - aMass) > kEkinMinStep) {
878
879 // The energy loss according to Bethe Bloch
880 iPdg = TMath::Abs(gMC->TrackPid());
881 if ((iPdg != kPdgElectron) ||
882 ((iPdg == kPdgElectron) &&
883 (pTot < kPTotMaxEl))) {
884 gMC->TrackMomentum(mom);
885 pTot = mom.Rho();
886 betaGamma = pTot / aMass;
887 pp = kPrim * BetheBloch(betaGamma);
888 // Take charge > 1 into account
889 charge = gMC->TrackCharge();
890 if (TMath::Abs(charge) > 1) {
891 pp = pp * charge*charge;
f73816f5 892 }
ce0d6231 893 }
894 else {
895 // Electrons above 20 Mev/c are at the plateau
896 pp = kPrim * kPlateau;
897 }
f73816f5 898
ce0d6231 899 if (pp > 0.0) {
900 do {
901 gMC->GetRandom()->RndmArray(1,random);
030b4415 902 }
ce0d6231 903 while ((random[0] == 1.0) ||
904 (random[0] == 0.0));
905 stepSize = - TMath::Log(random[0]) / pp;
906 gMC->SetMaxStep(stepSize);
907 }
030b4415 908
5c7f4665 909 }
030b4415 910
d3f347ff 911 }
030b4415 912
5c7f4665 913 }
030b4415 914
5c7f4665 915}
916
a328fff9 917//_____________________________________________________________________________
918void AliTRDv1::StepManagerFixedStep()
919{
920 //
921 // Slow simulator. Every charged track produces electron cluster as hits
922 // along its path across the drift volume. The step size is fixed in
923 // this version of the step manager.
924 //
f2979d08 925 // Works for Xe/CO2 as well as Ar/CO2
926 //
a328fff9 927
ce0d6231 928 // PDG code electron
929 const Int_t kPdgElectron = 11;
930
053767a4 931 Int_t layer = 0;
932 Int_t stack = 0;
933 Int_t sector = 0;
934 Int_t det = 0;
a328fff9 935 Int_t qTot;
936
937 Float_t hits[3];
938 Double_t eDep;
939
940 Bool_t drRegion = kFALSE;
941 Bool_t amRegion = kFALSE;
942
2c8bf4aa 943 TString cIdPath;
944 Char_t cIdSector[3];
945 cIdSector[2] = 0;
946
a328fff9 947 TString cIdCurrent;
948 TString cIdSensDr = "J";
949 TString cIdSensAm = "K";
950 Char_t cIdChamber[3];
2c8bf4aa 951 cIdChamber[2] = 0;
a328fff9 952
030b4415 953 TLorentzVector pos;
954 TLorentzVector mom;
a328fff9 955
053767a4 956 const Int_t kNlayer = AliTRDgeometry::Nlayer();
957 const Int_t kNstack = AliTRDgeometry::Nstack();
958 const Int_t kNdetsec = kNlayer * kNstack;
a328fff9 959
030b4415 960 const Double_t kBig = 1.0e+12;
a328fff9 961 const Float_t kEkinMinStep = 1.0e-5; // Minimum energy for the step size adjustment
962
963 // Set the maximum step size to a very large number for all
964 // neutral particles and those outside the driftvolume
965 gMC->SetMaxStep(kBig);
966
967 // If not charged track or already stopped or disappeared, just return.
968 if ((!gMC->TrackCharge()) ||
ce0d6231 969 gMC->IsTrackDisappeared()) {
970 return;
971 }
a328fff9 972
973 // Inside a sensitive volume?
974 cIdCurrent = gMC->CurrentVolName();
975
ce0d6231 976 if (cIdSensDr == cIdCurrent[1]) {
977 drRegion = kTRUE;
978 }
979 if (cIdSensAm == cIdCurrent[1]) {
980 amRegion = kTRUE;
981 }
a328fff9 982
030b4415 983 if ((!drRegion) &&
984 (!amRegion)) {
985 return;
986 }
a328fff9 987
988 // The hit coordinates and charge
989 gMC->TrackPosition(pos);
990 hits[0] = pos[0];
991 hits[1] = pos[1];
992 hits[2] = pos[2];
993
2c8bf4aa 994 // The sector number (0 - 17), according to standard coordinate system
995 cIdPath = gGeoManager->GetPath();
996 cIdSector[0] = cIdPath[21];
997 cIdSector[1] = cIdPath[22];
053767a4 998 sector = atoi(cIdSector);
a328fff9 999
1000 // The plane and chamber number
030b4415 1001 cIdChamber[0] = cIdCurrent[2];
1002 cIdChamber[1] = cIdCurrent[3];
a328fff9 1003 Int_t idChamber = (atoi(cIdChamber) % kNdetsec);
053767a4 1004 stack = ((Int_t) idChamber / kNlayer);
1005 layer = ((Int_t) idChamber % kNlayer);
e0d47c25 1006
030b4415 1007 // The detector number
053767a4 1008 det = fGeometry->GetDetector(layer,stack,sector);
030b4415 1009
25ca55ce 1010 // 0: InFlight 1:Entering 2:Exiting
030b4415 1011 Int_t trkStat = 0;
a328fff9 1012
1013 // Special hits only in the drift region
ce0d6231 1014 if ((drRegion) &&
1015 (gMC->IsTrackEntering())) {
a328fff9 1016
ce0d6231 1017 // Create a track reference at the entrance of each
1018 // chamber that contains the momentum components of the particle
1019 gMC->TrackMomentum(mom);
e6add757 1020 AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber(), AliTrackReference::kTRD);
ce0d6231 1021 trkStat = 1;
a328fff9 1022
ce0d6231 1023 // Create the hits from TR photons if electron/positron is
1024 // entering the drift volume
1025 if ((fTR) &&
1026 (TMath::Abs(gMC->TrackPid()) == kPdgElectron)) {
f2979d08 1027 CreateTRhit(det);
030b4415 1028 }
a328fff9 1029
ce0d6231 1030 }
1031 else if ((amRegion) &&
1032 (gMC->IsTrackExiting())) {
1033
1034 // Create a track reference at the exit of each
1035 // chamber that contains the momentum components of the particle
1036 gMC->TrackMomentum(mom);
e6add757 1037 AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber(), AliTrackReference::kTRD);
ce0d6231 1038 trkStat = 2;
1039
a328fff9 1040 }
1041
1042 // Calculate the charge according to GEANT Edep
1043 // Create a new dEdx hit
1044 eDep = TMath::Max(gMC->Edep(),0.0) * 1.0e+09;
f2979d08 1045 qTot = (Int_t) (eDep / fWion);
ce0d6231 1046 if ((qTot) ||
1047 (trkStat)) {
1048 AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber()
1049 ,det
1050 ,hits
1051 ,qTot
25ca55ce 1052 ,gMC->TrackTime()*1.0e06
ce0d6231 1053 ,drRegion);
1054 }
a328fff9 1055
1056 // Set Maximum Step Size
1057 // Produce only one hit if Ekin is below cutoff
030b4415 1058 if ((gMC->Etot() - gMC->TrackMass()) < kEkinMinStep) {
1059 return;
1060 }
a328fff9 1061 gMC->SetMaxStep(fStepSize);
1062
1063}
1064
5c7f4665 1065//_____________________________________________________________________________
1066Double_t AliTRDv1::BetheBloch(Double_t bg)
1067{
1068 //
1069 // Parametrization of the Bethe-Bloch-curve
1070 // The parametrization is the same as for the TPC and is taken from Lehrhaus.
1071 //
1072
1073 // This parameters have been adjusted to averaged values from GEANT
f57bb418 1074 const Double_t kP1 = 7.17960e-02;
1075 const Double_t kP2 = 8.54196;
1076 const Double_t kP3 = 1.38065e-06;
1077 const Double_t kP4 = 5.30972;
1078 const Double_t kP5 = 2.83798;
5c7f4665 1079
f73816f5 1080 // Lower cutoff of the Bethe-Bloch-curve to limit step sizes
1081 const Double_t kBgMin = 0.8;
1082 const Double_t kBBMax = 6.83298;
f73816f5 1083
1084 if (bg > kBgMin) {
030b4415 1085 Double_t yy = bg / TMath::Sqrt(1.0 + bg*bg);
5c7f4665 1086 Double_t aa = TMath::Power(yy,kP4);
030b4415 1087 Double_t bb = TMath::Power((1.0/bg),kP5);
5c7f4665 1088 bb = TMath::Log(kP3 + bb);
030b4415 1089 return ((kP2 - aa - bb) * kP1 / aa);
5c7f4665 1090 }
f73816f5 1091 else {
1092 return kBBMax;
1093 }
d3f347ff 1094
fe4da5cc 1095}
5c7f4665 1096
a328fff9 1097//_____________________________________________________________________________
c4214bc0 1098Double_t AliTRDv1::BetheBlochGeant(Double_t bg)
a328fff9 1099{
1100 //
1101 // Return dN/dx (number of primary collisions per centimeter)
1102 // for given beta*gamma factor.
1103 //
1104 // Implemented by K.Oyama according to GEANT 3 parametrization shown in
1105 // A.Andronic's webpage: http://www-alice.gsi.de/trd/papers/dedx/dedx.html
1106 // This must be used as a set with IntSpecGeant.
1107 //
1108
030b4415 1109 Int_t i = 0;
a328fff9 1110
030b4415 1111 Double_t arrG[20] = { 1.100000, 1.200000, 1.300000, 1.500000
1112 , 1.800000, 2.000000, 2.500000, 3.000000
1113 , 4.000000, 7.000000, 10.000000, 20.000000
1114 , 40.000000, 70.000000, 100.000000, 300.000000
1115 , 600.000000, 1000.000000, 3000.000000, 10000.000000 };
a328fff9 1116
030b4415 1117 Double_t arrNC[20] = { 75.009056, 45.508083, 35.299252, 27.116327
1118 , 22.734999, 21.411915, 19.934095, 19.449375
1119 , 19.344431, 20.185553, 21.027925, 22.912676
1120 , 24.933352, 26.504053, 27.387468, 29.566597
1121 , 30.353779, 30.787134, 31.129285, 31.157350 };
1122
1123 // Betagamma to gamma
1124 Double_t g = TMath::Sqrt(1.0 + bg*bg);
a328fff9 1125
1126 // Find the index just before the point we need.
030b4415 1127 for (i = 0; i < 18; i++) {
1128 if ((arrG[i] < g) &&
1129 (arrG[i+1] > g)) {
a328fff9 1130 break;
030b4415 1131 }
1132 }
a328fff9 1133
1134 // Simple interpolation.
030b4415 1135 Double_t pp = ((arrNC[i+1] - arrNC[i]) / (arrG[i+1] - arrG[i]))
1136 * (g - arrG[i]) + arrNC[i];
a328fff9 1137
030b4415 1138 return pp;
a328fff9 1139
1140}
1141
5c7f4665 1142//_____________________________________________________________________________
1143Double_t Ermilova(Double_t *x, Double_t *)
1144{
1145 //
1146 // Calculates the delta-ray energy distribution according to Ermilova.
1147 // Logarithmic scale !
1148 //
1149
1150 Double_t energy;
1151 Double_t dpos;
1152 Double_t dnde;
1153
030b4415 1154 Int_t pos1;
1155 Int_t pos2;
5c7f4665 1156
8230f242 1157 const Int_t kNv = 31;
5c7f4665 1158
030b4415 1159 Float_t vxe[kNv] = { 2.3026, 2.9957, 3.4012, 3.6889, 3.9120
1160 , 4.0943, 4.2485, 4.3820, 4.4998, 4.6052
1161 , 4.7005, 5.0752, 5.2983, 5.7038, 5.9915
1162 , 6.2146, 6.5221, 6.9078, 7.3132, 7.6009
1163 , 8.0064, 8.5172, 8.6995, 8.9872, 9.2103
1164 , 9.4727, 9.9035, 10.3735, 10.5966, 10.8198
1165 , 11.5129 };
1166
1167 Float_t vye[kNv] = { 80.0, 31.0, 23.3, 21.1, 21.0
1168 , 20.9, 20.8, 20.0, 16.0, 11.0
1169 , 8.0, 6.0, 5.2, 4.6, 4.0
1170 , 3.5, 3.0, 1.4, 0.67, 0.44
1171 , 0.3, 0.18, 0.12, 0.08, 0.056
1172 , 0.04, 0.023, 0.015, 0.011, 0.01
1173 , 0.004 };
5c7f4665 1174
1175 energy = x[0];
1176
1177 // Find the position
030b4415 1178 pos1 = 0;
1179 pos2 = 0;
5c7f4665 1180 dpos = 0;
1181 do {
1182 dpos = energy - vxe[pos2++];
1183 }
1184 while (dpos > 0);
1185 pos2--;
030b4415 1186 if (pos2 > kNv) {
1187 pos2 = kNv - 1;
1188 }
5c7f4665 1189 pos1 = pos2 - 1;
1190
1191 // Differentiate between the sampling points
1192 dnde = (vye[pos1] - vye[pos2]) / (vxe[pos2] - vxe[pos1]);
1193
1194 return dnde;
1195
1196}
a328fff9 1197
1198//_____________________________________________________________________________
1199Double_t IntSpecGeant(Double_t *x, Double_t *)
1200{
1201 //
1202 // Integrated spectrum from Geant3
1203 //
1204
96efaf83 1205 const Int_t npts = 83;
030b4415 1206 Double_t arre[npts] = { 2.421257, 2.483278, 2.534301, 2.592230
1207 , 2.672067, 2.813299, 3.015059, 3.216819
1208 , 3.418579, 3.620338, 3.868209, 3.920198
1209 , 3.978284, 4.063923, 4.186264, 4.308605
1210 , 4.430946, 4.553288, 4.724261, 4.837736
1211 , 4.999842, 5.161949, 5.324056, 5.486163
1212 , 5.679688, 5.752998, 5.857728, 5.962457
1213 , 6.067185, 6.171914, 6.315653, 6.393674
1214 , 6.471694, 6.539689, 6.597658, 6.655627
1215 , 6.710957, 6.763648, 6.816338, 6.876198
1216 , 6.943227, 7.010257, 7.106285, 7.252151
1217 , 7.460531, 7.668911, 7.877290, 8.085670
1218 , 8.302979, 8.353585, 8.413120, 8.483500
1219 , 8.541030, 8.592857, 8.668865, 8.820485
1220 , 9.037086, 9.253686, 9.470286, 9.686887
1221 , 9.930838, 9.994655, 10.085822, 10.176990
1222 , 10.268158, 10.359325, 10.503614, 10.627565
1223 , 10.804637, 10.981709, 11.158781, 11.335854
1224 , 11.593397, 11.781165, 12.049404, 12.317644
1225 , 12.585884, 12.854123, 14.278421, 16.975889
1226 , 20.829416, 24.682943, 28.536469 };
1227
1228 Double_t arrdnde[npts] = { 10.960000, 10.960000, 10.359500, 9.811340
1229 , 9.1601500, 8.206670, 6.919630, 5.655430
1230 , 4.6221300, 3.777610, 3.019560, 2.591950
1231 , 2.5414600, 2.712920, 3.327460, 4.928240
1232 , 7.6185300, 10.966700, 12.225800, 8.094750
1233 , 3.3586900, 1.553650, 1.209600, 1.263840
1234 , 1.3241100, 1.312140, 1.255130, 1.165770
1235 , 1.0594500, 0.945450, 0.813231, 0.699837
1236 , 0.6235580, 2.260990, 2.968350, 2.240320
1237 , 1.7988300, 1.553300, 1.432070, 1.535520
1238 , 1.4429900, 1.247990, 1.050750, 0.829549
1239 , 0.5900280, 0.395897, 0.268741, 0.185320
1240 , 0.1292120, 0.103545, 0.0949525, 0.101535
1241 , 0.1276380, 0.134216, 0.123816, 0.104557
1242 , 0.0751843, 0.0521745, 0.0373546, 0.0275391
1243 , 0.0204713, 0.0169234, 0.0154552, 0.0139194
1244 , 0.0125592, 0.0113638, 0.0107354, 0.0102137
1245 , 0.00845984, 0.00683338, 0.00556836, 0.00456874
1246 , 0.0036227, 0.00285991, 0.00226664, 0.00172234
1247 , 0.00131226, 0.00100284, 0.000465492, 7.26607e-05
1248 , 3.63304e-06, 0.0000000, 0.0000000 };
1249
1250 Int_t i;
a328fff9 1251 Double_t energy = x[0];
a328fff9 1252
225211c8 1253 if (energy >= arre[npts-1]) {
1254 return 0.0;
1255 }
1256
030b4415 1257 for (i = 0; i < npts; i++) {
1258 if (energy < arre[i]) {
1259 break;
1260 }
1261 }
a328fff9 1262
030b4415 1263 if (i == 0) {
1264 AliErrorGeneral("AliTRDv1::IntSpecGeant","Given energy value is too small or zero");
1265 }
a328fff9 1266
f57bb418 1267 return arrdnde[i];
a328fff9 1268
1269}