CbmRoot
Loading...
Searching...
No Matches
CbmLitCheckEnergyLossMuons.cxx
Go to the documentation of this file.
1/* Copyright (C) 2009-2020 GSI/JINR-LIT, Darmstadt/Dubna
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Andrey Lebedev [committer] */
4
6
7#include "CbmDrawHist.h"
8#include "CbmUtils.h"
9#include "TCanvas.h"
13
14#include <TAxis.h>
15#include <TGraph.h>
16
17#include <boost/assign/list_of.hpp>
18
19#include <cstdlib>
20
21using boost::assign::list_of;
22
23CbmLitCheckEnergyLossMuons::CbmLitCheckEnergyLossMuons() : fMat("iron"), fOutputDir("./test/")
24{
25 fMom[0] = 47.04;
26 fMom[1] = 56.16;
27 fMom[2] = 68.02;
28 fMom[3] = 85.09;
29 fMom[4] = 100.3;
30 fMom[5] = 152.7;
31 fMom[6] = 176.4;
32 fMom[7] = 221.8;
33 fMom[8] = 286.8;
34 fMom[9] = 364.2;
35 fMom[10] = 391.7;
36 fMom[11] = 494.5;
37 fMom[12] = 899.5;
38 fMom[13] = 1101.;
39 fMom[14] = 1502.;
40 fMom[15] = 2103.;
41 fMom[16] = 3104.;
42 fMom[17] = 4104.;
43 fMom[18] = 8105.;
44 fMom[19] = 10110.;
45 fMom[20] = 14110.;
46 fMom[21] = 20110.;
47 fMom[22] = 30110.;
48 fMom[23] = 40110.;
49 fMom[24] = 80110.;
50 fMom[25] = 100100.;
51}
52
54
56{
58
60
61 if (fMat == "iron") {
63 }
64 else if (fMat == "tungsten") {
66 }
67 else if (fMat == "carbon") {
69 }
70
71 CalcEloss();
72
73 DrawGraphs();
74}
75
77{
78 TCanvas* c1 = new TCanvas("mean_energy_loss_muons", "mean_energy_loss_muons", 700, 500);
79
80 fTable[0]->GetXaxis()->SetLimits(45, 100102);
81 fTable[0]->SetMinimum(0.1);
82 fTable[0]->SetMaximum(10.);
83
84 DrawGraph(list_of(fTable[0])(fTable[1])(fTable[2])(fTable[3])(fCalc[0])(fCalc[1])(fCalc[2])(fCalc[3]),
85 list_of("total (table)")("ionization (table)")("bremsstrahlung (table)")("pair production (table)")(
86 "total (calculation)")("ionization (calculation)")("bremsstrahlung (calculation)")(
87 "pair production (calculation)"),
88 kLog, kLog, true, 0.20, 0.97, 0.9, 0.7);
89
91}
92
94{
95 CbmLitMaterialInfo material;
96 if (fMat == "tungsten") {
97 material.SetA(183.84);
98 material.SetZ(74.0);
99 material.SetRho(19.3);
100 material.SetRL(0.35);
101 }
102 else if (fMat == "iron") {
103 material.SetA(55.85);
104 material.SetZ(26.0);
105 material.SetRho(7.87);
106 material.SetRL(1.757622);
107 }
108 else if (fMat == "carbon") {
109 material.SetA(12.0107);
110 material.SetZ(6);
111 material.SetRho(2.265);
112 material.SetRL(18.8);
113 }
114 else {
115 exit(0);
116 }
117
120 for (Int_t i = 0; i < 26; i++) {
121 par.SetQp(1. / (fMom[i] * 1e-3));
122 fCalc[1]->SetPoint(i, fMom[i], me.BetheBloch(&par, &material) * 1e3);
123 fCalc[2]->SetPoint(i, fMom[i], me.BetheHeitler(&par, &material) * 1e3);
124 fCalc[3]->SetPoint(i, fMom[i], me.PairProduction(&par, &material) * 1e3);
125 fCalc[0]->SetPoint(i, fMom[i], me.dEdx(&par, &material) * 1e3);
126 }
127}
128
130{
131 for (int i = 0; i < 4; i++) {
132 fTable[i] = new TGraph();
133 fTable[i]->GetXaxis()->SetTitle("Momentum [MeV/c]");
134 fTable[i]->GetYaxis()->SetTitle("Energy loss [Mev*cm^2/g]");
135 fCalc[i] = new TGraph();
136 fCalc[i]->GetXaxis()->SetTitle("Momentum [MeV/c]");
137 fCalc[i]->GetYaxis()->SetTitle("Energy loss [Mev*cm^2/g]");
138 }
139}
140
142{
143 Double_t table_iron[26][4] = {
144 //ion //brems //pair //total
145 {5.494, 0.000, 0.000, 5.494}, {4.321, 0.000, 0.000, 4.321}, {3.399, 0.000, 0.000, 3.399},
146 {2.654, 0.000, 0.000, 2.654}, {2.274, 0.000, 0.000, 2.274}, {1.717, 0.000, 0.000, 1.717},
147
148 {1.616, 0.000, 0.000, 1.616}, {1.516, 0.000, 0.000, 1.516}, {1.463, 0.000, 0.000, 1.463},
149 {1.451, 0.000, 0.000, 1.451}, {1.452, 0.000, 0.000, 1.453}, {1.467, 0.000, 0.000, 1.467},
150 {1.548, 0.000, 0.000, 1.548},
151
152 {1.581, 0.001, 0.000, 1.582}, {1.635, 0.001, 0.000, 1.637}, {1.694, 0.002, 0.001, 1.697},
153 {1.760, 0.003, 0.002, 1.767}, {1.806, 0.004, 0.004, 1.816}, {1.911, 0.010, 0.011, 1.936},
154
155 {1.942, 0.014, 0.015, 1.975}, {1.987, 0.021, 0.024, 2.039}, {2.032, 0.033, 0.039, 2.113},
156 {2.080, 0.054, 0.068, 2.214}, {2.112, 0.076, 0.099, 2.303}, {2.184, 0.171, 0.236, 2.623},
157 {2.207, 0.221, 0.310, 2.777}};
158
159 for (Int_t i = 0; i < 26; i++) {
160 fTable[1]->SetPoint(i, fMom[i], table_iron[i][0]);
161 fTable[2]->SetPoint(i, fMom[i], table_iron[i][1]);
162 fTable[3]->SetPoint(i, fMom[i], table_iron[i][2]);
163 fTable[0]->SetPoint(i, fMom[i], table_iron[i][3]);
164 }
165}
166
168{
169 Double_t table_tungsten[26][4] = {
170 //ion //brems //pair //total
171 {4.000, 0.000, 0.000, 4.00}, {3.185, 0.000, 0.000, 3.185}, {2.534, 0.000, 0.000, 2.534},
172 {2.000, 0.000, 0.000, 2.000}, {1.726, 0.000, 0.000, 1.726}, {1.323, 0.000, 0.000, 1.323},
173
174 {1.251, 0.000, 0.000, 1.251}, {1.182, 0.000, 0.000, 1.182}, {1.149, 0.000, 0.000, 1.149},
175 {1.145, 0.000, 0.000, 1.145}, {1.150, 0.000, 0.000, 1.150}, {1.168, 0.000, 0.000, 1.168},
176 {1.247, 0.001, 0.000, 1.249},
177
178 {1.279, 0.001, 0.000, 1.281}, {1.329, 0.002, 0.001, 1.333}, {1.384, 0.004, 0.002, 1.391},
179 {1.446, 0.006, 0.005, 1.459}, {1.489, 0.009, 0.009, 1.509}, {1.587, 0.023, 0.026, 1.640},
180
181 {1.616, 0.031, 0.036, 1.687}, {1.658, 0.047, 0.057, 1.768}, {1.700, 0.073, 0.092, 1.874},
182 {1.745, 0.120, 0.159, 2.035}, {1.774, 0.169, 0.231, 2.190}, {1.840, 0.383, 0.548, 2.800},
183 {1.860, 0.496, 0.717, 3.110}};
184
185 for (Int_t i = 0; i < 26; i++) {
186 fTable[1]->SetPoint(i, fMom[i], table_tungsten[i][0]);
187 fTable[2]->SetPoint(i, fMom[i], table_tungsten[i][1]);
188 fTable[3]->SetPoint(i, fMom[i], table_tungsten[i][2]);
189 fTable[0]->SetPoint(i, fMom[i], table_tungsten[i][3]);
190 }
191}
192
194{
195 Double_t table_carbon[26][4] = {
196 //ion //brems //pair //total
197 {7.119, 0.000, 0.000, 7.119}, {5.551, 0.000, 0.000, 5.551}, {4.332, 0.000, 0.000, 4.332},
198 {3.356, 0.000, 0.000, 3.356}, {2.861, 0.000, 0.000, 2.861}, {2.127, 0.000, 0.000, 2.127},
199
200 {1.992, 0.000, 0.000, 1.992}, {1.854, 0.000, 0.000, 1.854}, {1.775, 0.000, 0.000, 1.775},
201 {1.745, 0.000, 0.000, 1.745}, {1.745, 0.000, 0.000, 1.745}, {1.751, 0.000, 0.000, 1.751},
202 {1.819, 0.000, 0.000, 1.820},
203
204 {1.850, 0.000, 0.000, 1.851}, {1.900, 0.000, 0.000, 1.901}, {1.955, 0.000, 0.000, 1.957},
205 {2.018, 0.001, 0.001, 2.021}, {2.062, 0.001, 0.001, 2.066}, {2.161, 0.003, 0.003, 2.171},
206
207 {2.191, 0.004, 0.004, 2.204}, {2.234, 0.006, 0.007, 2.254}, {2.278, 0.010, 0.011, 2.308},
208 {2.326, 0.016, 0.019, 2.374}, {2.359, 0.022, 0.028, 2.427}, {2.434, 0.050, 0.068, 2.587},
209 {2.458, 0.065, 0.089, 2.655}};
210
211 for (Int_t i = 0; i < 26; i++) {
212 fTable[1]->SetPoint(i, fMom[i], table_carbon[i][0]);
213 fTable[2]->SetPoint(i, fMom[i], table_carbon[i][1]);
214 fTable[3]->SetPoint(i, fMom[i], table_carbon[i][2]);
215 fTable[0]->SetPoint(i, fMom[i], table_carbon[i][3]);
216 }
217}
218
void SetDefaultDrawStyle()
void DrawGraph(TGraph *graph, HistScale logx, HistScale logy, const string &drawOpt, Int_t color, Int_t lineWidth, Int_t lineStyle, Int_t markerSize, Int_t markerStyle)
Helper functions for drawing 1D and 2D histograms and graphs.
@ kLog
Definition CbmDrawHist.h:68
ClassImp(CbmLitCheckEnergyLossMuons)
Calculation of multiple scattering and energy loss.
Data class for track parameters.
Calculation of multiple scattering and energy loss.
litfloat BetheHeitler(const CbmLitTrackParam *par, const CbmLitMaterialInfo *mat) const
litfloat dEdx(const CbmLitTrackParam *par, const CbmLitMaterialInfo *mat) const
litfloat PairProduction(const CbmLitTrackParam *par, const CbmLitMaterialInfo *mat) const
litfloat BetheBloch(const CbmLitTrackParam *par, const CbmLitMaterialInfo *mat) const
void SetA(litfloat A)
void SetRL(litfloat rl)
void SetRho(litfloat rho)
void SetZ(litfloat Z)
Data class for track parameters.
void SetQp(litfloat qp)
void SaveCanvasAsImage(TCanvas *c, const string &dir, const string &option)
Definition CbmUtils.cxx:36