CbmRoot
Loading...
Searching...
No Matches
CbmTrdRadiator.cxx
Go to the documentation of this file.
1/* Copyright (C) 2004-2021 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Matus Kalisky [committer], Cyrano Bergmann, Adrian Meyer-Ahrens, David Emschermann, Florian Uhlig */
4
5// ------------------------------------------------------------------------
6// ----- CbmTrdRadiator source file -----
7// ----- Created 10/11/04 by M.Kalisky -----
8// ------------------------------------------------------------------------
9#include "CbmTrdRadiator.h"
10
11#include "CbmTrdGas.h" // for CbmTrdGas
12#include "CbmTrdPoint.h" // for CbmTrdPoint
13
14#include <Logger.h> // for Logger, LOG
15
16#include <TDirectory.h> // for TDirectory
17#include <TFile.h> // for TFile, gFile
18#include <TGeoManager.h> // for TGeoManager, gGeoManager
19#include <TH1.h> // for TH1D
20#include <TMath.h> // for Exp, Sqrt, Cos, Pi
21#include <TMathBase.h> // for Abs
22#include <TObjString.h> // for TObjString
23#include <TRandom.h> // for TRandom
24#include <TRandom3.h> // for gRandom
25#include <TVector3.h> // for TVector3
26
27#include <cstdio> // for printf, sprintf
28#include <cstdlib> // for exit
29#include <iomanip> // for setprecision, __iom_t5
30
31using std::setprecision;
32
33// ----- Default constructor ------------------------------------------------
34CbmTrdRadiator::CbmTrdRadiator(Bool_t SimpleTR, Int_t Nfoils, Float_t FoilThick, Float_t GapThick, TString material,
35 TString window) // If "Kapton" is set, 0.05 mu aluminum will be added
36 : fWindowFoil(window)
37 , fRadType("")
38 , fDetType(-1)
39 , fFirstPass(kTRUE)
40 , fSimpleTR(SimpleTR)
41 , fNFoils(Nfoils)
42 , fFoilThick(FoilThick)
43 , fGapThick(GapThick)
44 , fFoilMaterial(material)
45 , fGasThick(-1)
46 , fFoilDens(-1.)
47 , fGapDens(-1.)
48 , fFoilOmega(-1.)
49 , fGapOmega(-1.)
50 , fnPhotonCorr(0.65)
51 , fFoilThickCorr(1.)
52 , fGapThickCorr(1.)
53 , fGasThickCorr(1.)
54 , fnTRprod(-1)
55 , fWinDens(-1.)
56 , fWinThick(-1.)
57 , WINDOW()
58 , fCom1(-1.)
59 , fCom2(-1.)
61 , fSigma(nullptr)
62 , fSigmaWin(nullptr)
63 , fSigmaDet(nullptr)
64 , fSpectrum(nullptr)
65 , fWinSpectrum(nullptr)
66 , fDetSpectrumA(nullptr)
67 , fDetSpectrum(nullptr)
68 , fTrackMomentum(nullptr)
69 , fFinal()
70 , fnTRabs()
71 , fnTRab(-1.)
72 , fELoss(-1.)
73 , fMom(-1., -1., -1.)
74{
75 for (Int_t i = 0; i < fNMom; i++) {
76 fFinal[i] = nullptr;
77 }
78}
79//-----------------------------------------------------------------------------
80
81// ----- Constructor using predefined radiator prototype------------------------------------
82CbmTrdRadiator::CbmTrdRadiator(Bool_t SimpleTR, TString prototype,
83 TString window) // If "Kapton" is set, 0.05 mu aluminum will be added
84 : fWindowFoil(window)
85 , fRadType(prototype)
86 , fDetType(-1)
87 , fFirstPass(kTRUE)
88 , fSimpleTR(SimpleTR)
89 , fNFoils(337)
90 , fFoilThick(0.0012)
91 , fGapThick(0.09)
92 , fFoilMaterial("")
93 , fGasThick(-1)
94 , fFoilDens(-1.)
95 , fGapDens(-1.)
96 , fFoilOmega(-1.)
97 , fGapOmega(-1.)
98 , fnPhotonCorr(0.65)
99 , fFoilThickCorr(1.)
100 , fGapThickCorr(1.)
101 , fGasThickCorr(1.)
102 , fnTRprod(-1)
103 , fWinDens(-1.)
104 , fWinThick(-1.)
105 , WINDOW()
106 , fCom1(-1.)
107 , fCom2(-1.)
109 , fSigma(nullptr)
110 , fSigmaWin(nullptr)
111 , fSigmaDet(nullptr)
112 , fSpectrum(nullptr)
113 , fWinSpectrum(nullptr)
114 , fDetSpectrumA(nullptr)
115 , fDetSpectrum(nullptr)
116 , fTrackMomentum(nullptr)
117 , fFinal()
118 , fnTRabs()
119 , fnTRab(-1.)
120 , fELoss(-1.)
121 , fMom(-1., -1., -1.)
122{
123 for (Int_t i = 0; i < fNMom; i++) {
124 fFinal[i] = nullptr;
125 }
126}
127//-----------------------------------------------------------------------------
128
129// ----- Destructor ----------------------------------------------------
131{
132 // FairRootManager *fManager = FairRootManager::Instance();
133 // fManager->Write();
134 if (fSpectrum) {
135 LOG(info) << " -I DELETING fSpectrum ";
136 delete fSpectrum;
137 }
138 if (fWinSpectrum) delete fWinSpectrum;
139 if (fDetSpectrum) delete fDetSpectrum;
140 if (fDetSpectrumA) delete fDetSpectrumA;
141
142 for (Int_t i = 0; i < fNMom; i++) {
143 if (fFinal[i]) delete fFinal[i];
144 }
145}
146
147// ----- Create Histogramms --------------------------------------------------
149{
150
151 // Create the needed histograms
152
154
155 Float_t SpLower = 1.0 - 0.5 * fSpBinWidth;
156 Float_t SpUpper = SpLower + (Float_t) fSpRange;
157
158 size_t buf_size = 50;
159 Char_t name[buf_size];
160
161 if (fSpectrum) delete fSpectrum;
162 fSpectrum = new TH1D("fSpectrum", "TR spectrum", fSpNBins, SpLower, SpUpper);
163
164 if (fWinSpectrum) delete fWinSpectrum;
165 fWinSpectrum = new TH1D("fWinSpectrum", "TR spectrum in Mylar", fSpNBins, SpLower, SpUpper);
166
167 if (fDetSpectrum) delete fDetSpectrum;
168 fDetSpectrum = new TH1D("fDetSpectrum", "TR spectrum escaped from detector", fSpNBins, SpLower, SpUpper);
169
170 if (fDetSpectrumA) delete fDetSpectrumA;
171 fDetSpectrumA = new TH1D("fDetSpectrumA", "TR spectrum absorbed in detector", fSpNBins, SpLower, SpUpper);
172
173 for (Int_t i = 0; i < fNMom; i++) {
174 snprintf(name, buf_size - 1, "fFinal%d", i + 1);
175 //LOG(info) <<"name : "<<name;
176 if (fFinal[i]) delete fFinal[i];
177 fFinal[i] = new TH1D(name, name, fSpNBins, SpLower, SpUpper);
178 }
179}
180
181
182// ----- Init function ----------------------------------------------------
184{
185 LOG(info) << "================CbmTrdRadiator===============";
187
188 //Check if a radiator type has been set
189 if (fRadType == "") {
190 //Set fnPhotonCorr also if no prototype has been specified
191 fnPhotonCorr = 0.65;
192 //Check specified radiator types
193 }
194 else if (fRadType == "tdr18") { //Default TDR18 radiator 30cm Box + 1.5 cm in carbon grid
195 fFoilMaterial = "pefoam20";
196 fNFoils = 337;
197 fGapThick = 900 * 0.0001;
198 fFoilThick = 12 * 0.0001;
199 fnPhotonCorr = 0.65;
200 }
201 else if (fRadType == "B++" || fRadType == "K++") {
202 fFoilMaterial = "pokalon";
203 fNFoils = 350;
204 fGapThick = 0.07;
205 fFoilThick = 0.0024;
206 fnPhotonCorr = 0.65;
207 }
208 else if (fRadType == "G30") {
209 fFoilMaterial = "pefiber";
210 fNFoils = 1791;
211 fGapThick = 50 * 0.0001;
212 fFoilThick = 17 * 0.0001;
213 fnPhotonCorr = 0.65;
214 }
215 else if (fRadType == "tdr18_25cm") { //Like tdr18 but 5 cm shorter (reduced layer spacing version)
216 fFoilMaterial = "pefoam20";
217 fNFoils = 281;
218 fGapThick = 900 * 0.0001;
219 fFoilThick = 12 * 0.0001;
220 fnPhotonCorr = 0.65;
221 }
222 else { // If radiator typ is neither empty nor one of the above (because of typo or else) set defaults
223
224 LOG(warn) << "CbmTrdRadiator::Init : *********** Radiator prototype not known ***********";
225 LOG(warn) << "CbmTrdRadiator::Init : switch to default parameters ";
226 fFoilMaterial = "pefoam20";
227 fNFoils = 337;
228 fGapThick = 900 * 0.0001;
229 fFoilThick = 12 * 0.0001;
230 fnPhotonCorr = 0.65;
231 }
232
233
234 //Material has been specified now, either in the constructor or set by radiator type above. Now, go through materials and set properties accordingly
235
236
237 // material dependent parameters for the gas between the foils of the radiator
238 //Properties of the air gaps between foils
239 fGapDens = 0.001205; // [g/cm3]
240 fGapOmega = 0.7; //plasma frequency ( from Anton )
241 //material = fFoilMaterial;
242 if (fFoilMaterial == "pefoam20") {
243 // (polyethylen foam foil)
244 fFoilDens = 0.90; // [g/cm3 ]
245 fFoilOmega = 20.64; //plasma frequency ( from Cyrano )
246 }
247 else if (fFoilMaterial == "polyethylen") {
248 // (polyethylen)
249 fFoilDens = 0.92; // [g/cm3 ]
250 fFoilOmega = 20.9; //plasma frequency ( from Anton )
251 }
252 else if (fFoilMaterial == "pefiber") {
253 // (polyethylen)
254 fFoilDens = 0.90; // [g/cm3 ]
255 fFoilOmega = 20.64; //plasma frequency ( from Cyrano )
256 }
257 else if (fFoilMaterial == "mylar") {
258 // (Mylar)
259 fFoilDens = 1.393; // [g/cm3 ]
260 fFoilOmega = 24.53; //plasma frequency ( from Cyrano )
261 }
262 else if (fFoilMaterial == "pokalon") {
263 // (PokalonN470)
264 fFoilDens = 1.150; // [g/cm3 ]
265 fFoilOmega = 22.43; //plasma frequency ( from Cyrano )
266 }
267 else { // If material is empty or one of the above (because of typo or else) set default
268 LOG(warn) << "CbmTrdRadiator::Init : *********** Radiator material not known ***********";
269 LOG(warn) << "CbmTrdRadiator::Init : switch to default material ";
270 // (polyethylen foam foil)
271 fFoilMaterial = "pefoam20";
272 fFoilDens = 0.90; // [g/cm3 ]
273 fFoilOmega = 20.64; //plasma frequency ( from Cyrano )
274 }
275
276 // foil between radiator and gas chamber
277 if (fWindowFoil == "Kapton") { // 0.05 mu Aluminum is added in SigmaWin() for Kapton
278 fWinDens = 1.42; //[g/cm3] ?? define values ?? default mylar
279 fWinThick = 0.0025; // [cm] ?? define values ?? default mylar
280 }
281 else if (fWindowFoil != "Mylar") {
282 if (!WINDOW.Init(fWindowFoil)) fWindowFoil = "Mylar";
283 }
284 if (fWindowFoil == "Mylar") {
285 fWinDens = 1.39; //[g/cm3]
286 fWinThick = 0.0100; //[cm]
287 }
288 if (!WINDOW.fN) {
289 LOG(info) << "CbmTrdRadiator::Init : Window type : " << fWindowFoil;
290 LOG(info) << "CbmTrdRadiator::Init : Window density : " << fWinDens << " g/cm^3";
291 LOG(info) << "CbmTrdRadiator::Init : Window thickness : " << fWinThick << " cm";
292 }
293
294
295 LOG(info) << "CbmTrdRadiator::Init : Prototype : " << fRadType;
296 LOG(info) << "CbmTrdRadiator::Init : Scaling factor : " << fnPhotonCorr;
297 LOG(info) << "CbmTrdRadiator::Init : Foil material : " << fFoilMaterial;
298 LOG(info) << "CbmTrdRadiator::Init : Nr. of foils : " << fNFoils;
299 LOG(info) << "CbmTrdRadiator::Init : Foil thickness : " << setprecision(4) << fFoilThick << " cm";
300 LOG(info) << "CbmTrdRadiator::Init : Gap thickness : " << fGapThick << " cm";
301 LOG(info) << "CbmTrdRadiator::Init : Simple TR production: " << setprecision(2) << fSimpleTR;
302 LOG(info) << "CbmTrdRadiator::Init : Foil plasm. frequ. : " << fFoilOmega << " keV";
303 LOG(info) << "CbmTrdRadiator::Init : Gap plasm. frequ. : " << fGapOmega << " keV";
304
305
306 // Get all the gas properties from CbmTrdGas
307 CbmTrdGas* fTrdGas = CbmTrdGas::Instance();
308 if (fTrdGas == 0) {
309 fTrdGas = new CbmTrdGas();
310 fTrdGas->Init();
311 }
312
313 fCom2 = fTrdGas->GetNobleGas();
314 fCom1 = fTrdGas->GetCO2();
315 fDetType = fTrdGas->GetDetType();
316 fGasThick = fTrdGas->GetGasThick();
317
318 LOG(info) << "CbmTrdRadiator::Init : Detector noble gas : " << fTrdGas->GetNobleGas();
319 LOG(info) << "CbmTrdRadiator::Init : Detector quenching gas: " << fTrdGas->GetCO2();
320 LOG(info) << "CbmTrdRadiator::Init : Detector type : " << fTrdGas->GetDetType();
321 LOG(info) << "CbmTrdRadiator::Init : Detector gas thick. : " << fTrdGas->GetGasThick() << " cm";
322
323 //LOG(info) << "************* End of Trd Radiator Init **************";
324 // If simplified version is used one didn't calculate the TR
325 // for each CbmTrdPoint in full glory, but create at startup
326 // some histograms for several momenta which are later used to
327 // get the TR much faster.
328 if (fSimpleTR == kTRUE) {
329
331
332 TFile* oldFile = gFile;
333 TDirectory* oldDir = gDirectory;
334
335 TFile* f1 = new TFile("TRhistos.root", "recreate");
336
337 for (Int_t i = 0; i < fNMom; i++) {
338 fFinal[i]->Write();
339 }
340 f1->Close();
341 f1->Delete();
342
343 gFile = oldFile;
344 gDirectory = oldDir;
345 }
346}
347
348//----------------------------------------------------------------------------
350{
351 memset(fThick, 0, NCOMPONENTS * sizeof(Float_t));
352 memset(fDens, 0, NCOMPONENTS * sizeof(Float_t));
353 memset(GetMu, 0, NCOMPONENTS * sizeof(GetMuPtr));
354 Init(def);
355}
356
357//----------------------------------------------------------------------------
359{
360 TObjArray* mats = def.Tokenize(";");
361 for (fN = 0; fN < mats->GetEntriesFast(); fN++) {
362 TObjString* mat = (TObjString*) (*mats)[fN];
363 TString mname = mat->String();
364 if (mname.CompareTo("Al") == 0) {
365 fMat[fN] = "Al";
366 fDens[fN] = 2.7;
368 }
369 else if (mname.CompareTo("C") == 0) {
370 fMat[fN] = "C";
371 fDens[fN] = 2.267;
373 }
374 else if (mname.CompareTo("HC") == 0) {
375 fMat[fN] = "HC";
376 fDens[fN] = 4.8e-2; // 48 Kg/m^3
378 }
379 else if (mname.CompareTo("Po") == 0) {
380 fMat[fN] = "Po";
381 fDens[fN] = 0.92;
383 }
384 else if (mname.CompareTo("Pok") == 0) {
385 fMat[fN] = "Pok";
386 fDens[fN] = 1.150;
388 }
389 else if (mname.CompareTo("Ka") == 0) {
390 fMat[fN] = "Ka";
391 fDens[fN] = 1.39;
393 }
394 else if (mname.CompareTo("Air") == 0) {
395 fMat[fN] = "Air";
396 fDens[fN] = 1.205e-3;
398 }
399 else if (mname.CompareTo("Xe") == 0) {
400 fMat[fN] = "Xe";
401 fDens[fN] = 5.9e-3;
403 }
404 else if (mname.CompareTo("CO2") == 0) {
405 fMat[fN] = "CO2";
406 fDens[fN] = 1.9767e-3;
408 }
409 else if (mname.CompareTo("My") == 0) {
410 fMat[fN] = "My";
411 fDens[fN] = 1.393;
413 }
414 else {
415 LOG(warn) << "CbmTrdEntranceWindow::Init : material " << mname.Data() << " not in the DB. Can't create "
416 << def.Data() << " entrance window structure. Default to Mylar.";
417 fN = 0;
418 return kFALSE;
419 }
420 LOG(info) << "CbmTrdEntranceWindow::Init : C" << fN << " type : " << fMat[fN];
421 LOG(info) << "CbmTrdEntranceWindow::Init : density : " << fDens[fN] << " g/cm^3";
422 LOG(info) << "CbmTrdEntranceWindow::Init : thickness : " << fThick[fN] << " cm";
423 }
424 mats->Delete();
425 delete mats;
426 return kTRUE;
427}
428
429//----------------------------------------------------------------------------
431{
432 if (n > NCOMPONENTS) {
433 LOG(warn) << "CbmTrdEntranceWindow : can support only up to " << NCOMPONENTS
434 << " plane-parallel elements. Please check your simulations";
435 n = NCOMPONENTS;
436 }
437 memcpy(WINDOW.fThick, w, n * sizeof(Float_t));
438}
439
440//---- Lattice Hit------------------------------------------------------------
442{
443 //printf("---------------------------------------------------\n");
444 Double_t point_in[3];
445 Double_t point_out[3];
446 if (nullptr != point) {
447 point_in[0] = point->GetXIn();
448 point_in[1] = point->GetYIn();
449 point_in[2] = point->GetZIn();
450
451 point_out[0] = point->GetXOut();
452 point_out[1] = point->GetYOut();
453 point_out[2] = point->GetZOut();
454 Double_t back_direction[3] = {(point_in[0] - point_out[0]), (point_in[1] - point_out[1]),
455 (point_in[2] - point_out[2])};
456 if (back_direction[2] > 0) {
457 LOG(debug2) << "CbmTrdRadiator::LatticeHit: MC-track points towards target!";
458 //for(Int_t i = 0; i < 3; i++)
459 //printf("%i: in:%8.2f out:%8.2f direction:%8.2f\n",i,point_in[i],point_out[i],back_direction[i]);
460 return false;
461 }
462 Double_t trackLength = TMath::Sqrt(back_direction[0] * back_direction[0] + back_direction[1] * back_direction[1]
463 + back_direction[2] * back_direction[2]);
464
465 trackLength *= 10.; // cm -> mm to get a step width of 1mm
466 //printf("track length:%7.2fmm\n",trackLength);
467 gGeoManager->FindNode((point_out[0] + point_in[0]) / 2, (point_out[1] + point_in[1]) / 2,
468 (point_out[2] + point_in[2]) / 2);
469
470 Double_t pos[3] = {point_in[0], point_in[1], point_in[2]}; // start at entrance point
471
472 for (Int_t i = 0; i < 3; i++) {
473 back_direction[i] /= trackLength;
474 //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]);
475 }
476 if (TString(gGeoManager->GetPath()).Contains("gas")) {
477 Int_t stepCount = 0;
478 while (
479 /*(!TString(gGeoManager->GetPath()).Contains("lattice") || !TString(gGeoManager->GetPath()).Contains("radiator")) &&*/
480 pos[2] >= point_in[2] - 5 && stepCount < 50) {
481 stepCount++;
482 //printf("step%i\n",stepCount);
483 for (Int_t i = 0; i < 3; i++) {
484 pos[i] += back_direction[i];
485 //printf("%8.2f ",pos[i]);
486 }
487 //printf(" %s\n",TString(gGeoManager->GetPath()).Data());
488 gGeoManager->FindNode(pos[0], pos[1], pos[2]);
489 if (TString(gGeoManager->GetPath()).Contains("radiator")) {
490 //printf ("%s false\n",TString(gGeoManager->GetPath()).Data());
491 return false;
492 }
493 else if (TString(gGeoManager->GetPath()).Contains("lattice")) {
494 //printf ("%s true <----------------------------------------\n",TString(gGeoManager->GetPath()).Data());
495 return true;
496 }
497 else if (TString(gGeoManager->GetPath()).Contains("frame")) {
498 //printf ("%s true <----------------------------------------\n",TString(gGeoManager->GetPath()).Data());
499 return true;
500 }
501 else {
502 //printf ("%s\n",TString(gGeoManager->GetPath()).Data());
503 }
504 }
505 }
506 else {
507 LOG(error) << "CbmTrdRadiator::LatticeHit: MC-track not in TRD! Node:" << TString(gGeoManager->GetPath()).Data()
508 << " gGeoManager->MasterToLocal() failed!";
509 return false;
510 }
511 }
512 else {
513 LOG(error) << "CbmTrdRadiator::LatticeHit: CbmTrdPoint == nullptr!";
514 return false;
515 }
516 return kTRUE;
517}
518//----------------------------------------------------------------------------
519// ---- Spectra Production ---------------------------------------------------
521{
522
523 fTrackMomentum = new Double_t[fNMom];
524
525 Double_t trackMomentum[fNMom] = {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};
526
527 for (Int_t imom = 0; imom < fNMom; imom++) {
528 fTrackMomentum[imom] = trackMomentum[imom];
529 }
530
531 for (Int_t i = 0; i < fNMom; i++) {
532
533 fMom.SetXYZ(0.0, 0.0, fTrackMomentum[i]);
534
535 ProcessTR();
536 fnTRabs[i] = fnTRab;
537
538 // copy the content of the current fDetSpectrumA into appropriate
539 // fFinal
540
541 Float_t tmp = 0;
542 for (Int_t j = 0; j < fSpNBins; j++) {
543 tmp = fDetSpectrumA->GetBinContent(j + 1);
544 fFinal[i]->SetBinContent(j + 1, tmp);
545 }
546 fFinal[i]->SetTitle(Form("%s-momentum%.2f", fFinal[i]->GetTitle(), fTrackMomentum[i]));
547 }
548}
549
550//----------------------------------------------------------------------------
551
552// ---- GetTR ---------------------------------------------------------------
554{
555
556 /*
557 * Simplified simulation of the TR
558 */
559
560 if (fSimpleTR == kTRUE) {
561
562 Float_t m = Mom.Mag();
563
564 // Find correct spectra
565 Int_t i = 0;
566 while (m > fTrackMomentum[i]) {
567 i++;
568 if (i == fNMom) {
569 i--;
570 break;
571 }
572 }
573
574 ELoss(i);
575 }
576
577 /*
578 * Full simulation of the TR production
579 */
580 else {
581
582 fMom = Mom;
583 ProcessTR();
584 }
585
586 return fELoss;
587}
588
589//----------------------------------------------------------------------------
590
591//----- main TR processing function ------------------------------
593{
594
595 // Compute the angle normalization factor -
596 // for different angles the detector thicknesses are different
597
598 Float_t fNormPhi = fMom.Mag() / fMom.Pz();
599
600 // Correct the thickness according to this factor
601
602 fFoilThickCorr = TMath::Abs(fFoilThick * fNormPhi);
603 fGapThickCorr = TMath::Abs(fGapThick * fNormPhi);
604 fGasThickCorr = TMath::Abs(fGasThick * fNormPhi);
605
606 fFirstPass = true;
607 // compute the TR spectra in the radiator
608 TRspectrum();
609 // compute the TR spectra in the mylar foil
611 // compute the TR spectra in the detector
613 // Loop over photons and compute the E losses
614 ELoss(-1);
615
616
617 if (fDetType == 1) { // if detector is MB type = dual-sided MWPC
618 fFirstPass = false;
619 Float_t energiefirst = fELoss;
620 Float_t fTRAbsfirst = fnTRab;
621
622 // compute the TR spectra in the mylar foil between the two gas layers
624 // compute the TR spectra in the second halve of the detector
626 // Loop over photons and compute the E losses
627 ELoss(-1);
628 fELoss += energiefirst;
629 fnTRab += fTRAbsfirst;
630 }
631}
632//----------------------------------------------------------------------------
633
634// ---- Compute the energy losses --------------------------------------------
636{
637
638 // calculate the energy loss of the TR photons in the
639 // detector gas. fNTRab is the number of TR photons
640 // absorbed in the gas. Since this number is low the
641 // actual number comes from a poisson distribution.
642
643 Int_t nPhoton;
644 Int_t maxPhoton = 50;
645 Double32_t eLoss = 0;
646 Float_t eTR = 0;
647
648 if (-1 == index) {
649 if (0 == fDetSpectrumA) {
650 LOG(error) << " -Err :: No input spectra ";
651 fELoss = -1;
652 return 1;
653 }
654 nPhoton = gRandom->Poisson(fnTRab * fnPhotonCorr);
655 if (nPhoton > maxPhoton) nPhoton = maxPhoton;
656 for (Int_t i = 0; i < nPhoton; i++) {
657 // energy of the photon
658 eTR = fDetSpectrumA->GetRandom();
659 eLoss += eTR;
660 }
661 }
662 else {
663 if (0 == fFinal[index]) {
664 LOG(error) << " -Err :: No input spectra ";
665 fELoss = -1;
666 return 1;
667 }
668 nPhoton = gRandom->Poisson(fnTRabs[index] * fnPhotonCorr);
669 if (nPhoton > maxPhoton) nPhoton = maxPhoton;
670 for (Int_t i = 0; i < nPhoton; i++) {
671 // energy of the photon
672 eTR = fFinal[index]->GetRandom();
673 eLoss += eTR;
674 }
675 }
676
677 // keV->GeV & set the result
678 fELoss = eLoss / 1000000.0;
679
680 return 0;
681}
682
683//----------------------------------------------------------------------------
684
685//----- TR spectrum ---------------------------------------------------
687{
688
689 // calculate the number of TR photons in the radiator which are not
690 // absorbed in the radiator.
691 // Formulas are tacken from M. Castellano et al.
692 // Computer Physics Communications 61 (1990)
693
694
695 // Where does this values come from, put them in the header file
696 const Float_t kAlpha = 0.0072973; // 1/137
697 const Int_t kSumMax = 30;
698
699 // calculate the gamma of the particle
700 Float_t gamma = GammaF();
701
702 fSpectrum->Reset();
703 fDetSpectrum->Reset();
704 fDetSpectrumA->Reset();
705
707
708 SetSigma(1);
709
710 //Float_t stemp = 0;
711
712 // Loop over energy
713 // stemp is the total energy of all created photons
714 for (Int_t iBin = 0; iBin < fSpNBins; iBin++) {
715
716 // keV -> eV
717 Float_t energyeV = (fSpBinWidth * iBin + 1.0) * 1e3;
718
719 Float_t csFoil = fFoilOmega / energyeV;
720 Float_t csGap = fGapOmega / energyeV;
721
722 Float_t rho1 = energyeV * fFoilThickCorr * 1e4 * 2.5 * (1.0 / (gamma * gamma) + csFoil * csFoil);
723 Float_t rho2 = energyeV * fFoilThickCorr * 1e4 * 2.5 * (1.0 / (gamma * gamma) + csGap * csGap);
724
725 // Calculate the sum over the production angle tetan
726 Float_t sum = 0;
727 for (Int_t iSum = 0; iSum < kSumMax; iSum++) {
728 Float_t tetan = (TMath::Pi() * 2.0 * (iSum + 1) - (rho1 + kappa * rho2)) / (kappa + 1.0);
729 if (tetan < 0.0) { tetan = 0.0; }
730 Float_t aux = 1.0 / (rho1 + tetan) - 1.0 / (rho2 + tetan);
731 sum += tetan * (aux * aux) * (1.0 - TMath::Cos(rho1 + tetan));
732 }
733 Float_t conv = 1.0 - TMath::Exp(-fNFoils * fSigma[iBin]);
734
735 // eV -> keV
736 Float_t energykeV = energyeV * 0.001;
737
738 // dN / domega
739 Float_t wn = kAlpha * 4.0 / (fSigma[iBin] * (kappa + 1.0)) * conv * sum / energykeV;
740 /*
741 wn = 4.0 * kAlpha
742 / (energykeV * (kappa + 1.0))
743 * (1.0 - TMath::Exp(Double_t(-fNFoils) * fSigma[iBin])) / (1.0 - TMath::Exp(-fSigma[iBin]))
744 * sum;
745 */
746 // save the result
747 fSpectrum->SetBinContent(iBin + 1, wn);
748 // compute the integral
749 //stemp += wn;
750 }
751
752 // <nTR> (binsize corr.)
753 //Float_t nTR0 = stemp * fSpBinWidth;
754 //LOG(info) << " No. of photons after radiator" << nTR0;
755
756 // compute the spectra behind the mylar foil (absorption)
757 // WinTRspectrum();
758
759 return 1;
760}
761
762//------------------------------------------------------------------
763
764//----- WinTRspectrum -----------------------------------------------
766{
767 //
768 // Computes the spectrum after absorption in Mylar foil
769 //
770
771 fWinSpectrum->Reset();
772
773 SetSigma(2);
774
775 //Float_t stemp = 0;
776 for (Int_t iBin = 0; iBin < fSpNBins; iBin++) {
777
778 Float_t sp = 0;
779 if (fFirstPass == true) { sp = fSpectrum->GetBinContent(iBin + 1); }
780 else {
781 sp = fDetSpectrum->GetBinContent(iBin + 1);
782 }
783 // compute the absorption coefficient
784 Float_t conv = TMath::Exp(-fSigmaWin[iBin]);
785 Float_t wn = sp * conv;
786
787 fWinSpectrum->SetBinContent(iBin + 1, wn);
788
789 //stemp += wn;
790 }
791
792 //fnTRprod = stemp * fSpBinWidth;
793 //LOG(info) << "No. of photons after My = "<< fnTRprod;
794
795 // compute the spectrum absorbed in the D & escaped from the D
796 // DetTRspectrum();
797 return 1;
798}
799//----------------------------------------------------------------------------
800
801//----- DetTRspectrum ------------------------------------------------
803{
804 //
805 // Computes the spectrum absorbed and passed in the gas mixture
806 //
807
808 SetSigma(3);
809 //Float_t stemp = 0;
810 Float_t stemp2 = 0;
811 for (Int_t iBin = 0; iBin < fSpNBins; iBin++) {
812
813 //passed spectrum
814 Float_t sp = 0;
815 sp = fWinSpectrum->GetBinContent(iBin + 1);
816 Float_t conv = TMath::Exp(-fSigmaDet[iBin]);
817 Float_t wn = sp * conv;
818
819 fDetSpectrum->SetBinContent(iBin + 1, wn);
820 //stemp += wn;
821
822 // absorbed spectrum
823 Float_t conv2 = 1 - TMath::Exp(-fSigmaDet[iBin]);
824 Float_t wn2 = sp * conv2;
825
826 fDetSpectrumA->SetBinContent(iBin + 1, wn2);
827 stemp2 += wn2;
828 }
829
830 //Float_t nTR1 = stemp * fSpBinWidth;
831 Float_t nTR2 = stemp2 * fSpBinWidth;
832
833 // Save the number of photons absorbed
834 fnTRab = nTR2;
835
836 //LOG(info) << " No. of photons escaped: " << nTR1;
837 //LOG(info) << " No. of photnos absorbed in the gas: " <<nTR2;
838
839 return 1;
840}
841//----------------------------------------------------------------------------
842
843//----- Gamma factor--------------------------------------------------
845{
846
847 // put this in the header ?
848 Float_t massE = 5.1099907e-4; // GeV/c2
849
850 Float_t p2 = fMom.Mag2();
851 Float_t result = TMath::Sqrt(p2 + massE * massE) / massE;
852
853 return result;
854}
855//----------------------------------------------------------------------------
856
857//----- SetSigma---------------------------------------------------------
859{
860 //
861 // Sets the absorbtion crosssection for the energies of the TR spectrum
862 //
863
864
865 if (SigmaT == 1) {
866 if (fSigma) delete[] fSigma;
867 fSigma = new Float_t[fSpNBins];
868 for (Int_t iBin = 0; iBin < fSpNBins; iBin++) {
869 Float_t energykeV = iBin * fSpBinWidth + 1.0;
870 fSigma[iBin] = Sigma(energykeV);
871 }
872 }
873
874 if (SigmaT == 2) {
875 if (fSigmaWin) delete[] fSigmaWin;
877 for (Int_t iBin = 0; iBin < fSpNBins; iBin++) {
878 Float_t energykeV = iBin * fSpBinWidth + 1.0;
879 fSigmaWin[iBin] = SigmaWin(energykeV);
880 }
881 }
882
883 if (SigmaT == 3) {
884 if (fSigmaDet) delete[] fSigmaDet;
886 for (Int_t iBin = 0; iBin < fSpNBins; iBin++) {
887 Float_t energykeV = iBin * fSpBinWidth + 1.0;
888 fSigmaDet[iBin] = SigmaDet(energykeV);
889 }
890 }
891}
892//----------------------------------------------------------------------------
893
894//----- Sigma-------------------------------------------------------------
896{
897 //
898 // Calculates the absorbtion crosssection for a one-foil-one-gap-radiator
899 //
900
901 // keV -> MeV
902 Float_t energyMeV = energykeV * 0.001;
903 Float_t foil = 0.0;
904
905 if (fFoilMaterial == "polyethylen") foil = GetMuPo(energyMeV) * fFoilDens * fFoilThickCorr;
906 else if (fFoilMaterial == "pefoam20")
907 foil = GetMuPo(energyMeV) * fFoilDens * fFoilThickCorr;
908 else if (fFoilMaterial == "pefiber")
909 foil = GetMuPo(energyMeV) * fFoilDens * fFoilThickCorr;
910 else if (fFoilMaterial == "mylar")
911 foil = GetMuMy(energyMeV) * fFoilDens * fFoilThickCorr;
912 else if (fFoilMaterial == "pokalon")
913 foil = GetMuPok(energyMeV) * fFoilDens * fFoilThickCorr;
914 else
915 LOG(error) << "ERROR:: unknown radiator material";
916
917
918 if (energyMeV >= 0.001) {
919 Float_t result = (foil + (GetMuAir(energyMeV) * fGapDens * fGapThickCorr));
920 return result;
921 }
922 else {
923 return 1e6;
924 }
925}
926//----------------------------------------------------------------------------
927
928//----- SigmaWin -------------------------------------------------------
930{
931 //
932 // Calculates the absorbtion crosssection for a one-foil
933 //
934
935 // keV -> MeV
936 Float_t energyMeV = energykeV * 0.001;
937 if (energyMeV >= 0.001) {
938 if (fWindowFoil == "Kapton") // aluminized kapton
939 return GetMuKa(energyMeV) * fWinDens * fWinThick + GetMuAl(energyMeV) * 2.70 /*[g/cm^3]*/ * 5E-6 /*[cm]*/;
940 else if (fWindowFoil == "Mylar") // mylar foil
941 return GetMuMy(energyMeV) * fWinDens * fWinThick;
942 else { // general sandwich structure
943 Float_t sigma(0.);
944 for (Int_t iwin(0); iwin < WINDOW.fN; iwin++)
945 sigma += WINDOW.GetMu[iwin](energyMeV) * WINDOW.fDens[iwin] * WINDOW.fThick[iwin];
946 return sigma;
947 }
948 }
949 else
950 return 1e6;
951}
952//----------------------------------------------------------------------------
953
954//----- SigmaDet --------------------------------------------------------
956{
957 //
958 //Calculates the absorbtion crosssection for a choosed gas
959 //
960
961 // keV -> MeV
962 Float_t energyMeV = energykeV * 0.001;
963
964 if (energyMeV >= 0.001) {
965 // densities for CO2 and Xe changed. think density for CO2 is just
966 // a typo. where the difference for Xe comes from i don't know
967 // Values are from http://pdg.lbl.gov/AtomicNuclearProperties/
968 // return(GetMuCO2(energyMeV) * 0.001806 * fGasThick * fCom1 +
969 // GetMuXe(energyMeV) * 0.00589 * fGasThick * fCom2);
970 return (GetMuCO2(energyMeV) * 0.00184 * fGasThick * fCom1 + GetMuXe(energyMeV) * 0.00549 * fGasThick * fCom2);
971 }
972 else {
973 return 1e6;
974 }
975}
976//----------------------------------------------------------------------------
977
978//----- GetMuPok --------------------------------------------------------
980{
981 //
982 // Returns the photon absorbtion cross section for pokalon N470
983 //
984 return Interpolate(energyMeV, fEn_pok, fMu_pok, fKN_pok);
985}
986//----------------------------------------------------------------------------
987//----- GetMuKa --------------------------------------------------------
989{
990 //
991 // Returns the photon absorbtion cross section for kapton
992 //
993 return Interpolate(energyMeV, fEn_ka, fMu_ka, fKN_ka);
994}
995//----------------------------------------------------------------------------
996
997//----- GetMuAl --------------------------------------------------------
999{
1000 //
1001 // Returns the photon absorbtion cross section for Al
1002 //
1003 return Interpolate(energyMeV, fEn_al, fMu_al, fKN_al);
1004}
1005//----------------------------------------------------------------------------
1006
1007//----- GetMuC --------------------------------------------------------
1009{
1010 //
1011 // Returns the photon absorbtion cross section for Carbon
1012 //
1013 return Interpolate(energyMeV, fEn_c, fMu_c, fKN_c);
1014}
1015//----------------------------------------------------------------------------
1016
1017//----- GetMuPo --------------------------------------------------------
1019{
1020 //
1021 // Returns the photon absorbtion cross section for polypropylene
1022 //
1023 return Interpolate(energyMeV, fEn_po, fMu_po, fKN_po);
1024}
1025//----------------------------------------------------------------------------
1026
1027//----- GetMuAir --------------------------------------------------------
1029{
1030 //
1031 // Returns the photon absorbtion cross section for Air
1032 //
1033 return Interpolate(energyMeV, fEn_air, fMu_air, fKN_air);
1034}
1035//----------------------------------------------------------------------------
1036
1037//----- GetMuXe --------------------------------------------------------
1039{
1040 //
1041 // Returns the photon absorbtion cross section for xenon
1042 //
1043 return Interpolate(energyMeV, fEn_xe, fMu_xe, fKN_xe);
1044}
1045//----------------------------------------------------------------------------
1046
1047//----- GetMuCO2 ------------------------------------------------------
1049{
1050 //
1051 // Returns the photon absorbtion cross section for CO2
1052 //
1053 return Interpolate(energyMeV, fEn_co2, fMu_co2, fKN_co2);
1054}
1055//----------------------------------------------------------------------------
1056
1057//----- GetMuMy -------------------------------------------------------
1059{
1060 //
1061 // Returns the photon absorbtion cross section for mylar
1062 //
1063 return Interpolate(energyMeV, fEn_my, fMu_my, fKN_my);
1064}
1065//----------------------------------------------------------------------------
1066
1067//----- Interpolate ------------------------------------------------------
1069{
1070 //
1071 // Interpolates the photon absorbtion cross section
1072 // for a given energy <energyMeV>.
1073 //
1074
1075 Float_t de = 0;
1076 Int_t index = 0;
1077 Int_t istat = Locate(en, n, energyMeV, index, de);
1078 if (istat == 0) {
1079 // Float_t result = (mu[index] - de * (mu[index] - mu[index+1])
1080 // / (en[index+1] - en[index] ));
1081 // return result;
1082 Float_t result =
1083 (TMath::Log(mu[index]) - de * (TMath::Log(mu[index]) - TMath::Log(mu[index + 1])) / (en[index + 1] - en[index]));
1084 return TMath::Exp(result);
1085 }
1086 else {
1087 return 0.0;
1088 }
1089}
1090//----------------------------------------------------------------------------
1091
1092//----- Locate ------------------------------------------------------------
1094{
1095 //
1096 // Locates a point (xval) in a 1-dim grid (xv(n))
1097 //
1098
1099 if (xval >= xv[n - 1]) return 1;
1100 if (xval < xv[0]) return -1;
1101
1102 Int_t km;
1103 Int_t kh = n - 1;
1104
1105 kl = 0;
1106 while (kh - kl > 1) {
1107 if (xval < xv[km = (kl + kh) / 2]) kh = km;
1108 else
1109 kl = km;
1110 }
1111 if (xval < xv[kl] || xval > xv[kl + 1] || kl >= n - 1) {
1112 printf("Locate failed xv[%d] %f xval %f xv[%d] %f!!!\n", kl, xv[kl], xval, kl + 1, xv[kl + 1]);
1113 exit(1);
1114 }
1115
1116 dx = xval - xv[kl];
1117 if (xval == 0.001) LOG(info) << "Locat = 0";
1118 return 0;
1119}
1120//----------------------------------------------------------------------------
1121// Photon absorption cross sections
1122
1123// Pokalon N470
1124constexpr Float_t CbmTrdRadiator::fEn_pok[46];
1125constexpr Float_t CbmTrdRadiator::fMu_pok[46];
1126
1127// Kapton
1128constexpr Float_t CbmTrdRadiator::fEn_ka[46];
1129constexpr Float_t CbmTrdRadiator::fMu_ka[46];
1130
1131// Aluminum
1132constexpr Float_t CbmTrdRadiator::fEn_al[48];
1133constexpr Float_t CbmTrdRadiator::fMu_al[48];
1134
1135// Polypropylene/Polyethylene
1136constexpr Float_t CbmTrdRadiator::fEn_po[36];
1137constexpr Float_t CbmTrdRadiator::fMu_po[36];
1138
1139// Carbon
1140constexpr Float_t CbmTrdRadiator::fEn_c[46];
1141constexpr Float_t CbmTrdRadiator::fMu_c[46];
1142
1143// Air
1144constexpr Float_t CbmTrdRadiator::fEn_air[38];
1145constexpr Float_t CbmTrdRadiator::fMu_air[38];
1146
1147// Xenon
1148constexpr Float_t CbmTrdRadiator::fEn_xe[48];
1149constexpr Float_t CbmTrdRadiator::fMu_xe[48];
1150
1151// CO2
1152constexpr Float_t CbmTrdRadiator::fEn_co2[36];
1153constexpr Float_t CbmTrdRadiator::fMu_co2[36];
1154
1155// Mylar
1156constexpr Float_t CbmTrdRadiator::fEn_my[36];
1157constexpr Float_t CbmTrdRadiator::fMu_my[36];
1158
1159
ClassImp(CbmConverterManager)
Container for gas properties of TRD.
Float_t(* GetMuPtr)(Float_t)
#define NCOMPONENTS
float Float_t
int Int_t
bool Bool_t
Double_t GetNobleGas() const
Definition CbmTrdGas.h:28
Double_t GetGasThick() const
Definition CbmTrdGas.h:27
Int_t GetDetType() const
Definition CbmTrdGas.h:26
void Init()
Definition CbmTrdGas.cxx:43
static CbmTrdGas * Instance()
Definition CbmTrdGas.h:33
Double_t GetCO2() const
Definition CbmTrdGas.h:29
double GetYOut() const
Definition CbmTrdPoint.h:69
double GetXIn() const
Definition CbmTrdPoint.h:65
double GetZIn() const
Definition CbmTrdPoint.h:67
double GetXOut() const
Definition CbmTrdPoint.h:68
double GetYIn() const
Definition CbmTrdPoint.h:66
double GetZOut() const
Definition CbmTrdPoint.h:70
static constexpr Float_t fMu_c[46]
static const Int_t fKN_xe
static constexpr Float_t fMu_ka[46]
static constexpr Float_t fEn_co2[36]
static Int_t Locate(const Float_t *xv, Int_t n, Float_t xval, Int_t &kl, Float_t &dx)
Bool_t LatticeHit(const CbmTrdPoint *point)
Float_t SigmaDet(Float_t energykeV)
static const Int_t fNMom
TR passed through Detector.
static Float_t GetMuKa(Float_t energyMeV)
CbmTrdRadiator(Bool_t SimpleTR=true, Int_t Nfoils=337, Float_t FoilThick=0.0012, Float_t GapThick=0.09, TString material="pefoam20", TString window="Kapton")
Float_t Sigma(Float_t energykeV)
static constexpr Float_t fEn_c[46]
static Float_t GetMuC(Float_t energyMeV)
static Float_t GetMuPo(Float_t energyMeV)
void SetSigma(Int_t SigmaT)
static constexpr Float_t fEn_my[36]
Float_t fnTRabs[fNMom]
Absorption spectra for different momenta.
Float_t * fSigmaDet
[fSpNBins] Array of sigma values for the entrance window of detector
static const Int_t fKN_pok
static Float_t GetMuXe(Float_t energyMeV)
static const Int_t fKN_co2
static Float_t GetMuAir(Float_t energyMeV)
static constexpr Float_t fMu_my[36]
static constexpr Float_t fEn_pok[46]
CbmTrdEntranceWindow WINDOW
Float_t * fSigmaWin
[fSpNBins] Array of sigma values for the foil of the radiator
static const Int_t fKN_po
static constexpr Float_t fMu_al[48]
static constexpr Float_t fMu_pok[46]
static Float_t Interpolate(Float_t energyMeV, const Float_t *en, const Float_t *mu, Int_t n)
Float_t SigmaWin(Float_t energykeV)
static constexpr Float_t fMu_co2[36]
TH1D * fDetSpectrumA
TR spectra in gas-window foil.
TH1D * fDetSpectrum
TR absorbed in Detector.
static const Int_t fKN_ka
Float_t GetTR(TVector3 mom)
Double_t * fTrackMomentum
static constexpr Float_t fEn_ka[46]
static const Int_t fSpNBins
TH1D * fSpectrum
[fSpNBins] Array of sigma values for the active gas
Int_t ELoss(Int_t index)
TH1D * fWinSpectrum
TR photon energy spectrum.
static constexpr Float_t fMu_po[36]
static constexpr Float_t fMu_air[38]
void SetEWwidths(Int_t n, Float_t *w)
static const Int_t fSpRange
static Float_t GetMuAl(Float_t energyMeV)
static constexpr Float_t fEn_al[48]
static const Int_t fKN_c
static constexpr Float_t fEn_po[36]
static constexpr Float_t fMu_xe[48]
static const Int_t fKN_air
static Float_t GetMuMy(Float_t energyMeV)
virtual ~CbmTrdRadiator()
static constexpr Float_t fEn_xe[48]
static Float_t GetMuCO2(Float_t energyMeV)
static Float_t GetMuPok(Float_t energyMeV)
static const Int_t fKN_al
static constexpr Float_t fEn_air[38]
static const Int_t fKN_my
TH1D * fFinal[fNMom]
[fNMom] Track momenta for which spectra