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