CbmRoot
CbmTrdRadiator.cxx
Go to the documentation of this file.
1 // ------------------------------------------------------------------------
2 // ----- CbmTrdRadiator source file -----
3 // ----- Created 10/11/04 by M.Kalisky -----
4 // ------------------------------------------------------------------------
5 #include "CbmTrdRadiator.h"
6 
7 #include "CbmTrdGas.h" // for CbmTrdGas
8 #include "CbmTrdPoint.h" // for CbmTrdPoint
9 
10 #include <FairLogger.h> // for Logger, LOG
11 
12 #include <TFile.h> // for TFile, gFile
13 #include <TGeoManager.h> // for TGeoManager, gGeoManager
14 #include <TH1.h> // for TH1D
15 #include <TMath.h> // for Exp, Sqrt, Cos, Pi
16 #include <TMathBase.h> // for Abs
17 #include <TRandom.h> // for TRandom
18 #include <TRandom3.h> // for gRandom
19 #include <TVector3.h> // for TVector3
20 
21 #include <iomanip> // for setprecision, __iom_t5
22 #include <stdio.h> // for printf, sprintf
23 #include <stdlib.h> // for exit
24 
25 using std::setprecision;
26 
27 
28 // ----- Default constructor ---------------------------------------
29 CbmTrdRadiator::CbmTrdRadiator() : CbmTrdRadiator(kTRUE, 130, 0.0013, 0.02) {}
30 //-----------------------------------------------------------------------------
31 
32 // ----- Constructor --------------------------------------------------
34  Int_t Nfoils,
35  Float_t FoilThick,
36  Float_t GapThick)
37  : CbmTrdRadiator(SimpleTR, Nfoils, FoilThick, GapThick, "polyethylen", "") {}
38 //-----------------------------------------------------------------------------
39 
40 // ----- Constructor --------------------------------------------------
42  Int_t Nfoils,
43  Float_t FoilThick,
44  Float_t GapThick,
45  TString material,
46  TString prototype)
47  : fWindowFoil("")
48  , fRadType(prototype)
49  , fDetType(-1)
50  , fFirstPass(kTRUE)
51  , fSimpleTR(SimpleTR)
52  , fNFoils(Nfoils)
53  , fFoilThick(FoilThick)
54  , fGapThick(GapThick)
55  , fFoilMaterial(material)
56  , fGasThick(-1)
57  , fFoilDens(-1.)
58  , fGapDens(-1.)
59  , fFoilOmega(-1.)
60  , fGapOmega(-1.)
61  , fnPhotonCorr(1.0)
62  , fFoilThickCorr(1.)
63  , fGapThickCorr(1.)
64  , fGasThickCorr(1.)
65  , fnTRprod(-1)
66  , fWinDens(-1.)
67  , fWinThick(-1.)
68  , fCom1(-1.)
69  , fCom2(-1.)
70  , fSpBinWidth((Float_t) fSpRange / (Float_t) fSpNBins)
71  , fSigma(nullptr)
72  , fSigmaWin(nullptr)
73  , fSigmaDet(nullptr)
74  , fSpectrum(nullptr)
75  , fWinSpectrum(nullptr)
76  , fDetSpectrumA(nullptr)
77  , fDetSpectrum(nullptr)
78  , fTrackMomentum(nullptr)
79  , fFinal()
80  , fnTRabs()
81  , fnTRab(-1.)
82  , fELoss(-1.)
83  , fMom(-1., -1., -1.) {
84  for (Int_t i = 0; i < fNMom; i++) {
85  fFinal[i] = nullptr;
86  }
87  //Init();
89 }
90 //-----------------------------------------------------------------------------
91 
92 // ----- Constructor --------------------------------------------------
93 CbmTrdRadiator::CbmTrdRadiator(Bool_t SimpleTR, TString prototype)
94  : fWindowFoil("Kapton")
95  , //which is the gas window of MS2012 prototypes. Anything else will result in a mylar gas window -> MuBu default
96  fRadType(prototype)
97  , fDetType(-1)
98  , fFirstPass(kTRUE)
99  , fSimpleTR(SimpleTR)
100  , fNFoils(130)
101  , fFoilThick(0.0013)
102  , fGapThick(0.02)
103  , fFoilMaterial("")
104  , fGasThick(-1)
105  , fFoilDens(-1.)
106  , fGapDens(-1.)
107  , fFoilOmega(-1.)
108  , fGapOmega(-1.)
109  , fnPhotonCorr(1.0)
110  , fFoilThickCorr(1.)
111  , fGapThickCorr(1.)
112  , fGasThickCorr(1.)
113  , fnTRprod(-1)
114  , fWinDens(-1.)
115  , fWinThick(-1.)
116  , fCom1(-1.)
117  , fCom2(-1.)
118  , fSpBinWidth((Float_t) fSpRange / (Float_t) fSpNBins)
119  , fSigma(nullptr)
120  , fSigmaWin(nullptr)
121  , fSigmaDet(nullptr)
122  , fSpectrum(nullptr)
123  , fWinSpectrum(nullptr)
124  , fDetSpectrumA(nullptr)
125  , fDetSpectrum(nullptr)
126  , fTrackMomentum(nullptr)
127  , fFinal()
128  , fnTRabs()
129  , fnTRab(-1.)
130  , fELoss(-1.)
131  , fMom(-1., -1., -1.) {
132  for (Int_t i = 0; i < fNMom; i++) {
133  fFinal[i] = nullptr;
134  }
135 
136  //Init(SimpleTR, prototype);
137  //SetRadPrototype(prototype);
139 }
140 //-----------------------------------------------------------------------------
141 
142 // ----- Destructor ----------------------------------------------------
144  // FairRootManager *fManager = FairRootManager::Instance();
145  // fManager->Write();
146  if (fSpectrum) {
147  LOG(info) << " -I DELETING fSpectrum ";
148  delete fSpectrum;
149  }
150  if (fWinSpectrum) delete fWinSpectrum;
151  if (fDetSpectrum) delete fDetSpectrum;
152  if (fDetSpectrumA) delete fDetSpectrumA;
153 
154  for (Int_t i = 0; i < fNMom; i++) {
155  if (fFinal[i]) delete fFinal[i];
156  }
157 }
158 
159 // ----- Create Histogramms --------------------------------------------------
161 
162  // Create the needed histograms
163 
164  fSpBinWidth = (Float_t) fSpRange / (Float_t) fSpNBins;
165 
166  Float_t SpLower = 1.0 - 0.5 * fSpBinWidth;
167  Float_t SpUpper = SpLower + (Float_t) fSpRange;
168 
169  Char_t name[50];
170 
171  if (fSpectrum) delete fSpectrum;
172  fSpectrum = new TH1D("fSpectrum", "TR spectrum", fSpNBins, SpLower, SpUpper);
173 
174  if (fWinSpectrum) delete fWinSpectrum;
175  fWinSpectrum = new TH1D(
176  "fWinSpectrum", "TR spectrum in Mylar", fSpNBins, SpLower, SpUpper);
177 
178  if (fDetSpectrum) delete fDetSpectrum;
179  fDetSpectrum = new TH1D("fDetSpectrum",
180  "TR spectrum escaped from detector",
181  fSpNBins,
182  SpLower,
183  SpUpper);
184 
185  if (fDetSpectrumA) delete fDetSpectrumA;
186  fDetSpectrumA = new TH1D("fDetSpectrumA",
187  "TR spectrum absorbed in detector",
188  fSpNBins,
189  SpLower,
190  SpUpper);
191 
192  for (Int_t i = 0; i < fNMom; i++) {
193  sprintf(name, "fFinal%d", i + 1);
194  //LOG(info) <<"name : "<<name;
195  if (fFinal[i]) delete fFinal[i];
196  fFinal[i] = new TH1D(name, name, fSpNBins, SpLower, SpUpper);
197  }
198 }
199 
200 // ----- Init function ----------------------------------------------------
201 void CbmTrdRadiator::Init(Bool_t SimpleTR,
202  Int_t Nfoils,
203  Float_t FoilThick,
204  Float_t GapThick,
205  TString material = "polyethylen") {
206 
207  LOG(info) << "================CbmTrdRadiator===============";
208  // Set initial parameters defining the radiator
209 
210  fNFoils = Nfoils;
211  fGapThick = GapThick;
212  fFoilThick = FoilThick;
213  fSimpleTR = SimpleTR;
214  fnPhotonCorr = 1.0;
215 
216 
218 
219  LOG(info) << "CbmTrdRadiator::Init : Foil material : " << material;
220  LOG(info) << "CbmTrdRadiator::Init : Nr. of foils : " << Nfoils;
221  LOG(info) << "CbmTrdRadiator::Init : Foil thickness : "
222  << setprecision(4) << FoilThick << " cm";
223  LOG(info) << "CbmTrdRadiator::Init : Gap thickness : " << GapThick
224  << " cm";
225  LOG(info) << "CbmTrdRadiator::Init : Simple TR production: "
226  << setprecision(2) << SimpleTR;
227 
228  // material dependent parameters for the radiator material
229  // (polyethylen)
230  if (material == "" || material == "polyethylen") {
231  fFoilDens = 0.92; // [g/cm3 ]
232  fFoilOmega = 20.9; //plasma frequency ( from Anton )
233  } else if (material == "pefoam20") {
234  // (polyethylen)
235  fFoilDens = 0.90; // [g/cm3 ]
236  fFoilOmega = 20.64; //plasma frequency ( from Cyrano )
237  } else if (material == "pefiber") {
238  // (polyethylen)
239  fFoilDens = 0.90; // [g/cm3 ]
240  fFoilOmega = 20.64; //plasma frequency ( from Cyrano )
241  } else if (material == "mylar") {
242  // (Mylar)
243  fFoilDens = 1.393; // [g/cm3 ]
244  fFoilOmega = 24.53; //plasma frequency ( from Cyrano )
245  } else if (material == "pokalon") {
246  // (Pokalon)
247  fFoilDens = 1.150; // [g/cm3 ]
248  fFoilOmega = 22.50; //plasma frequency ( from Cyrano )
249  } else {
250  // (polyethylen)
251  material = "polyethylen";
252  fFoilDens = 0.92; // [g/cm3 ]
253  fFoilOmega = 20.9; //plasma frequency ( from Anton )
254  }
255  fFoilMaterial = material;
256  // material dependent parameters for the gas between the foils of the
257  // radiator
258  // changed gap density at 16.04.07 FU
259  // density of air is 0.001205, 0.00186
260  //fGapDens = 0.00186; // [g/cm3]
261  fGapDens = 0.001205; // [g/cm3]
262  fGapOmega = 0.7; //plasma frequency ( from Anton )
263 
264  // foil between radiator and gas chamber
265  // TODO: implement kapton foil also in trd geometry
266  fWinDens = 1.39; // [g/cm3]
267  fWinThick = 0.0025; // [cm]
268 
269 
270  // Get all the gas properties from CbmTrdGas
271  CbmTrdGas* fTrdGas = CbmTrdGas::Instance();
272  if (fTrdGas == 0) {
273  fTrdGas = new CbmTrdGas();
274  fTrdGas->Init();
275  }
276 
277  fCom2 = fTrdGas->GetNobleGas();
278  fCom1 = fTrdGas->GetCO2();
279  fDetType = fTrdGas->GetDetType();
280  fGasThick = fTrdGas->GetGasThick();
281 
282  LOG(info) << "CbmTrdRadiator::Init : Detector noble gas : "
283  << fTrdGas->GetNobleGas();
284  LOG(info) << "CbmTrdRadiator::Init : Detector cruncher gas: "
285  << fTrdGas->GetCO2();
286  LOG(info) << "CbmTrdRadiator::Init : Detector type : "
287  << fTrdGas->GetDetType();
288  LOG(info) << "CbmTrdRadiator::Init : Detector gas thick. : "
289  << fTrdGas->GetGasThick() << " cm";
290 
291  //LOG(info) << "************* End of Trd Radiator Init **************";
292 
293  // If simplified version is used one didn't calculate the TR
294  // for each CbmTrdPoint in full glory, but create at startup
295  // some histograms for several momenta which are later used to
296  // get the TR much faster.
297  if (fSimpleTR == kTRUE) {
298 
299  ProduceSpectra();
300 
301  TFile* oldfile = gFile;
302  TFile* f1 = new TFile("TRhistos.root", "recreate");
303  f1->cd();
304 
305  for (Int_t i = 0; i < fNMom; i++) {
306  fFinal[i]->Write();
307  }
308 
309  f1->Close();
310  f1->Delete();
311  gFile = oldfile;
312  }
313 }
314 //----------------------------------------------------------------------------
315 
316 // ----- Init function ----------------------------------------------------
317 void CbmTrdRadiator::Init(Bool_t SimpleTR,
318  Int_t Nfoils,
319  Float_t FoilThick,
320  Float_t GapThick) {
321 
322  LOG(info) << "================CbmTrdRadiator===============";
323  // Set initial parameters defining the radiator
324 
325  fNFoils = Nfoils;
326  fGapThick = GapThick;
327  fFoilThick = FoilThick;
328  fSimpleTR = SimpleTR;
329  fnPhotonCorr = 1.0;
330 
332 
333  //LOG(info) << "CbmTrdRadiator::Init : Foil material : " << material;
334  LOG(info) << "CbmTrdRadiator::Init : Nr. of foils : " << Nfoils;
335  LOG(info) << "CbmTrdRadiator::Init : Foil thickness : "
336  << setprecision(4) << FoilThick << " cm";
337  LOG(info) << "CbmTrdRadiator::Init : Gap thickness : " << GapThick
338  << " cm";
339  LOG(info) << "CbmTrdRadiator::Init : Simple TR production: "
340  << setprecision(2) << SimpleTR;
341 
342 
343  // material dependend parameters for the radiator material
344  // (polyethylen)
345  fFoilDens = 0.92; // [g/cm3 ]
346  fFoilOmega = 20.9; //plasma frequency ( from Anton )
347 
348  // material dependent parameters for the gas between the foils of the
349  // radiator
350  // changed gap density at 16.04.07 FU
351  // density of air is 0.001205, 0.00186
352  //fGapDens = 0.00186; // [g/cm3]
353  fGapDens = 0.001205; // [g/cm3]
354  fGapOmega = 0.7; //plasma frequency ( from Anton )
355 
356  // foil between radiator and gas chamber
357  // TODO: implement kapton foil also in trd geometry
358  fWinDens = 1.39; // [g/cm3]
359  fWinThick = 0.0025; // [cm]
360 
361 
362  // Get all the gas properties from CbmTrdGas
363  CbmTrdGas* fTrdGas = CbmTrdGas::Instance();
364  if (fTrdGas == 0) {
365  fTrdGas = new CbmTrdGas();
366  fTrdGas->Init();
367  }
368 
369  fCom2 = fTrdGas->GetNobleGas();
370  fCom1 = fTrdGas->GetCO2();
371  fDetType = fTrdGas->GetDetType();
372  fGasThick = fTrdGas->GetGasThick();
373  /*
374  LOG(info) << "Detector noble gas : "<< fTrdGas->GetNobleGas();
375  LOG(info) << "Detector cruncher gas: "<< fTrdGas->GetCO2();
376  LOG(info) << "Detector type : "<< fTrdGas->GetDetType();
377  LOG(info) << "Detector gas thick. : "<< fTrdGas->GetGasThick() <<" cm";
378  */
379  //LOG(info) << "************* End of Trd Radiator Init **************";
380 
381  // If simplified version is used one didn't calculate the TR
382  // for each CbmTrdPoint in full glory, but create at startup
383  // some histograms for several momenta which are later used to
384  // get the TR much faster.
385  if (fSimpleTR == kTRUE) {
386 
387  ProduceSpectra();
388 
389  TFile* oldfile = gFile;
390  TFile* f1 = new TFile("TRhistos.root", "recreate");
391  f1->cd();
392 
393  for (Int_t i = 0; i < fNMom; i++) {
394  fFinal[i]->Write();
395  }
396 
397  f1->Close();
398  f1->Delete();
399  gFile = oldfile;
400  }
401 }
402 //----------------------------------------------------------------------------
403 
404 // ----- Init function ----------------------------------------------------
406  LOG(info) << "================CbmTrdRadiator===============";
407  TString material;
409  if (fRadType == "") {
410  fFoilMaterial = "polyethylen";
411  fnPhotonCorr = 1.0;
412 
413  /*
414  LOG(info) << "CbmTrdRadiator::Init : Foil material : " << fFoilMaterial;
415  //LOG(info) << "CbmTrdRadiator::Init : Nr. of foils : " << fNfoils;
416  LOG(info) << "CbmTrdRadiator::Init : Foil thickness : " << setprecision(4) << fFoilThick <<" cm";
417  LOG(info) << "CbmTrdRadiator::Init : Gap thickness : " << fGapThick << " cm";
418  LOG(info) << "CbmTrdRadiator::Init : Simple TR production: " << setprecision(2)<< fSimpleTR;
419  */
420  // material dependend parameters for the radiator material
421  // (polyethylen)
422  fFoilDens = 0.92; // [g/cm3 ]
423  fFoilOmega = 20.9; //plasma frequency ( from Anton )
424 
425  // material dependent parameters for the gas between the foils of the
426  // radiator
427  // changed gap density at 16.04.07 FU
428  // density of air is 0.001205, 0.00186
429  //fGapDens = 0.00186; // [g/cm3]
430  fGapDens = 0.001205; // [g/cm3]
431  fGapOmega = 0.7; //plasma frequency ( from Anton )
432 
433  // foil between radiator and gas chamber
434  // TODO: implement kapton foil also in trd geometry
435  fWinDens = 1.39; // [g/cm3]
436  fWinThick = 0.0025; // [cm]
437  } else {
438  if (fRadType == "default") {
439  material = "polyethylen";
440  fNFoils = 130;
441  fGapThick = 0.02;
442  fFoilThick = 0.0013;
443  fnPhotonCorr = 1.0;
444  } else if (fRadType == "A") {
445  material = "polyethylen";
446  fNFoils = 200;
447  fGapThick = 700 * 0.0001;
448  fFoilThick = 15 * 0.0001;
449  fnPhotonCorr = 0.40;
450  } else if (fRadType == "Bshort" || fRadType == "Kshort") {
451  material = "pokalon";
452  fNFoils = 100;
453  fGapThick = 0.07;
454  fFoilThick = 0.0024;
455  fnPhotonCorr = 0.65;
456  } else if (fRadType == "B" || fRadType == "K") {
457  material = "pokalon";
458  fNFoils = 250;
459  fGapThick = 0.07;
460  fFoilThick = 0.0024;
461  fnPhotonCorr = 0.65;
462  } else if (fRadType == "B++" || fRadType == "K++") {
463  material = "pokalon";
464  fNFoils = 350;
465  fGapThick = 0.07;
466  fFoilThick = 0.0024;
467  fnPhotonCorr = 0.65;
468  } else if (fRadType == "C") {
469  material = "polyethylen";
470  fNFoils = 200;
471  fGapThick = 700 * 0.0001;
472  fFoilThick = 15 * 0.0001;
473  fnPhotonCorr = 0.30;
474  } else if (fRadType == "D") {
475  material = "polyethylen";
476  fNFoils = 100;
477  fGapThick = 500 * 0.0001;
478  fFoilThick = 15 * 0.0001;
479  fnPhotonCorr = 0.45;
480  } else if (fRadType == "E") {
481  material = "polyethylen";
482  fNFoils = 120;
483  fGapThick = 500 * 0.0001;
484  fFoilThick = 20 * 0.0001;
485  fnPhotonCorr = 0.60;
486  } else if (fRadType == "F") {
487  material = "polyethylen";
488  fNFoils = 220;
489  fGapThick = 250 * 0.0001;
490  fFoilThick = 20 * 0.0001;
491  fnPhotonCorr = 0.65;
492  } else if (fRadType == "G30") {
493  material = "pefiber";
494  fNFoils = 200;
495  fGapThick = 700 * 0.0001;
496  fFoilThick = 15 * 0.0001;
497  fnPhotonCorr = 0.65;
498  } else if (fRadType == "H") {
499  material = "pefoam20";
500  fNFoils = 200;
501  fGapThick = 700 * 0.0001;
502  fFoilThick = 15 * 0.0001;
503  fnPhotonCorr = 0.65;
504  } else if (fRadType == "H++") {
505  material = "pefoam20";
506  fNFoils = 200;
507  fGapThick = 700 * 0.0001;
508  fFoilThick = 15 * 0.0001;
509  fnPhotonCorr = 0.65;
510  } else {
511  LOG(warn) << "CbmTrdRadiator::Init : *********** Radiator prototype not "
512  "known ***********";
513  LOG(warn) << "CbmTrdRadiator::Init : switch to default "
514  "parameters ";
515  material = "polyethylen";
516  fNFoils = 130;
517  fGapThick = 0.02;
518  fFoilThick = 0.0013;
519  fnPhotonCorr = 1.0;
520  }
521  fFoilMaterial = material;
522  if (material == "" || material == "polyethylen") {
523  fFoilDens = 0.92; // [g/cm3 ]
524  fFoilOmega = 20.9; //plasma frequency ( from Anton )
525  } else if (material == "pefoam20") {
526  // (polyethylen)
527  fFoilDens = 0.90; // [g/cm3 ]
528  fFoilOmega = 20.64; //plasma frequency ( from Cyrano )
529  } else if (material == "pefiber") {
530  // (polyethylen)
531  fFoilDens = 0.90; // [g/cm3 ]
532  fFoilOmega = 20.64; //plasma frequency ( from Cyrano )
533  } else if (material == "mylar") {
534  // (Mylar)
535  fFoilDens = 1.393; // [g/cm3 ]
536  fFoilOmega = 24.53; //plasma frequency ( from Cyrano )
537  } else if (material == "pokalon") {
538  // (Pokalon)
539  fFoilDens = 1.150; // [g/cm3 ]
540  fFoilOmega = 22.50; //plasma frequency ( from Cyrano )
541  } else {
542  // (polyethylen)
543  fFoilDens = 0.92; // [g/cm3 ]
544  fFoilOmega = 20.9; //plasma frequency ( from Anton )
545  }
546  fGapDens = 0.001205; // [g/cm3]
547  fGapOmega = 0.7; //plasma frequency ( from Anton )
548  // foil between radiator and gas chamber
549  fWinDens =
550  1.42; //[g/cm3] for Kapton -> Alu is added only in case of kapton !!!
551  fWinThick = 0.0025; // [cm]
552 
553  LOG(info) << "CbmTrdRadiator::Init : Prototype : " << fRadType;
554  LOG(info) << "CbmTrdRadiator::Init : Scaling factor : "
555  << fnPhotonCorr;
556  LOG(info) << "CbmTrdRadiator::Init : Foil material : "
557  << fFoilMaterial;
558  LOG(info) << "CbmTrdRadiator::Init : Nr. of foils : " << fNFoils;
559  LOG(info) << "CbmTrdRadiator::Init : Foil thickness : "
560  << setprecision(4) << fFoilThick << " cm";
561  LOG(info) << "CbmTrdRadiator::Init : Gap thickness : " << fGapThick
562  << " cm";
563  LOG(info) << "CbmTrdRadiator::Init : Simple TR production: "
564  << setprecision(2) << fSimpleTR;
565  LOG(info) << "CbmTrdRadiator::Init : Foil plasm. frequ. : " << fFoilOmega
566  << " keV";
567  LOG(info) << "CbmTrdRadiator::Init : Gap plasm. frequ. : " << fGapOmega
568  << " keV";
569  LOG(info) << "CbmTrdRadiator::Init : Window density : " << fWinDens
570  << " g/cm^3";
571  LOG(info) << "CbmTrdRadiator::Init : Window thickness : " << fWinThick
572  << " cm";
573  }
574 
575 
576  // Get all the gas properties from CbmTrdGas
577  CbmTrdGas* fTrdGas = CbmTrdGas::Instance();
578  if (fTrdGas == 0) {
579  fTrdGas = new CbmTrdGas();
580  fTrdGas->Init();
581  }
582 
583  fCom2 = fTrdGas->GetNobleGas();
584  fCom1 = fTrdGas->GetCO2();
585  fDetType = fTrdGas->GetDetType();
586  fGasThick = fTrdGas->GetGasThick();
587 
588  LOG(info) << "CbmTrdRadiator::Init : Detector noble gas : "
589  << fTrdGas->GetNobleGas();
590  LOG(info) << "CbmTrdRadiator::Init : Detector cruncher gas: "
591  << fTrdGas->GetCO2();
592  LOG(info) << "CbmTrdRadiator::Init : Detector type : "
593  << fTrdGas->GetDetType();
594  LOG(info) << "CbmTrdRadiator::Init : Detector gas thick. : "
595  << fTrdGas->GetGasThick() << " cm";
596 
597  //LOG(info) << "************* End of Trd Radiator Init **************";
598  // If simplified version is used one didn't calculate the TR
599  // for each CbmTrdPoint in full glory, but create at startup
600  // some histograms for several momenta which are later used to
601  // get the TR much faster.
602  if (fSimpleTR == kTRUE) {
603 
604  ProduceSpectra();
605 
606  TFile* f1 = new TFile("TRhistos.root", "recreate");
607 
608  for (Int_t i = 0; i < fNMom; i++) {
609  fFinal[i]->Write();
610  }
611  f1->Close();
612  f1->Delete();
613  }
614 }
615 //----------------------------------------------------------------------------
616 
617 //---- Lattice Hit------------------------------------------------------------
619  //printf("---------------------------------------------------\n");
620  Double_t point_in[3];
621  Double_t point_out[3];
622  if (nullptr != point) {
623  point_in[0] = point->GetXIn();
624  point_in[1] = point->GetYIn();
625  point_in[2] = point->GetZIn();
626 
627  point_out[0] = point->GetXOut();
628  point_out[1] = point->GetYOut();
629  point_out[2] = point->GetZOut();
630  Double_t back_direction[3] = {(point_in[0] - point_out[0]),
631  (point_in[1] - point_out[1]),
632  (point_in[2] - point_out[2])};
633  if (back_direction[2] > 0) {
634  LOG(debug2)
635  << "CbmTrdRadiator::LatticeHit: MC-track points towards target!";
636  //for(Int_t i = 0; i < 3; i++)
637  //printf("%i: in:%8.2f out:%8.2f direction:%8.2f\n",i,point_in[i],point_out[i],back_direction[i]);
638  return false;
639  }
640  Double_t trackLength = TMath::Sqrt(back_direction[0] * back_direction[0]
641  + back_direction[1] * back_direction[1]
642  + back_direction[2] * back_direction[2]);
643 
644  trackLength *= 10.; // cm -> mm to get a step width of 1mm
645  //printf("track length:%7.2fmm\n",trackLength);
646  gGeoManager->FindNode((point_out[0] + point_in[0]) / 2,
647  (point_out[1] + point_in[1]) / 2,
648  (point_out[2] + point_in[2]) / 2);
649 
650  Double_t pos[3] = {
651  point_in[0], point_in[1], point_in[2]}; // start at entrance point
652 
653  for (Int_t i = 0; i < 3; i++) {
654  back_direction[i] /= trackLength;
655  //printf("%i: in:%8.2f out:%8.2f start:%8.2f direction:%8.2f\n",i,point_in[i],point_out[i],pos[i],back_direction[i]);
656  }
657  if (TString(gGeoManager->GetPath()).Contains("gas")) {
658  Int_t stepCount = 0;
659  while (
660  /*(!TString(gGeoManager->GetPath()).Contains("lattice") || !TString(gGeoManager->GetPath()).Contains("radiator")) &&*/
661  pos[2] >= point_in[2] - 5
662  && stepCount < 50) {
663  stepCount++;
664  //printf("step%i\n",stepCount);
665  for (Int_t i = 0; i < 3; i++) {
666  pos[i] += back_direction[i];
667  //printf("%8.2f ",pos[i]);
668  }
669  //printf(" %s\n",TString(gGeoManager->GetPath()).Data());
670  gGeoManager->FindNode(pos[0], pos[1], pos[2]);
671  if (TString(gGeoManager->GetPath()).Contains("radiator")) {
672  //printf ("%s false\n",TString(gGeoManager->GetPath()).Data());
673  return false;
674  } else if (TString(gGeoManager->GetPath()).Contains("lattice")) {
675  //printf ("%s true <----------------------------------------\n",TString(gGeoManager->GetPath()).Data());
676  return true;
677  } else if (TString(gGeoManager->GetPath()).Contains("frame")) {
678  //printf ("%s true <----------------------------------------\n",TString(gGeoManager->GetPath()).Data());
679  return true;
680  } else {
681  //printf ("%s\n",TString(gGeoManager->GetPath()).Data());
682  }
683  }
684  } else {
685  LOG(error) << "CbmTrdRadiator::LatticeHit: MC-track not in TRD! Node:"
686  << TString(gGeoManager->GetPath()).Data()
687  << " gGeoManager->MasterToLocal() failed!";
688  return false;
689  }
690  } else {
691  LOG(error) << "CbmTrdRadiator::LatticeHit: CbmTrdPoint == nullptr!";
692  return false;
693  }
694  return kTRUE;
695 }
696 //----------------------------------------------------------------------------
697 // ---- Spectra Production ---------------------------------------------------
699 
700  fTrackMomentum = new Double_t[fNMom];
701 
702  Double_t trackMomentum[fNMom] = {
703  0.1, 0.25, 0.5, 1.0, 1.5, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0, 9.0, 10.0};
704 
705  for (Int_t imom = 0; imom < fNMom; imom++) {
706  fTrackMomentum[imom] = trackMomentum[imom];
707  }
708 
709  for (Int_t i = 0; i < fNMom; i++) {
710 
711  fMom.SetXYZ(0.0, 0.0, fTrackMomentum[i]);
712 
713  ProcessTR();
714  fnTRabs[i] = fnTRab;
715 
716  // copy the content of the current fDetSpectrumA into appropriate
717  // fFinal
718 
719  Float_t tmp = 0;
720  for (Int_t j = 0; j < fSpNBins; j++) {
721  tmp = fDetSpectrumA->GetBinContent(j + 1);
722  fFinal[i]->SetBinContent(j + 1, tmp);
723  }
724  fFinal[i]->SetTitle(
725  Form("%s-momentum%.2f", fFinal[i]->GetTitle(), fTrackMomentum[i]));
726  }
727 }
728 
729 //----------------------------------------------------------------------------
730 
731 // ---- GetTR ---------------------------------------------------------------
732 Float_t CbmTrdRadiator::GetTR(TVector3 Mom) {
733 
734  /*
735  * Simplified simulation of the TR
736  */
737 
738  if (fSimpleTR == kTRUE) {
739 
740  Float_t m = Mom.Mag();
741 
742  // Find correct spectra
743  Int_t i = 0;
744  while (m > fTrackMomentum[i]) {
745  i++;
746  if (i == fNMom) {
747  i--;
748  break;
749  }
750  }
751 
752  ELoss(i);
753 
754  }
755 
756  /*
757  * Full simulation of the TR production
758  */
759  else {
760 
761  fMom = Mom;
762  ProcessTR();
763  }
764 
765  return fELoss;
766 }
767 
768 //----------------------------------------------------------------------------
769 
770 //----- main TR processing function ------------------------------
772 
773  // Compute the angle normalization factor -
774  // for different angles the detector thicknesses are different
775 
776  Float_t fNormPhi = fMom.Mag() / fMom.Pz();
777 
778  // Correct the thickness according to this factor
779 
780  fFoilThickCorr = TMath::Abs(fFoilThick * fNormPhi);
781  fGapThickCorr = TMath::Abs(fGapThick * fNormPhi);
782  fGasThickCorr = TMath::Abs(fGasThick * fNormPhi);
783 
784  fFirstPass = true;
785  // compute the TR spectra in the radiator
786  TRspectrum();
787  // compute the TR spectra in the mylar foil
788  WinTRspectrum();
789  // compute the TR spectra in the detector
790  DetTRspectrum();
791  // Loop over photons and compute the E losses
792  ELoss(-1);
793 
794 
795  if (fDetType == 1) { // if detector is MB type = dual-sided MWPC
796  fFirstPass = false;
797  Float_t energiefirst = fELoss;
798  Float_t fTRAbsfirst = fnTRab;
799 
800  // compute the TR spectra in the mylar foil between the two gas layers
801  WinTRspectrum();
802  // compute the TR spectra in the second halve of the detector
803  DetTRspectrum();
804  // Loop over photons and compute the E losses
805  ELoss(-1);
806  fELoss += energiefirst;
807  fnTRab += fTRAbsfirst;
808  }
809 }
810 //----------------------------------------------------------------------------
811 
812 // ---- Compute the energy losses --------------------------------------------
813 Int_t CbmTrdRadiator::ELoss(Int_t index) {
814 
815  // calculate the energy loss of the TR photons in the
816  // detector gas. fNTRab is the number of TR photons
817  // absorbed in the gas. Since this number is low the
818  // actual number comes from a poisson distribution.
819 
820  Int_t nPhoton;
821  Int_t maxPhoton = 50;
822  Double32_t eLoss = 0;
823  Float_t eTR = 0;
824 
825  if (-1 == index) {
826  if (0 == fDetSpectrumA) {
827  LOG(error) << " -Err :: No input spectra ";
828  fELoss = -1;
829  return 1;
830  }
831  nPhoton = gRandom->Poisson(fnTRab * fnPhotonCorr);
832  if (nPhoton > maxPhoton) nPhoton = maxPhoton;
833  for (Int_t i = 0; i < nPhoton; i++) {
834  // energy of the photon
835  eTR = fDetSpectrumA->GetRandom();
836  eLoss += eTR;
837  }
838  } else {
839  if (0 == fFinal[index]) {
840  LOG(error) << " -Err :: No input spectra ";
841  fELoss = -1;
842  return 1;
843  }
844  nPhoton = gRandom->Poisson(fnTRabs[index] * fnPhotonCorr);
845  if (nPhoton > maxPhoton) nPhoton = maxPhoton;
846  for (Int_t i = 0; i < nPhoton; i++) {
847  // energy of the photon
848  eTR = fFinal[index]->GetRandom();
849  eLoss += eTR;
850  }
851  }
852 
853  // keV->GeV & set the result
854  fELoss = eLoss / 1000000.0;
855 
856  return 0;
857 }
858 
859 //----------------------------------------------------------------------------
860 
861 //----- TR spectrum ---------------------------------------------------
863 
864  // calculate the number of TR photons in the radiator which are not
865  // absorbed in the radiator.
866  // Formulas are tacken from M. Castellano et al.
867  // Computer Physics Communications 61 (1990)
868 
869 
870  // Where does this values come from, put them in the header file
871  const Float_t kAlpha = 0.0072973; // 1/137
872  const Int_t kSumMax = 30;
873 
874  // calculate the gamma of the particle
875  Float_t gamma = GammaF();
876 
877  fSpectrum->Reset();
878  fDetSpectrum->Reset();
879  fDetSpectrumA->Reset();
880 
881  Float_t kappa = fGapThickCorr / fFoilThickCorr;
882 
883  SetSigma(1);
884 
885  Float_t stemp = 0;
886 
887  // Loop over energy
888  // stemp is the total energy of all created photons
889  for (Int_t iBin = 0; iBin < fSpNBins; iBin++) {
890 
891  // keV -> eV
892  Float_t energyeV = (fSpBinWidth * iBin + 1.0) * 1e3;
893 
894  Float_t csFoil = fFoilOmega / energyeV;
895  Float_t csGap = fGapOmega / energyeV;
896 
897  Float_t rho1 = energyeV * fFoilThickCorr * 1e4 * 2.5
898  * (1.0 / (gamma * gamma) + csFoil * csFoil);
899  Float_t rho2 = energyeV * fFoilThickCorr * 1e4 * 2.5
900  * (1.0 / (gamma * gamma) + csGap * csGap);
901 
902  // Calculate the sum over the production angle tetan
903  Float_t sum = 0;
904  for (Int_t iSum = 0; iSum < kSumMax; iSum++) {
905  Float_t tetan = (TMath::Pi() * 2.0 * (iSum + 1) - (rho1 + kappa * rho2))
906  / (kappa + 1.0);
907  if (tetan < 0.0) { tetan = 0.0; }
908  Float_t aux = 1.0 / (rho1 + tetan) - 1.0 / (rho2 + tetan);
909  sum += tetan * (aux * aux) * (1.0 - TMath::Cos(rho1 + tetan));
910  }
911  Float_t conv = 1.0 - TMath::Exp(-fNFoils * fSigma[iBin]);
912 
913  // eV -> keV
914  Float_t energykeV = energyeV * 0.001;
915 
916  // dN / domega
917  Float_t wn =
918  kAlpha * 4.0 / (fSigma[iBin] * (kappa + 1.0)) * conv * sum / energykeV;
919  /*
920  wn = 4.0 * kAlpha
921  / (energykeV * (kappa + 1.0))
922  * (1.0 - TMath::Exp(Double_t(-fNFoils) * fSigma[iBin])) / (1.0 - TMath::Exp(-fSigma[iBin]))
923  * sum;
924  */
925  // save the result
926  fSpectrum->SetBinContent(iBin + 1, wn);
927  // compute the integral
928  stemp += wn;
929  }
930 
931  // <nTR> (binsize corr.)
932  //Float_t nTR0 = stemp * fSpBinWidth;
933  //LOG(info) << " No. of photons after radiator" << nTR0;
934 
935  // compute the spectra behind the mylar foil (absorption)
936  // WinTRspectrum();
937 
938  return 1;
939 }
940 
941 //------------------------------------------------------------------
942 
943 //----- WinTRspectrum -----------------------------------------------
945  //
946  // Computes the spectrum after absorption in Mylar foil
947  //
948 
949  fWinSpectrum->Reset();
950 
951  SetSigma(2);
952 
953  Float_t stemp = 0;
954  for (Int_t iBin = 0; iBin < fSpNBins; iBin++) {
955 
956  Float_t sp = 0;
957  if (fFirstPass == true) {
958  sp = fSpectrum->GetBinContent(iBin + 1);
959  } else {
960  sp = fDetSpectrum->GetBinContent(iBin + 1);
961  }
962  // compute the absorption coefficient
963  Float_t conv = TMath::Exp(-fSigmaWin[iBin]);
964  Float_t wn = sp * conv;
965 
966  fWinSpectrum->SetBinContent(iBin + 1, wn);
967 
968  stemp += wn;
969  }
970 
971  //fnTRprod = stemp * fSpBinWidth;
972  //LOG(info) << "No. of photons after My = "<< fnTRprod;
973 
974  // compute the spectrum absorbed in the D & escaped from the D
975  // DetTRspectrum();
976  return 1;
977 }
978 //----------------------------------------------------------------------------
979 
980 //----- DetTRspectrum ------------------------------------------------
982  //
983  // Computes the spectrum absorbed and passed in the gas mixture
984  //
985 
986  SetSigma(3);
987  Float_t stemp = 0;
988  Float_t stemp2 = 0;
989  for (Int_t iBin = 0; iBin < fSpNBins; iBin++) {
990 
991  //passed spectrum
992  Float_t sp = 0;
993  sp = fWinSpectrum->GetBinContent(iBin + 1);
994  Float_t conv = TMath::Exp(-fSigmaDet[iBin]);
995  Float_t wn = sp * conv;
996 
997  fDetSpectrum->SetBinContent(iBin + 1, wn);
998  stemp += wn;
999 
1000  // absorbed spectrum
1001  Float_t conv2 = 1 - TMath::Exp(-fSigmaDet[iBin]);
1002  Float_t wn2 = sp * conv2;
1003 
1004  fDetSpectrumA->SetBinContent(iBin + 1, wn2);
1005  stemp2 += wn2;
1006  }
1007 
1008  //Float_t nTR1 = stemp * fSpBinWidth;
1009  Float_t nTR2 = stemp2 * fSpBinWidth;
1010 
1011  // Save the number of photons absorbed
1012  fnTRab = nTR2;
1013 
1014  //LOG(info) << " No. of photons escaped: " << nTR1;
1015  //LOG(info) << " No. of photnos absorbed in the gas: " <<nTR2;
1016 
1017  return 1;
1018 }
1019 //----------------------------------------------------------------------------
1020 
1021 //----- Gamma factor--------------------------------------------------
1023 
1024  // put this in the header ?
1025  Float_t massE = 5.1099907e-4; // GeV/c2
1026 
1027  Float_t p2 = fMom.Mag2();
1028  Float_t result = TMath::Sqrt(p2 + massE * massE) / massE;
1029 
1030  return result;
1031 }
1032 //----------------------------------------------------------------------------
1033 
1034 //----- SetSigma---------------------------------------------------------
1035 void CbmTrdRadiator::SetSigma(Int_t SigmaT) {
1036  //
1037  // Sets the absorbtion crosssection for the energies of the TR spectrum
1038  //
1039 
1040 
1041  if (SigmaT == 1) {
1042  if (fSigma) delete[] fSigma;
1043  fSigma = new Float_t[fSpNBins];
1044  for (Int_t iBin = 0; iBin < fSpNBins; iBin++) {
1045  Float_t energykeV = iBin * fSpBinWidth + 1.0;
1046  fSigma[iBin] = Sigma(energykeV);
1047  }
1048  }
1049 
1050  if (SigmaT == 2) {
1051  if (fSigmaWin) delete[] fSigmaWin;
1052  fSigmaWin = new Float_t[fSpNBins];
1053  for (Int_t iBin = 0; iBin < fSpNBins; iBin++) {
1054  Float_t energykeV = iBin * fSpBinWidth + 1.0;
1055  fSigmaWin[iBin] = SigmaWin(energykeV);
1056  }
1057  }
1058 
1059  if (SigmaT == 3) {
1060  if (fSigmaDet) delete[] fSigmaDet;
1061  fSigmaDet = new Float_t[fSpNBins];
1062  for (Int_t iBin = 0; iBin < fSpNBins; iBin++) {
1063  Float_t energykeV = iBin * fSpBinWidth + 1.0;
1064  fSigmaDet[iBin] = SigmaDet(energykeV);
1065  }
1066  }
1067 }
1068 //----------------------------------------------------------------------------
1069 
1070 //----- Sigma-------------------------------------------------------------
1071 Float_t CbmTrdRadiator::Sigma(Float_t energykeV) {
1072  //
1073  // Calculates the absorbtion crosssection for a one-foil-one-gap-radiator
1074  //
1075 
1076  // keV -> MeV
1077  Float_t energyMeV = energykeV * 0.001;
1078  Float_t foil = 0.0;
1079 
1080  if (fFoilMaterial == "polyethylen")
1081  foil = GetMuPo(energyMeV) * fFoilDens * fFoilThickCorr;
1082  else if (fFoilMaterial == "pefoam20")
1083  foil = GetMuPo(energyMeV) * fFoilDens * fFoilThickCorr;
1084  else if (fFoilMaterial == "pefiber")
1085  foil = GetMuPo(energyMeV) * fFoilDens * fFoilThickCorr;
1086  else if (fFoilMaterial == "mylar")
1087  foil = GetMuMy(energyMeV) * fFoilDens * fFoilThickCorr;
1088  else if (fFoilMaterial == "pokalon")
1089  foil = GetMuPok(energyMeV) * fFoilDens * fFoilThickCorr;
1090  else
1091  LOG(error) << "ERROR:: unknown radiator material";
1092 
1093 
1094  if (energyMeV >= 0.001) {
1095  Float_t result = (foil + (GetMuAir(energyMeV) * fGapDens * fGapThickCorr));
1096  return result;
1097  } else {
1098  return 1e6;
1099  }
1100 }
1101 //----------------------------------------------------------------------------
1102 
1103 //----- SigmaWin -------------------------------------------------------
1104 Float_t CbmTrdRadiator::SigmaWin(Float_t energykeV) {
1105  //
1106  // Calculates the absorbtion crosssection for a one-foil
1107  //
1108 
1109  // keV -> MeV
1110  Float_t energyMeV = energykeV * 0.001;
1111  //if (fWindowFoil=="Kapton"){
1112  //if (energyMeV >= 0.001) {
1113  //return((GetMuKa(energyMeV) * fWinDens * fWinThick) + (GetMuAl(energyMeV) * 2.70/*[g/cm^3]*/ * 5E-6/*[cm]*/));
1114  //}
1115  // else {
1116  // return 1e6;
1117  // }
1118  // }
1119  // else {
1120  if (energyMeV >= 0.001) {
1121  return (GetMuMy(energyMeV) * fWinDens * fWinThick);
1122  } else {
1123  return 1e6;
1124  }
1125  // }
1126 }
1127 //----------------------------------------------------------------------------
1128 
1129 //----- SigmaDet --------------------------------------------------------
1130 Float_t CbmTrdRadiator::SigmaDet(Float_t energykeV) {
1131  //
1132  //Calculates the absorbtion crosssection for a choosed gas
1133  //
1134 
1135  // keV -> MeV
1136  Float_t energyMeV = energykeV * 0.001;
1137 
1138  if (energyMeV >= 0.001) {
1139  // densities for CO2 and Xe changed. think density for CO2 is just
1140  // a typo. where the difference for Xe comes from i don't know
1141  // Values are from http://pdg.lbl.gov/AtomicNuclearProperties/
1142  // return(GetMuCO2(energyMeV) * 0.001806 * fGasThick * fCom1 +
1143  // GetMuXe(energyMeV) * 0.00589 * fGasThick * fCom2);
1144  return (GetMuCO2(energyMeV) * 0.00184 * fGasThick * fCom1
1145  + GetMuXe(energyMeV) * 0.00549 * fGasThick * fCom2);
1146  } else {
1147  return 1e6;
1148  }
1149 }
1150 //----------------------------------------------------------------------------
1151 
1152 //----- GetMuPok --------------------------------------------------------
1153 Float_t CbmTrdRadiator::GetMuPok(Float_t energyMeV) {
1154  //
1155  // Returns the photon absorbtion cross section for pokalon N470
1156  //
1157  const Int_t kN = 46;
1158  Float_t en[kN] = {
1159  1.000E-03, 1.500E-03, 2.000E-03, 3.000E-03, 4.000E-03, 5.000E-03, 6.000E-03,
1160  8.000E-03, 1.000E-02, 1.500E-02, 2.000E-02, 3.000E-02, 4.000E-02, 5.000E-02,
1161  6.000E-02, 8.000E-02, 1.000E-01, 1.500E-01, 2.000E-01, 3.000E-01, 4.000E-01,
1162  5.000E-01, 6.000E-01, 8.000E-01, 1.000E+00, 1.022E+00, 1.250E+00, 1.500E+00,
1163  2.000E+00, 2.044E+00, 3.000E+00, 4.000E+00, 5.000E+00, 6.000E+00, 7.000E+00,
1164  8.000E+00, 9.000E+00, 1.000E+01, 1.100E+01, 1.200E+01, 1.300E+01, 1.400E+01,
1165  1.500E+01, 1.600E+01, 1.800E+01, 2.000E+01};
1166  Float_t mu[kN] = {
1167  2.537E+03, 8.218E+02, 3.599E+02, 1.093E+02, 4.616E+01, 2.351E+01, 1.352E+01,
1168  5.675E+00, 2.938E+00, 9.776E-01, 5.179E-01, 2.847E-01, 2.249E-01, 2.003E-01,
1169  1.866E-01, 1.705E-01, 1.600E-01, 1.422E-01, 1.297E-01, 1.125E-01, 1.007E-01,
1170  9.193E-02, 8.501E-02, 7.465E-02, 6.711E-02, 6.642E-02, 6.002E-02, 5.463E-02,
1171  4.686E-02, 4.631E-02, 3.755E-02, 3.210E-02, 2.851E-02, 2.597E-02, 2.409E-02,
1172  2.263E-02, 2.149E-02, 2.056E-02, 1.979E-02, 1.916E-02, 1.862E-02, 1.816E-02,
1173  1.777E-02, 1.743E-02, 1.687E-02, 1.644E-02};
1174  return Interpolate(energyMeV, en, mu, kN);
1175 }
1176 //----------------------------------------------------------------------------
1177 //----- GetMuKa --------------------------------------------------------
1178 Float_t CbmTrdRadiator::GetMuKa(Float_t energyMeV) {
1179  //
1180  // Returns the photon absorbtion cross section for kapton
1181  //
1182  const Int_t kN = 46;
1183  Float_t en[kN] = {
1184  1.000E-03, 1.500E-03, 2.000E-03, 3.000E-03, 4.000E-03, 5.000E-03, 6.000E-03,
1185  8.000E-03, 1.000E-02, 1.500E-02, 2.000E-02, 3.000E-02, 4.000E-02, 5.000E-02,
1186  6.000E-02, 8.000E-02, 1.000E-01, 1.500E-01, 2.000E-01, 3.000E-01, 4.000E-01,
1187  5.000E-01, 6.000E-01, 8.000E-01, 1.000E+00, 1.022E+00, 1.250E+00, 1.500E+00,
1188  2.000E+00, 2.044E+00, 3.000E+00, 4.000E+00, 5.000E+00, 6.000E+00, 7.000E+00,
1189  8.000E+00, 9.000E+00, 1.000E+01, 1.100E+01, 1.200E+01, 1.300E+01, 1.400E+01,
1190  1.500E+01, 1.600E+01, 1.800E+01, 2.000E+01};
1191  Float_t mu[kN] = {
1192  2.731E+03, 8.875E+02, 3.895E+02, 1.185E+02, 5.013E+01, 2.555E+01, 1.470E+01,
1193  6.160E+00, 3.180E+00, 1.043E+00, 5.415E-01, 2.880E-01, 2.235E-01, 1.973E-01,
1194  1.830E-01, 1.666E-01, 1.560E-01, 1.385E-01, 1.263E-01, 1.095E-01, 9.799E-02,
1195  8.944E-02, 8.270E-02, 7.262E-02, 6.529E-02, 6.462E-02, 5.839E-02, 5.315E-02,
1196  4.561E-02, 4.507E-02, 3.659E-02, 3.133E-02, 2.787E-02, 2.543E-02, 2.362E-02,
1197  2.223E-02, 2.114E-02, 2.025E-02, 1.953E-02, 1.893E-02, 1.842E-02, 1.799E-02,
1198  1.762E-02, 1.730E-02, 1.678E-02, 1.638E-02};
1199  return Interpolate(energyMeV, en, mu, kN);
1200 }
1201 //----------------------------------------------------------------------------
1202 
1203 //----- GetMuAl --------------------------------------------------------
1204 Float_t CbmTrdRadiator::GetMuAl(Float_t energyMeV) {
1205  //
1206  // Returns the photon absorbtion cross section for Al
1207  //
1208 
1209  const Int_t kN = 48;
1210 
1211  Float_t en[kN] = {
1212  1.000E-03, 1.500E-03, 1.560E-03, 1.560E-03, 2.000E-03, 3.000E-03, 4.000E-03,
1213  5.000E-03, 6.000E-03, 8.000E-03, 1.000E-02, 1.500E-02, 2.000E-02, 3.000E-02,
1214  4.000E-02, 5.000E-02, 6.000E-02, 8.000E-02, 1.000E-01, 1.500E-01, 2.000E-01,
1215  3.000E-01, 4.000E-01, 5.000E-01, 6.000E-01, 8.000E-01, 1.000E+00, 1.022E+00,
1216  1.250E+00, 1.500E+00, 2.000E+00, 2.044E+00, 3.000E+00, 4.000E+00, 5.000E+00,
1217  6.000E+00, 7.000E+00, 8.000E+00, 9.000E+00, 1.000E+01, 1.100E+01, 1.200E+01,
1218  1.300E+01, 1.400E+01, 1.500E+01, 1.600E+01, 1.800E+01, 2.000E+01};
1219 
1220  Float_t mu[kN] = {
1221  1.185E+03, 4.023E+02, 3.621E+02, 3.957E+03, 2.263E+03, 7.881E+02, 3.605E+02,
1222  1.934E+02, 1.153E+02, 5.032E+01, 2.621E+01, 7.955E+00, 3.442E+00, 1.128E+00,
1223  5.684E-01, 3.681E-01, 2.778E-01, 2.018E-01, 1.704E-01, 1.378E-01, 1.223E-01,
1224  1.042E-01, 9.276E-02, 8.445E-02, 7.802E-02, 6.841E-02, 6.146E-02, 6.080E-02,
1225  5.496E-02, 5.006E-02, 4.324E-02, 4.277E-02, 3.541E-02, 3.106E-02, 2.836E-02,
1226  2.655E-02, 2.529E-02, 2.437E-02, 2.369E-02, 2.318E-02, 2.279E-02, 2.249E-02,
1227  2.226E-02, 2.208E-02, 2.195E-02, 2.185E-02, 2.173E-02, 2.168E-02};
1228 
1229  return Interpolate(energyMeV, en, mu, kN);
1230 }
1231 //----------------------------------------------------------------------------
1232 
1233 
1234 //----- GetMuPo --------------------------------------------------------
1235 Float_t CbmTrdRadiator::GetMuPo(Float_t energyMeV) {
1236  //
1237  // Returns the photon absorbtion cross section for polypropylene
1238  //
1239 
1240  const Int_t kN = 36;
1241 
1242  Float_t mu[kN] = {
1243  1.894E+03, 5.999E+02, 2.593E+02, 7.743E+01, 3.242E+01, 1.643E+01,
1244  9.432E+00, 3.975E+00, 2.088E+00, 7.452E-01, 4.315E-01, 2.706E-01,
1245  2.275E-01, 2.084E-01, 1.970E-01, 1.823E-01, 1.719E-01, 1.534E-01,
1246  1.402E-01, 1.217E-01, 1.089E-01, 9.947E-02, 9.198E-02, 8.078E-02,
1247  7.262E-02, 6.495E-02, 5.910E-02, 5.064E-02, 4.045E-02, 3.444E-02,
1248  3.045E-02, 2.760E-02, 2.383E-02, 2.145E-02, 1.819E-02, 1.658E-02};
1249 
1250  Float_t en[kN] = {
1251  1.000E-03, 1.500E-03, 2.000E-03, 3.000E-03, 4.000E-03, 5.000E-03,
1252  6.000E-03, 8.000E-03, 1.000E-02, 1.500E-02, 2.000E-02, 3.000E-02,
1253  4.000E-02, 5.000E-02, 6.000E-02, 8.000E-02, 1.000E-01, 1.500E-01,
1254  2.000E-01, 3.000E-01, 4.000E-01, 5.000E-01, 6.000E-01, 8.000E-01,
1255  1.000E+00, 1.250E+00, 1.500E+00, 2.000E+00, 3.000E+00, 4.000E+00,
1256  5.000E+00, 6.000E+00, 8.000E+00, 1.000E+01, 1.500E+01, 2.000E+01};
1257 
1258  return Interpolate(energyMeV, en, mu, kN);
1259 }
1260 //----------------------------------------------------------------------------
1261 
1262 //----- GetMuAir --------------------------------------------------------
1263 Float_t CbmTrdRadiator::GetMuAir(Float_t energyMeV) {
1264  //
1265  // Returns the photon absorbtion cross section for Air
1266  //
1267 
1268  const Int_t kN = 38;
1269 
1270 
1271  Float_t mu[kN] = {
1272  0.35854E+04, 0.11841E+04, 0.52458E+03, 0.16143E+03, 0.15722E+03,
1273  0.14250E+03, 0.77538E+02, 0.40099E+02, 0.23313E+02, 0.98816E+01,
1274  0.51000E+01, 0.16079E+01, 0.77536E+00, 0.35282E+00, 0.24790E+00,
1275  0.20750E+00, 0.18703E+00, 0.16589E+00, 0.15375E+00, 0.13530E+00,
1276  0.12311E+00, 0.10654E+00, 0.95297E-01, 0.86939E-01, 0.80390E-01,
1277  0.70596E-01, 0.63452E-01, 0.56754E-01, 0.51644E-01, 0.44382E-01,
1278  0.35733E-01, 0.30721E-01, 0.27450E-01, 0.25171E-01, 0.22205E-01,
1279  0.20399E-01, 0.18053E-01, 0.18057E-01};
1280 
1281  Float_t en[kN] = {
1282  0.10000E-02, 0.15000E-02, 0.20000E-02, 0.30000E-02, 0.32029E-02,
1283  0.32029E-02, 0.40000E-02, 0.50000E-02, 0.60000E-02, 0.80000E-02,
1284  0.10000E-01, 0.15000E-01, 0.20000E-01, 0.30000E-01, 0.40000E-01,
1285  0.50000E-01, 0.60000E-01, 0.80000E-01, 0.10000E+00, 0.15000E+00,
1286  0.20000E+00, 0.30000E+00, 0.40000E+00, 0.50000E+00, 0.60000E+00,
1287  0.80000E+00, 0.10000E+01, 0.12500E+01, 0.15000E+01, 0.20000E+01,
1288  0.30000E+01, 0.40000E+01, 0.50000E+01, 0.60000E+01, 0.80000E+01,
1289  0.10000E+02, 0.15000E+02, 0.20000E+02};
1290 
1291  return Interpolate(energyMeV, en, mu, kN);
1292 }
1293 //----------------------------------------------------------------------------
1294 
1295 //----- GetMuXe --------------------------------------------------------
1296 Float_t CbmTrdRadiator::GetMuXe(Float_t energyMeV) {
1297  //
1298  // Returns the photon absorbtion cross section for xenon
1299  //
1300 
1301  const Int_t kN = 48;
1302 
1303  Float_t mu[kN] = {
1304  9.413E+03, 8.151E+03, 7.035E+03, 7.338E+03, 4.085E+03, 2.088E+03, 7.780E+02,
1305  3.787E+02, 2.408E+02, 6.941E+02, 6.392E+02, 6.044E+02, 8.181E+02, 7.579E+02,
1306  6.991E+02, 8.064E+02, 6.376E+02, 3.032E+02, 1.690E+02, 5.743E+01, 2.652E+01,
1307  8.930E+00, 6.129E+00, 3.316E+01, 2.270E+01, 1.272E+01, 7.825E+00, 3.633E+00,
1308  2.011E+00, 7.202E-01, 3.760E-01, 1.797E-01, 1.223E-01, 9.699E-02, 8.281E-02,
1309  6.696E-02, 5.785E-02, 5.054E-02, 4.594E-02, 4.078E-02, 3.681E-02, 3.577E-02,
1310  3.583E-02, 3.634E-02, 3.797E-02, 3.987E-02, 4.445E-02, 4.815E-02};
1311 
1312  Float_t en[kN] = {
1313  1.00000E-03, 1.07191E-03, 1.14900E-03, 1.14900E-03, 1.50000E-03,
1314  2.00000E-03, 3.00000E-03, 4.00000E-03, 4.78220E-03, 4.78220E-03,
1315  5.00000E-03, 5.10370E-03, 5.10370E-03, 5.27536E-03, 5.45280E-03,
1316  5.45280E-03, 6.00000E-03, 8.00000E-03, 1.00000E-02, 1.50000E-02,
1317  2.00000E-02, 3.00000E-02, 3.45614E-02, 3.45614E-02, 4.00000E-02,
1318  5.00000E-02, 6.00000E-02, 8.00000E-02, 1.00000E-01, 1.50000E-01,
1319  2.00000E-01, 3.00000E-01, 4.00000E-01, 5.00000E-01, 6.00000E-01,
1320  8.00000E-01, 1.00000E+00, 1.25000E+00, 1.50000E+00, 2.00000E+00,
1321  3.00000E+00, 4.00000E+00, 5.00000E+00, 6.00000E+00, 8.00000E+00,
1322  1.00000E+01, 1.50000E+01, 2.00000E+01};
1323 
1324  return Interpolate(energyMeV, en, mu, kN);
1325 }
1326 //----------------------------------------------------------------------------
1327 
1328 //----- GetMuCO2 ------------------------------------------------------
1329 Float_t CbmTrdRadiator::GetMuCO2(Float_t energyMeV) {
1330  //
1331  // Returns the photon absorbtion cross section for CO2
1332  //
1333 
1334  const Int_t kN = 36;
1335 
1336  Float_t mu[kN] = {0.39383E+04, 0.13166E+04, 0.58750E+03, 0.18240E+03,
1337  0.77996E+02, 0.40024E+02, 0.23116E+02, 0.96997E+01,
1338  0.49726E+01, 0.15543E+01, 0.74915E+00, 0.34442E+00,
1339  0.24440E+00, 0.20589E+00, 0.18632E+00, 0.16578E+00,
1340  0.15394E+00, 0.13558E+00, 0.12336E+00, 0.10678E+00,
1341  0.95510E-01, 0.87165E-01, 0.80587E-01, 0.70769E-01,
1342  0.63626E-01, 0.56894E-01, 0.51782E-01, 0.44499E-01,
1343  0.35839E-01, 0.30825E-01, 0.27555E-01, 0.25269E-01,
1344  0.22311E-01, 0.20516E-01, 0.18184E-01, 0.17152E-01};
1345 
1346  Float_t en[kN] = {0.10000E-02, 0.15000E-02, 0.20000E-02, 0.30000E-02,
1347  0.40000E-02, 0.50000E-02, 0.60000E-02, 0.80000E-02,
1348  0.10000E-01, 0.15000E-01, 0.20000E-01, 0.30000E-01,
1349  0.40000E-01, 0.50000E-01, 0.60000E-01, 0.80000E-01,
1350  0.10000E+00, 0.15000E+00, 0.20000E+00, 0.30000E+00,
1351  0.40000E+00, 0.50000E+00, 0.60000E+00, 0.80000E+00,
1352  0.10000E+01, 0.12500E+01, 0.15000E+01, 0.20000E+01,
1353  0.30000E+01, 0.40000E+01, 0.50000E+01, 0.60000E+01,
1354  0.80000E+01, 0.10000E+02, 0.15000E+02, 0.20000E+02};
1355 
1356  return Interpolate(energyMeV, en, mu, kN);
1357 }
1358 //----------------------------------------------------------------------------
1359 
1360 //----- GetMuMy -------------------------------------------------------
1361 Float_t CbmTrdRadiator::GetMuMy(Float_t energyMeV) {
1362  //
1363  // Returns the photon absorbtion cross section for mylar
1364  //
1365 
1366  const Int_t kN = 36;
1367 
1368  Float_t mu[kN] = {
1369  2.911E+03, 9.536E+02, 4.206E+02, 1.288E+02, 5.466E+01, 2.792E+01,
1370  1.608E+01, 6.750E+00, 3.481E+00, 1.132E+00, 5.798E-01, 3.009E-01,
1371  2.304E-01, 2.020E-01, 1.868E-01, 1.695E-01, 1.586E-01, 1.406E-01,
1372  1.282E-01, 1.111E-01, 9.947E-02, 9.079E-02, 8.395E-02, 7.372E-02,
1373  6.628E-02, 5.927E-02, 5.395E-02, 4.630E-02, 3.715E-02, 3.181E-02,
1374  2.829E-02, 2.582E-02, 2.257E-02, 2.057E-02, 1.789E-02, 1.664E-02};
1375 
1376  Float_t en[kN] = {1.00000E-03, 1.50000E-03, 2.00000E-03, 3.00000E-03,
1377  4.00000E-03, 5.00000E-03, 6.00000E-03, 8.00000E-03,
1378  1.00000E-02, 1.50000E-02, 2.00000E-02, 3.00000E-02,
1379  4.00000E-02, 5.00000E-02, 6.00000E-02, 8.00000E-02,
1380  1.00000E-01, 1.50000E-01, 2.00000E-01, 3.00000E-01,
1381  4.00000E-01, 5.00000E-01, 6.00000E-01, 8.00000E-01,
1382  1.00000E+00, 1.25000E+00, 1.50000E+00, 2.00000E+00,
1383  3.00000E+00, 4.00000E+00, 5.00000E+00, 6.00000E+00,
1384  8.00000E+00, 1.00000E+01, 1.50000E+01, 2.00000E+01};
1385 
1386  return Interpolate(energyMeV, en, mu, kN);
1387 }
1388 //----------------------------------------------------------------------------
1389 
1390 //----- Interpolate ------------------------------------------------------
1391 Float_t CbmTrdRadiator::Interpolate(Float_t energyMeV,
1392  Float_t* en,
1393  Float_t* mu,
1394  Int_t n) {
1395  //
1396  // Interpolates the photon absorbtion cross section
1397  // for a given energy <energyMeV>.
1398  //
1399 
1400  Float_t de = 0;
1401  Int_t index = 0;
1402  Int_t istat = Locate(en, n, energyMeV, index, de);
1403  if (istat == 0) {
1404  Float_t result =
1405  (mu[index]
1406  - de * (mu[index] - mu[index + 1]) / (en[index + 1] - en[index]));
1407  return result;
1408  } else {
1409  return 0.0;
1410  }
1411 }
1412 //----------------------------------------------------------------------------
1413 
1414 //----- Locate ------------------------------------------------------------
1415 Int_t CbmTrdRadiator::Locate(Float_t* xv,
1416  Int_t n,
1417  Float_t xval,
1418  Int_t& kl,
1419  Float_t& dx) {
1420  //
1421  // Locates a point (xval) in a 1-dim grid (xv(n))
1422  //
1423 
1424  if (xval >= xv[n - 1]) return 1;
1425  if (xval < xv[0]) return -1;
1426 
1427  Int_t km;
1428  Int_t kh = n - 1;
1429 
1430  kl = 0;
1431  while (kh - kl > 1) {
1432  if (xval < xv[km = (kl + kh) / 2])
1433  kh = km;
1434  else
1435  kl = km;
1436  }
1437  if (xval < xv[kl] || xval > xv[kl + 1] || kl >= n - 1) {
1438  printf("Locate failed xv[%d] %f xval %f xv[%d] %f!!!\n",
1439  kl,
1440  xv[kl],
1441  xval,
1442  kl + 1,
1443  xv[kl + 1]);
1444  exit(1);
1445  }
1446 
1447  dx = xval - xv[kl];
1448  if (xval == 0.001) LOG(info) << "Locat = 0";
1449  return 0;
1450 }
1451 //----------------------------------------------------------------------------
1452 
1453 
CbmTrdGas::GetDetType
Int_t GetDetType() const
Definition: CbmTrdGas.h:22
CbmTrdPoint::GetYOut
Double_t GetYOut() const
Definition: CbmTrdPoint.h:67
CbmTrdRadiator::DetTRspectrum
Int_t DetTRspectrum()
Definition: CbmTrdRadiator.cxx:981
CbmTrdPoint::GetZIn
Double_t GetZIn() const
Definition: CbmTrdPoint.h:65
CbmTrdRadiator::fFoilThick
Float_t fFoilThick
Definition: CbmTrdRadiator.h:152
CbmTrdRadiator::GetMuCO2
Float_t GetMuCO2(Float_t energyMeV)
Definition: CbmTrdRadiator.cxx:1329
CbmTrdRadiator::fGasThick
Float_t fGasThick
Definition: CbmTrdRadiator.h:155
CbmTrdGas.h
Container for gas properties of TRD.
CbmTrdRadiator::fSpRange
static const Int_t fSpRange
Definition: CbmTrdRadiator.h:183
CbmTrdRadiator::CbmTrdRadiator
CbmTrdRadiator()
Definition: CbmTrdRadiator.cxx:29
CbmTrdRadiator::fDetSpectrum
TH1D * fDetSpectrum
TR absorbed in Detector.
Definition: CbmTrdRadiator.h:195
CbmTrdRadiator::fSimpleTR
Bool_t fSimpleTR
Definition: CbmTrdRadiator.h:150
i
int i
Definition: L1/vectors/P4_F32vec4.h:25
CbmTrdRadiator::~CbmTrdRadiator
virtual ~CbmTrdRadiator()
Definition: CbmTrdRadiator.cxx:143
CbmTrdRadiator::SigmaDet
Float_t SigmaDet(Float_t energykeV)
Definition: CbmTrdRadiator.cxx:1130
CbmTrdRadiator::fDetType
Int_t fDetType
Definition: CbmTrdRadiator.h:148
CbmTrdRadiator::ProduceSpectra
void ProduceSpectra()
Definition: CbmTrdRadiator.cxx:698
CbmTrdPoint::GetXOut
Double_t GetXOut() const
Definition: CbmTrdPoint.h:66
CbmTrdRadiator::GammaF
Float_t GammaF()
Definition: CbmTrdRadiator.cxx:1022
CbmTrdRadiator::fSigmaWin
Float_t * fSigmaWin
[fSpNBins] Array of sigma values for the foil of the radiator
Definition: CbmTrdRadiator.h:189
CbmTrdRadiator::Interpolate
Float_t Interpolate(Float_t energyMeV, Float_t *en, Float_t *mu, Int_t n)
Definition: CbmTrdRadiator.cxx:1391
CbmTrdRadiator::fGapDens
Float_t fGapDens
Definition: CbmTrdRadiator.h:161
CbmTrdRadiator::Sigma
Float_t Sigma(Float_t energykeV)
Definition: CbmTrdRadiator.cxx:1071
CbmTrdRadiator::fWinThick
Float_t fWinThick
Definition: CbmTrdRadiator.h:177
CbmTrdRadiator.h
CbmTrdGas::GetNobleGas
Double_t GetNobleGas() const
Definition: CbmTrdGas.h:24
CbmTrdRadiator::GetMuMy
Float_t GetMuMy(Float_t energyMeV)
Definition: CbmTrdRadiator.cxx:1361
CbmTrdRadiator::fFoilMaterial
TString fFoilMaterial
Definition: CbmTrdRadiator.h:154
CbmTrdRadiator::ELoss
Int_t ELoss(Int_t index)
Definition: CbmTrdRadiator.cxx:813
CbmTrdRadiator::LatticeHit
Bool_t LatticeHit(const CbmTrdPoint *point)
Definition: CbmTrdRadiator.cxx:618
CbmTrdRadiator::fGapThick
Float_t fGapThick
Definition: CbmTrdRadiator.h:153
CbmTrdRadiator::fSpBinWidth
Float_t fSpBinWidth
Definition: CbmTrdRadiator.h:184
CbmTrdRadiator::fNFoils
Int_t fNFoils
Definition: CbmTrdRadiator.h:151
CbmTrdRadiator::fWinSpectrum
TH1D * fWinSpectrum
TR photon energy spectrum.
Definition: CbmTrdRadiator.h:193
CbmTrdRadiator::fnPhotonCorr
Float_t fnPhotonCorr
Definition: CbmTrdRadiator.h:164
CbmTrdRadiator::fFoilOmega
Float_t fFoilOmega
Definition: CbmTrdRadiator.h:162
CbmTrdGas::GetCO2
Double_t GetCO2() const
Definition: CbmTrdGas.h:25
CbmTrdRadiator::fnTRabs
Float_t fnTRabs[fNMom]
Absorption spectra for different momenta.
Definition: CbmTrdRadiator.h:202
CbmTrdRadiator::SetSigma
void SetSigma(Int_t SigmaT)
Definition: CbmTrdRadiator.cxx:1035
CbmTrdRadiator::fTrackMomentum
Double_t * fTrackMomentum
Definition: CbmTrdRadiator.h:198
CbmTrdRadiator::fGapThickCorr
Float_t fGapThickCorr
Definition: CbmTrdRadiator.h:169
CbmTrdRadiator::fSigmaDet
Float_t * fSigmaDet
[fSpNBins] Array of sigma values for the entrance window of detector
Definition: CbmTrdRadiator.h:190
CbmTrdRadiator::fGapOmega
Float_t fGapOmega
Definition: CbmTrdRadiator.h:163
CbmTrdRadiator::fFirstPass
Bool_t fFirstPass
Definition: CbmTrdRadiator.h:149
CbmTrdRadiator::GetMuXe
Float_t GetMuXe(Float_t energyMeV)
Definition: CbmTrdRadiator.cxx:1296
CbmTrdGas::GetGasThick
Double_t GetGasThick() const
Definition: CbmTrdGas.h:23
CbmTrdRadiator::Locate
Int_t Locate(Float_t *xv, Int_t n, Float_t xval, Int_t &kl, Float_t &dx)
Definition: CbmTrdRadiator.cxx:1415
CbmTrdRadiator::fDetSpectrumA
TH1D * fDetSpectrumA
TR spectra in gas-window foil.
Definition: CbmTrdRadiator.h:194
ClassImp
ClassImp(CbmConverterManager) InitStatus CbmConverterManager
Definition: CbmConverterManager.cxx:12
CbmTrdRadiator::GetMuAl
Float_t GetMuAl(Float_t energyMeV)
Definition: CbmTrdRadiator.cxx:1204
CbmTrdRadiator::fNMom
static const Int_t fNMom
TR passed through Detector.
Definition: CbmTrdRadiator.h:197
CbmTrdRadiator::fFoilThickCorr
Float_t fFoilThickCorr
Definition: CbmTrdRadiator.h:168
xMath::Pi
double Pi()
Definition: xMath.h:5
CbmTrdPoint::GetZOut
Double_t GetZOut() const
Definition: CbmTrdPoint.h:68
CbmTrdRadiator::ProcessTR
void ProcessTR()
Definition: CbmTrdRadiator.cxx:771
CbmTrdRadiator
Definition: CbmTrdRadiator.h:24
CbmTrdRadiator::SigmaWin
Float_t SigmaWin(Float_t energykeV)
Definition: CbmTrdRadiator.cxx:1104
CbmTrdPoint
Definition: CbmTrdPoint.h:23
CbmTrdRadiator::GetMuPo
Float_t GetMuPo(Float_t energyMeV)
Definition: CbmTrdRadiator.cxx:1235
CbmTrdRadiator::Init
void Init()
Definition: CbmTrdRadiator.cxx:405
CbmTrdRadiator::TRspectrum
Int_t TRspectrum()
Definition: CbmTrdRadiator.cxx:862
CbmTrdRadiator::fnTRab
Float_t fnTRab
Definition: CbmTrdRadiator.h:203
CbmTrdRadiator::fFinal
TH1D * fFinal[fNMom]
[fNMom] Track momenta for which spectra
Definition: CbmTrdRadiator.h:201
CbmTrdRadiator::fWinDens
Float_t fWinDens
Definition: CbmTrdRadiator.h:176
CbmTrdPoint::GetXIn
Double_t GetXIn() const
Definition: CbmTrdPoint.h:63
CbmTrdGas::Instance
static CbmTrdGas * Instance()
Definition: CbmTrdGas.h:29
CbmTrdRadiator::GetTR
Float_t GetTR(TVector3 mom)
Definition: CbmTrdRadiator.cxx:732
CbmTrdRadiator::fCom2
Float_t fCom2
Definition: CbmTrdRadiator.h:180
CbmTrdRadiator::WinTRspectrum
Int_t WinTRspectrum()
Definition: CbmTrdRadiator.cxx:944
CbmTrdRadiator::fMom
TVector3 fMom
Definition: CbmTrdRadiator.h:207
m
__m128 m
Definition: L1/vectors/P4_F32vec4.h:26
CbmTrdPoint.h
CbmTrdRadiator::GetMuKa
Float_t GetMuKa(Float_t energyMeV)
Definition: CbmTrdRadiator.cxx:1178
pos
TVector3 pos
Definition: CbmMvdSensorDigiToHitTask.cxx:60
CbmTrdGas
Definition: CbmTrdGas.h:14
CbmTrdRadiator::fGasThickCorr
Float_t fGasThickCorr
Definition: CbmTrdRadiator.h:170
CbmTrdRadiator::fSigma
Float_t * fSigma
Definition: CbmTrdRadiator.h:187
CbmTrdRadiator::GetMuPok
Float_t GetMuPok(Float_t energyMeV)
Definition: CbmTrdRadiator.cxx:1153
CbmTrdGas::Init
void Init()
Definition: CbmTrdGas.cxx:40
CbmTrdRadiator::fCom1
Float_t fCom1
Definition: CbmTrdRadiator.h:179
CbmTrdRadiator::fSpNBins
static const Int_t fSpNBins
Definition: CbmTrdRadiator.h:182
CbmTrdRadiator::fSpectrum
TH1D * fSpectrum
[fSpNBins] Array of sigma values for the active gas
Definition: CbmTrdRadiator.h:192
CbmTrdRadiator::CreateHistograms
void CreateHistograms()
Definition: CbmTrdRadiator.cxx:160
CbmTrdRadiator::fRadType
TString fRadType
Definition: CbmTrdRadiator.h:147
CbmTrdPoint::GetYIn
Double_t GetYIn() const
Definition: CbmTrdPoint.h:64
CbmTrdRadiator::fELoss
Float_t fELoss
Definition: CbmTrdRadiator.h:205
CbmTrdRadiator::GetMuAir
Float_t GetMuAir(Float_t energyMeV)
Definition: CbmTrdRadiator.cxx:1263
CbmTrdRadiator::fFoilDens
Float_t fFoilDens
Definition: CbmTrdRadiator.h:160