BmnRoot
Loading...
Searching...
No Matches
BmnTof700Qa.cxx
Go to the documentation of this file.
1#include "BmnTof700Qa.h"
2#include <TStopwatch.h>
3#include "BmnGlobalTrack.h"
4#include <TFile.h>
5#include "UniRun.h"
6#include "BmnFieldMap.h"
7#include "BmnNewFieldMap.h"
8#include "CbmStsTrack.h"
9#include "TMath.h"
10#include <fstream>
11
12using namespace TMath;
13
14static Double_t fWorkTime = 0.0;
15
17 fStsTracks = nullptr;
18 fGlobalTracks = nullptr;
19 fKalman = nullptr;
20 fTof700Hits = nullptr;
21 fHistoManager = new BmnQaHistoManager();
22 fKalmanType = 2; //1 - L1 KF, 2 - Lebedev's KF
23 fIsField = kTRUE;
24 fEventId = 0;
25}
26
27InitStatus BmnTof700Qa::Init() {
28
29 TStopwatch sw;
30 sw.Start();
31
32 FairRootManager* ioman = FairRootManager::Instance();
33 if (!ioman)
34 Fatal("Init", "FairRootManager is not instantiated");
35
36 /*CbmStsDigiScheme* digiScheme = */CbmStsDigiScheme::Instance(4); //for L1 KF
37 fKalman = new BmnKalmanFilter();
38
39 fTof700Hits = (TClonesArray*)ioman->GetObject("BmnTof700Hit");
40 fStsTracks = (TClonesArray*)ioman->GetObject("StsVector");
41 fStsHits = (TClonesArray*)ioman->GetObject("StsHit");
42 fGlobalTracks = (TClonesArray*)ioman->GetObject("BmnGlobalTrack");
43 fVertex = (CbmVertex*)ioman->GetObject("PrimaryVertex.");
44
45 CreateHistograms();
46
47
48 ifstream infile("tof700_strip_time_shifts.txt");
49
50 Int_t m, s, n;
51 Float_t sh;
52 while (infile >> m >> s >> n >> sh) {
53 timeShifts[m][s] = sh;
54 }
55
56 sw.Stop();
57 fWorkTime += sw.RealTime();
58 return kSUCCESS;
59}
60
61void BmnTof700Qa::Exec(Option_t* opt) {
62
63 TStopwatch sw;
64 sw.Start();
65
66 // Int_t modOrder[58] = {
67 // 0, 1, 4, 7, 10, 13, 6, 30, 24, 32, 26, 34, 28, 52,
68 // 14, 15, 16, 17, 18, 19, 36, 39, 37, 40, 38, 41,
69 // 3, 12, 2, 5, 8, 11, 20, 31, 25, 33, 27, 35, 29, 51,
70 // 42, 43, 44, 45, 46, 47,
71 // 53, 54, 55, 56, 57, 58,
72 // 22, 21, 23,
73 // 48, 49, 50
74 // };
75
76
77 Float_t dx[nModules] = {
78 +0.13, +0.15, -0.21, -0.28, +0.09, -0.29, -0.13, -0.28, -0.14, -0.08,
79 -0.07, -0.13, -0.30, -0.26, +0.36, +0.24, +0.33, +0.46, -0.06, +0.09,
80 -0.07, +0.36, -0.23, +0.18, +0.05, -0.38, +0.16, +0.26, +0.28, -0.21,
81 -0.06, -0.29, -0.03, -0.45, +0.04, -0.33, +0.16, +0.15, +0.38, +0.19,
82 +0.54, +0.68, -0.16, -0.04, +0.10, +0.23, -0.13, +0.00, -0.10, +0.56,
83 +0.50, -0.18, +0.18, +0.44, +0.50, +0.39, +0.29, +0.37, +0.01
84 };
85
86 Float_t dy[nModules] = {
87 -0.06, -0.08, -0.32, -0.25, -0.13, -0.31, -0.18, -0.19, -0.34, -0.45,
88 -0.21, -0.34, -0.22, -0.21, -0.09, -0.22, -0.29, -0.28, -0.02, -0.15,
89 -0.34, -0.17, -0.08, -0.01, -0.05, -0.09, +0.01, -0.06, +0.12, +0.00,
90 -0.08, -0.20, -0.05, -0.07, +0.13, +0.00, -0.03, +0.02, +0.09, -0.07,
91 +0.13, +0.21, -0.12, -0.16, -0.04, -0.02, +0.05, -0.04, -0.10, +0.19,
92 +0.10, +0.04, +0.16, -0.01, -0.25, -0.29, -0.16, -0.02, -0.20
93 };
94
95 if (!IsActive()) return;
96
97 cout << "Event ID: " << fEventId++ << endl;
98
99 if (fVertex->GetNTracks() < 2) return;
100 if (fVertex->GetZ() < -0.5 || fVertex->GetZ() > 0.5) return;
101 if (fVertex->GetX() < -1.0 || fVertex->GetX() > 1.5) return;
102 if (fVertex->GetY() < -1.0 || fVertex->GetY() > 1.2) return;
103
104 //for STS
105 for (Int_t iTr = 0; iTr < fStsTracks->GetEntriesFast(); ++iTr) {
106 CbmStsTrack* track = (CbmStsTrack*)fStsTracks->At(iTr);
107
108 if (track->GetNStsHits() < 5) continue;
109
110 FairTrackParam parF = *(track->GetParamFirst());
111 Int_t pdg = 2212;
112 Float_t mom = 0.0;
113 if (fIsField) {
114 mom = 1.0 / parF.GetQp();
115 pdg = (mom > 0) ? 2212 : -211;
116 }
117
118 // FairTrackParam parL1 = *(track->GetParamLast());
119 // FairTrackParam parL2 = *(track->GetParamLast());
120
121 // PropagateToZ(&parL1, 400, pdg, 1);
122 // PropagateToZ(&parL2, 400, pdg, 2);
123
124 // fHistoManager->H2("Compare_X")->Fill(parL1.GetX(), parL2.GetX());
125 // fHistoManager->H2("Compare_Y")->Fill(parL1.GetY(), parL2.GetY());
126 // fHistoManager->H2("Compare_Tx")->Fill(parL1.GetTx(), parL2.GetTx());
127 // fHistoManager->H2("Compare_Ty")->Fill(parL1.GetTy(), parL2.GetTy());
128 // fHistoManager->H2("Compare_P")->Fill(1.0 / parL1.GetQp(), 1.0 / parL2.GetQp());
129
130
131 Double_t length = 0.0;
132
133 if (PropagateToZ(&parF, fVertex->GetZ(), pdg, fKalmanType, &length) == kFALSE) continue;
134
135 //if (fKalman->TGeoTrackPropagate(&parF, fVertex->GetZ(), pdg, nullptr, nullptr, kTRUE) == kBMNERROR) continue;
136 if (TMath::Sqrt(TMath::Sq(parF.GetX() - fVertex->GetX()) + TMath::Sq(parF.GetY() - fVertex->GetY())) > 1.0) continue;
137
138 //ADD LOOP OVER HITS TO CALCULATE CBMSTSTRACK LENGTH ====>
139 CbmStsHit* cbmhit = (CbmStsHit*)fStsHits->At(track->GetStsHitIndex(0));
140 Double_t xPrev = cbmhit->GetX();
141 Double_t yPrev = cbmhit->GetY();
142 Double_t zPrev = cbmhit->GetZ();
143 for (Int_t iHit = 1; iHit < track->GetNStsHits(); ++iHit) {
144 cbmhit = (CbmStsHit*)fStsHits->At(track->GetStsHitIndex(iHit));
145 length += Sqrt(Sq(cbmhit->GetX() - xPrev) + Sq(cbmhit->GetY() - yPrev) + Sq(cbmhit->GetZ() - zPrev));
146 xPrev = cbmhit->GetX();
147 yPrev = cbmhit->GetY();
148 zPrev = cbmhit->GetZ();
149 }
150 //<====
151
152
153 FairTrackParam parL = *(track->GetParamLast());
154 Float_t minResX = 100.0;
155 Float_t minResY = 100.0;
156 Int_t minIdx = -1;
157 for (Int_t iH = 0; iH < fTof700Hits->GetEntriesFast(); ++iH) {
158 BmnTofHit* hit = (BmnTofHit*)fTof700Hits->At(iH);
159 Int_t mod = ((hit->GetDetectorID() & 0x0000FF00) >> 8) - 1;
160 //Double_t hitX = hit->GetX();
161 //Double_t hitY = hit->GetY();
162 Double_t hitZ = hit->GetZ();// -dz[mod];
163 PropagateToZ(&parL, hitZ, pdg, fKalmanType);
164 //if (fKalman->TGeoTrackPropagate(&parL, hitZ, pdg, nullptr, nullptr) == kBMNERROR) continue;
165 Float_t resX = parL.GetX() - hit->GetX() - GetDxMom(&parL) - dx[mod];
166 Float_t resY = parL.GetY() - hit->GetY() - GetDyMom(&parL) - dy[mod];
167 if (Abs(resX) > 3.0 * GetSigxMom(&parL) || Abs(resY) > 3.0 * GetSigxMom(&parL)) continue;
168 if (Abs(resX) < minResX && Abs(resY) < minResY) {
169 minResX = resX;
170 minResY = resY;
171 minIdx = iH;
172 }
173 }
174 if (minIdx == -1) continue;
175
176 FairTrackParam minParL = *(track->GetParamLast());
177 BmnTofHit* minHit = (BmnTofHit*)fTof700Hits->At(minIdx);
178 if (PropagateToZ(&minParL, minHit->GetZ(), pdg, fKalmanType, &length) == kFALSE) continue;
179 Int_t mod = ((minHit->GetDetectorID() & 0x0000FF00) >> 8) - 1;
180 Int_t strip = (minHit->GetDetectorID() & 0x000000FF) - 1;
181
182 // fHistoManager->H1(Form("TOF700_ResidX_%d", mod))->Fill(minResX);
183 // fHistoManager->H1(Form("TOF700_ResidY_%d", mod))->Fill(minResY);
184 // fHistoManager->H2(Form("TOF700_dx_tx_%d", mod))->Fill(minParL.GetTx(), minResX);
185
186 if (fIsField) {
187
188 mom = 1.0 / parF.GetQp();
189 Float_t time = minHit->GetTimeStamp() + timeShifts[mod][strip];
190 Float_t beta = length / time / (TMath::C() * 1e-7);
191 Float_t m2 = mom * mom * (1.0 / beta / beta - 1.0);
192
193 fHistoManager->H2("TOF700_Banan_total")->Fill(mom, beta);
194 fHistoManager->H2("TOF700_M2_total")->Fill(mom, m2);
195 fHistoManager->H2("TOF700_dx_mom_total")->Fill(mom, minResX);
196 fHistoManager->H2("TOF700_dy_mom_total")->Fill(mom, minResY);
197
198
199 //cout << minHit->GetLength() << " " << minHit->GetTimeStamp() << " " << TMath::C() * 1e-7 << endl;
200 //cout << minResX << " " << minResY << " " << mom << " " << m2 << " " << beta << endl;
201
202 if (mom > 0.0 && mom < 5.0 && m2 > 0.5 && m2 < 1.5) { //p
203 // fHistoManager->H1("TOF700_proton_ResidX_total")->Fill(resX);
204 // fHistoManager->H1("TOF700_proton_ResidY_total")->Fill(resY);
205 fHistoManager->H1(Form("TOF700_proton_ResidX_%d", mod))->Fill(minResX);
206 fHistoManager->H1(Form("TOF700_proton_ResidY_%d", mod))->Fill(minResY);
207 fHistoManager->H2("TOF700_proton_dx_mom_total")->Fill(mom, minResX);
208 fHistoManager->H2("TOF700_proton_dy_mom_total")->Fill(mom, minResY);
209 //fHistoManager->H2(Form("TOF700_proton_Banan_%d", mod))->Fill(mom, beta);
210 // fHistoManager->H2("TOF700_proton_Mass2_modules")->Fill(mod, m2);
211 // fHistoManager->H2(Form("TOF700_proton_Mass2_strips_%d", mod))->Fill(strip, m2);
212 //fHistoManager->H2("TOF700_proton_Banan_total")->Fill(mom, beta);
213 // fHistoManager->H2("TOF700_proton_Mass2_total")->Fill(mom, m2);
214 // fHistoManager->H2("TOF700_proton_dx_tx_total")->Fill(par.GetTx(), resX);
215 fHistoManager->H2(Form("TOF700_proton_dx_tx_%d", mod))->Fill(minParL.GetTx(), minResX);
216 }
217
218 if (mom > 0.0 && mom < 6.0 && m2 > 3.0 && m2 < 5.5) { //d
219 // fHistoManager->H1("TOF700_deutron_ResidX_total")->Fill(resX);
220 // fHistoManager->H1("TOF700_deutron_ResidY_total")->Fill(resY);
221 fHistoManager->H1(Form("TOF700_deutron_ResidX_%d", mod))->Fill(minResX);
222 fHistoManager->H1(Form("TOF700_deutron_ResidY_%d", mod))->Fill(minResY);
223 //fHistoManager->H2(Form("TOF700_deutron_Banan_%d", mod))->Fill(mom, beta);
224 fHistoManager->H2("TOF700_deutron_dx_mom_total")->Fill(mom, minResX);
225 fHistoManager->H2("TOF700_deutron_dy_mom_total")->Fill(mom, minResY);
226 // fHistoManager->H2("TOF700_deutron_Mass2_modules")->Fill(mod, m2);
227 // fHistoManager->H2(Form("TOF700_deutron_Mass2_strips_%d", mod))->Fill(strip, m2);
228 // fHistoManager->H2("TOF700_deutron_Banan_total")->Fill(mom, track->GetBeta(2));
229 // fHistoManager->H2("TOF700_deutron_Mass2_total")->Fill(mom, m2);
230 // fHistoManager->H2("TOF700_deutron_dx_tx_total")->Fill(par.GetTx(), resX);
231 fHistoManager->H2(Form("TOF700_deutron_dx_tx_%d", mod))->Fill(minParL.GetTx(), minResX);
232 }
233
234 if (mom > 0.0 && mom < 3.5 && m2 > -0.2 && m2 < 0.2) { //pi+
235 // fHistoManager->H1("TOF700_piPlus_ResidX_total")->Fill(resX);
236 // fHistoManager->H1("TOF700_piPlus_ResidY_total")->Fill(resY);
237 fHistoManager->H1(Form("TOF700_piPlus_ResidX_%d", mod))->Fill(minResX);
238 fHistoManager->H1(Form("TOF700_piPlus_ResidY_%d", mod))->Fill(minResY);
239 //fHistoManager->H2(Form("TOF700_piPlus_Banan_%d", mod))->Fill(mom, beta);
240 fHistoManager->H2("TOF700_piPlus_dx_mom_total")->Fill(mom, minResX);
241 fHistoManager->H2("TOF700_piPlus_dy_mom_total")->Fill(mom, minResY);
242 // fHistoManager->H2("TOF700_piPlus_Mass2_modules")->Fill(mod, m2);
243 // fHistoManager->H2(Form("TOF700_piPlus_Mass2_strips_%d", mod))->Fill(strip, m2);
244 // fHistoManager->H2("TOF700_piPlus_Banan_total")->Fill(mom, track->GetBeta(2));
245 // fHistoManager->H2("TOF700_piPlus_Mass2_total")->Fill(mom, m2);
246 // fHistoManager->H2("TOF700_piPlus_dx_tx_total")->Fill(par.GetTx(), resX);
247 fHistoManager->H2(Form("TOF700_piPlus_dx_tx_%d", mod))->Fill(minParL.GetTx(), minResX);
248 }
249
250 if (mom > -3.5 && mom < 0.0 && m2 > -0.2 && m2 < 0.2) { //pi-
251 // fHistoManager->H1("TOF700_piMinus_ResidX_total")->Fill(resX);
252 // fHistoManager->H1("TOF700_piMinus_ResidY_total")->Fill(resY);
253 fHistoManager->H1(Form("TOF700_piMinus_ResidX_%d", mod))->Fill(minResX);
254 fHistoManager->H1(Form("TOF700_piMinus_ResidY_%d", mod))->Fill(minResY);
255 //fHistoManager->H2(Form("TOF700_piMinus_Banan_%d", mod))->Fill(mom, beta);
256 fHistoManager->H2("TOF700_piMinus_dx_mom_total")->Fill(mom, minResX);
257 fHistoManager->H2("TOF700_piMinus_dy_mom_total")->Fill(mom, minResY);
258 // fHistoManager->H2("TOF700_piMinus_Mass2_modules")->Fill(mod, m2);
259 // fHistoManager->H2(Form("TOF700_piMinus_Mass2_strips_%d", mod))->Fill(strip, m2);
260 // fHistoManager->H2("TOF700_piMinus_Banan_total")->Fill(mom, track->GetBeta(2));
261 // fHistoManager->H2("TOF700_piMinus_Mass2_total")->Fill(mom, m2);
262 // fHistoManager->H2("TOF700_piMinus_dx_tx_total")->Fill(par.GetTx(), resX);
263 fHistoManager->H2(Form("TOF700_piMinus_dx_tx_%d", mod))->Fill(minParL.GetTx(), minResX);
264 }
265 }
266 }
267
268 sw.Stop();
269 fWorkTime += sw.RealTime();
270
271}
272
274
275 TStopwatch sw;
276 sw.Start();
277
278 FairRootManager::Instance()->GetOutFile()->cd();
279 TDirectory* tof700Dir = FairRootManager::Instance()->GetOutFile()->mkdir("TOF700");
280
281 // for (Int_t i = 0; i < nModules; ++i) {
282 // //fHistoManager->H1(Form("TOF700_dx_tx_%d", i))->Fit("pol1", "QSR", "", -0.5, 0.5);
283 // //fHistoManager->H1(Form("TOF700_ResidX_%d", i))->Fit("gaus", "QSR", "", -2, 2);
284 // //fHistoManager->H1(Form("TOF700_ResidY_%d", i))->Fit("gaus", "QSR", "", -2, 2);
285 // }
286
287 tof700Dir->cd();
288 fHistoManager->WriteToFile();
289
290 sw.Stop();
291 fWorkTime += sw.RealTime();
292 printf("Work time of BmnTof700Qa: %4.2f sec.\n", fWorkTime);
293}
294
295void BmnTof700Qa::PropagateToZ(FairTrackParam* par, Double_t zDst, Int_t pdg, Int_t kalmanType) {
296 if (kalmanType == 1) {//L1 KF
297 CbmKFTrack kftr = CbmKFTrack(*par);
298 kftr.SetPID(pdg);
299 kftr.Extrapolate(zDst);
300 kftr.GetTrackParam(*par);
301 } else if (kalmanType == 2) { //Lebedev's KF
302 fKalman->TGeoTrackPropagate(par, zDst, pdg, nullptr, nullptr);
303 }
304}
305
306Bool_t BmnTof700Qa::PropagateToZ(FairTrackParam* par, Double_t zDst, Int_t pdg, Int_t kalmanType, Double_t* len) {
307 Double_t lenPrev = *len;
308 CbmKFTrack kftr = CbmKFTrack(*par);
309 kftr.SetPID(pdg);
310 Int_t n = 10; //number of steps
311 Double_t zLength = zDst - par->GetZ();
312 Double_t step = zLength / n;
313 Double_t xPrev = par->GetX();
314 Double_t yPrev = par->GetY();
315 Double_t zPrev = par->GetZ();
316 Double_t z = zPrev;
317 for (Int_t i = 0; i < n; ++i) {
318 z += step;
319 if (kalmanType == 1) {//L1 KF
320 kftr.Extrapolate(z);
321 kftr.GetTrackParam(*par);
322 } else if (kalmanType == 2) { //Lebedev's KF
323 if (fKalman->TGeoTrackPropagate(par, z, pdg, nullptr, nullptr) == kBMNERROR)return kFALSE;
324 }
325 if (len)
326 *len += Sqrt(Sq(par->GetX() - xPrev) + Sq(par->GetY() - yPrev) + Sq(par->GetZ() - zPrev));
327 xPrev = par->GetX();
328 yPrev = par->GetY();
329 zPrev = par->GetZ();
330 }
331 if (Abs(*len - lenPrev) < 1.0) return kFALSE;
332 else return kTRUE;
333}
334
335void BmnTof700Qa::CreateHistograms() {
336
337 // fHistoManager->Create2 <TH2F>("Compare_X", "; X, L1; X, Leb", 400, -100.0, 100.0, 400, -100.0, 100.0);
338 // fHistoManager->Create2 <TH2F>("Compare_Y", "; Y, L1; Y, Leb", 400, -100.0, 100.0, 400, -100.0, 100.0);
339 // fHistoManager->Create2 <TH2F>("Compare_Tx", "; Tx, L1; Tx, Leb", 400, -1.0, 1.0, 400, -1.0, 1.0);
340 // fHistoManager->Create2 <TH2F>("Compare_Ty", "; Ty, L1; Ty, Leb", 400, -1.0, 1.0, 400, -1.0, 1.0);
341 // fHistoManager->Create2 <TH2F>("Compare_P", "; P, L1; P, Leb", 400, -5.0, 10.0, 400, -5.0, 10.0);
342
343 fHistoManager->Create2 <TH2F>("TOF700_Banan_total", "Total #beta vs rigidity; rigidity, GeV/c/Q; #beta", 800, -5.0, 10.0, 800, 0.2, 1.1);
344 fHistoManager->Create2 <TH2F>("TOF700_M2_total", "Total m^{2} vs rigidity; rigidity, GeV/c/Q; m^{2}, GeV/c^{2}", 800, -5.0, 10.0, 800, -0.5, 10);
345 fHistoManager->Create2 <TH2F>("TOF700_dx_mom_total", "Total dX vs rigidity; rigidity, GeV/c/Q; dX, cm", 800, -5.0, 10.0, 800, -5.0, 5.0);
346 fHistoManager->Create2 <TH2F>("TOF700_dy_mom_total", "Total dY vs rigidity; rigidity, GeV/c/Q; dY, cm", 800, -5.0, 10.0, 800, -5.0, 5.0);
347
348 fHistoManager->Create2 <TH2F>("TOF700_proton_dx_mom_total", "Total dX vs rigidity; rigidity, GeV/c/Q; dX, cm", 800, -5.0, 10.0, 800, -5.0, 5.0);
349 fHistoManager->Create2 <TH2F>("TOF700_proton_dy_mom_total", "Total dY vs rigidity; rigidity, GeV/c/Q; dY, cm", 800, -5.0, 10.0, 800, -5.0, 5.0);
350
351 fHistoManager->Create2 <TH2F>("TOF700_deutron_dx_mom_total", "Total dX vs rigidity; rigidity, GeV/c/Q; dX, cm", 800, -5.0, 10.0, 800, -5.0, 5.0);
352 fHistoManager->Create2 <TH2F>("TOF700_deutron_dy_mom_total", "Total dY vs rigidity; rigidity, GeV/c/Q; dY, cm", 800, -5.0, 10.0, 800, -5.0, 5.0);
353
354 fHistoManager->Create2 <TH2F>("TOF700_piMinus_dx_mom_total", "Total dX vs rigidity; rigidity, GeV/c/Q; dX, cm", 800, -5.0, 10.0, 800, -5.0, 5.0);
355 fHistoManager->Create2 <TH2F>("TOF700_piMinus_dy_mom_total", "Total dY vs rigidity; rigidity, GeV/c/Q; dY, cm", 800, -5.0, 10.0, 800, -5.0, 5.0);
356
357 fHistoManager->Create2 <TH2F>("TOF700_piPlus_dx_mom_total", "Total dX vs rigidity; rigidity, GeV/c/Q; dX, cm", 800, -5.0, 10.0, 800, -5.0, 5.0);
358 fHistoManager->Create2 <TH2F>("TOF700_piPlus_dy_mom_total", "Total dY vs rigidity; rigidity, GeV/c/Q; dY, cm", 800, -5.0, 10.0, 800, -5.0, 5.0);
359
360 // fHistoManager->Create1 <TH1F>("TOF700_proton_ResidX_total", "Total residuals in X coordinate; dX, cm; N", 100, -5.0, 5.0);
361 // fHistoManager->Create1 <TH1F>("TOF700_proton_ResidY_total", "Total residuals in Y coordinate; dY, cm; N", 100, -5.0, 5.0);
362 // fHistoManager->Create2 <TH2F>("TOF700_proton_Banan_total", "Total #beta vs rigidity; rigidity, GeV/c/Q; #beta", 400, -5.0, 10.0, 400, 0.2, 1.1);
363 // fHistoManager->Create2 <TH2F>("TOF700_proton_Mass2_total", "Total m^{2} vs rigidity; rigidity, GeV/c/Q; m^{2}", 400, -5.0, 10.0, 400, -0.5, 10.0);
364 // fHistoManager->Create2 <TH2F>("TOF700_proton_dx_tx_total", "Total dX vs Tx; t_{x}; dx, cm", 400, -0.5, 0.5, 400, -5.0, 5.0);
365 // fHistoManager->Create2 <TH2F>("TOF700_proton_Mass2_modules", "m^{2} vs module; module number; m^{2}", 60, 0, 60, 400, -0.5, 1.5);
366
367 // fHistoManager->Create1 <TH1F>("TOF700_deutron_ResidX_total", "Total residuals in X coordinate; dX, cm; N", 100, -5.0, 5.0);
368 // fHistoManager->Create1 <TH1F>("TOF700_deutron_ResidY_total", "Total residuals in Y coordinate; dY, cm; N", 100, -5.0, 5.0);
369 // fHistoManager->Create2 <TH2F>("TOF700_deutron_Banan_total", "Total #beta vs rigidity; rigidity, GeV/c/Q; #beta", 400, -5.0, 10.0, 400, 0.2, 1.1);
370 // fHistoManager->Create2 <TH2F>("TOF700_deutron_Mass2_total", "Total m^{2} vs rigidity; rigidity, GeV/c/Q; m^{2}", 400, -5.0, 10.0, 400, -0.5, 10.0);
371 // fHistoManager->Create2 <TH2F>("TOF700_deutron_dx_tx_total", "Total dX vs Tx; t_{x}; dx, cm", 400, -0.5, 0.5, 400, -5.0, 5.0);
372 // fHistoManager->Create2 <TH2F>("TOF700_deutron_Mass2_modules", "m^{2} vs module; module number; m^{2}", 60, 0, 60, 400, -0.5, 1.5);
373
374 // fHistoManager->Create1 <TH1F>("TOF700_piMinus_ResidX_total", "Total residuals in X coordinate; dX, cm; N", 100, -5.0, 5.0);
375 // fHistoManager->Create1 <TH1F>("TOF700_piMinus_ResidY_total", "Total residuals in Y coordinate; dY, cm; N", 100, -5.0, 5.0);
376 // fHistoManager->Create2 <TH2F>("TOF700_piMinus_Banan_total", "Total #beta vs rigidity; rigidity, GeV/c/Q; #beta", 400, -5.0, 10.0, 400, 0.2, 1.1);
377 // fHistoManager->Create2 <TH2F>("TOF700_piMinus_Mass2_total", "Total m^{2} vs rigidity; rigidity, GeV/c/Q; m^{2}", 400, -5.0, 10.0, 400, -0.5, 10.0);
378 // fHistoManager->Create2 <TH2F>("TOF700_piMinus_dx_tx_total", "Total dX vs Tx; t_{x}; dx, cm", 400, -0.5, 0.5, 400, -5.0, 5.0);
379 // fHistoManager->Create2 <TH2F>("TOF700_piMinus_Mass2_modules", "m^{2} vs module; module number; m^{2}", 60, 0, 60, 400, -0.5, 1.5);
380
381 // fHistoManager->Create1 <TH1F>("TOF700_piPlus_ResidX_total", "Total residuals in X coordinate; dX, cm; N", 100, -5.0, 5.0);
382 // fHistoManager->Create1 <TH1F>("TOF700_piPlus_ResidY_total", "Total residuals in Y coordinate; dY, cm; N", 100, -5.0, 5.0);
383 // fHistoManager->Create2 <TH2F>("TOF700_piPlus_Banan_total", "Total #beta vs rigidity; rigidity, GeV/c/Q; #beta", 400, -5.0, 10.0, 400, 0.2, 1.1);
384 // fHistoManager->Create2 <TH2F>("TOF700_piPlus_Mass2_total", "Total m^{2} vs rigidity; rigidity, GeV/c/Q; m^{2}", 400, -5.0, 10.0, 400, -0.5, 10.0);
385 // fHistoManager->Create2 <TH2F>("TOF700_piPlus_dx_tx_total", "Total dX vs Tx; t_{x}; dx, cm", 400, -0.5, 0.5, 400, -5.0, 5.0);
386 // fHistoManager->Create2 <TH2F>("TOF700_piPlus_Mass2_modules", "m^{2} vs module; module number; m^{2}", 60, 0, 60, 400, -0.5, 1.5);
387
388 for (Int_t i = 0; i < nModules; ++i) {
389
390 // fHistoManager->Create1 <TH1F>(Form("TOF700_ResidX_%d", i), Form("Residuals in X coordinate for module %d; dX, cm; N", i), 100, -5.0, 5.0);
391 // fHistoManager->Create1 <TH1F>(Form("TOF700_ResidY_%d", i), Form("Residuals in Y coordinate for module %d; dY, cm; N", i), 100, -5.0, 5.0);
392 // fHistoManager->Create2 <TH2F>(Form("TOF700_dx_tx_%d", i), Form("dX vs Tx for module %d; t_{x}; dx, cm", i), 400, -0.5, 0.5, 400, -5.0, 5.0);
393
394 // fHistoManager->Create2 <TH2F>(Form("TOF700_proton_Mass2_strips_%d", i), Form("m^{2} vs strips for module %d; strip number; m^{2}", i), 32, 0, 32, 100, -0.5, 1.5);
395 fHistoManager->Create1 <TH1F>(Form("TOF700_proton_ResidX_%d", i), Form("Residuals in X coordinate for module %d; dX, cm; N", i), 100, -5.0, 5.0);
396 fHistoManager->Create1 <TH1F>(Form("TOF700_proton_ResidY_%d", i), Form("Residuals in Y coordinate for module %d; dY, cm; N", i), 100, -5.0, 5.0);
397 fHistoManager->Create2 <TH2F>(Form("TOF700_proton_dx_tx_%d", i), Form("dX vs Tx for module %d; t_{x}; dx, cm", i), 400, -0.5, 0.5, 400, -5.0, 5.0);
398 //fHistoManager->Create2 <TH2F>(Form("TOF700_proton_Banan_%d", i), Form("#beta vs rigidity for module %d; rigidity, GeV/c/Q; #beta", i), 400, -5.0, 10.0, 400, 0.2, 1.1);
399
400 // fHistoManager->Create2 <TH2F>(Form("TOF700_deutron_Mass2_strips_%d", i), Form("m^{2} vs strips for module %d; strip number; m^{2}", i), 32, 0, 32, 100, -0.5, 1.5);
401 fHistoManager->Create1 <TH1F>(Form("TOF700_deutron_ResidX_%d", i), Form("Residuals in X coordinate for module %d; dX, cm; N", i), 100, -5.0, 5.0);
402 fHistoManager->Create1 <TH1F>(Form("TOF700_deutron_ResidY_%d", i), Form("Residuals in Y coordinate for module %d; dY, cm; N", i), 100, -5.0, 5.0);
403 fHistoManager->Create2 <TH2F>(Form("TOF700_deutron_dx_tx_%d", i), Form("dX vs Tx for module %d; t_{x}; dx, cm", i), 400, -0.5, 0.5, 400, -5.0, 5.0);
404 //fHistoManager->Create2 <TH2F>(Form("TOF700_deutron_Banan_%d", i), Form("#beta vs rigidity for module %d; rigidity, GeV/c/Q; #beta", i), 400, -5.0, 10.0, 400, 0.2, 1.1);
405
406 // fHistoManager->Create2 <TH2F>(Form("TOF700_piMinus_Mass2_strips_%d", i), Form("m^{2} vs strips for module %d; strip number; m^{2}", i), 32, 0, 32, 100, -0.5, 1.5);
407 fHistoManager->Create1 <TH1F>(Form("TOF700_piMinus_ResidX_%d", i), Form("Residuals in X coordinate for module %d; dX, cm; N", i), 100, -5.0, 5.0);
408 fHistoManager->Create1 <TH1F>(Form("TOF700_piMinus_ResidY_%d", i), Form("Residuals in Y coordinate for module %d; dY, cm; N", i), 100, -5.0, 5.0);
409 fHistoManager->Create2 <TH2F>(Form("TOF700_piMinus_dx_tx_%d", i), Form("dX vs Tx for module %d; t_{x}; dx, cm", i), 400, -0.5, 0.5, 400, -5.0, 5.0);
410 //fHistoManager->Create2 <TH2F>(Form("TOF700_piMinus_Banan_%d", i), Form("#beta vs rigidity for module %d; rigidity, GeV/c/Q; #beta", i), 400, -5.0, 10.0, 400, 0.2, 1.1);
411
412 // fHistoManager->Create2 <TH2F>(Form("TOF700_piPlus_Mass2_strips_%d", i), Form("m^{2} vs strips for module %d; strip number; m^{2}", i), 32, 0, 32, 100, -0.5, 1.5);
413 fHistoManager->Create1 <TH1F>(Form("TOF700_piPlus_ResidX_%d", i), Form("Residuals in X coordinate for module %d; dX, cm; N", i), 100, -5.0, 5.0);
414 fHistoManager->Create1 <TH1F>(Form("TOF700_piPlus_ResidY_%d", i), Form("Residuals in Y coordinate for module %d; dY, cm; N", i), 100, -5.0, 5.0);
415 fHistoManager->Create2 <TH2F>(Form("TOF700_piPlus_dx_tx_%d", i), Form("dX vs Tx for module %d; t_{x}; dx, cm", i), 400, -0.5, 0.5, 400, -5.0, 5.0);
416 //fHistoManager->Create2 <TH2F>(Form("TOF700_piPlus_Banan_%d", i), Form("#beta vs rigidity for module %d; rigidity, GeV/c/Q; #beta", i), 400, -5.0, 10.0, 400, 0.2, 1.1);
417 //fHistoManager->Create1 <TH1F>(Form("TOF700_Mass2_%d", i), Form("m^{2} for module %d; rigidity, GeV/c/Q; m^{2}", i), 400, -0.5, 10.0);
418 }
419}
420
421Float_t BmnTof700Qa::GetDxMom(FairTrackParam* par) {
422 // DX PROTONS
423 // 3.70/x^0.28-2.10
424 // DX PI-MINUS
425 //-7.61/(-x)^0.15+5.85
426
427 Float_t mom = 1.0 / par->GetQp();
428
429 return (mom > 0) ? 3.70 / TMath::Power(mom, 0.28) - 2.10 : -7.61 / TMath::Power(-1.0 * mom, 0.15) + 5.85;
430}
431
432Float_t BmnTof700Qa::GetDyMom(FairTrackParam* par) {
433 // DY PROTONS
434 // -1.56/x^0.94+0.18
435 // DY PI-MINUS
436 // +10.81/(-x)^0.05-9.53
437
438 Float_t mom = 1.0 / par->GetQp();
439
440 return (mom > 0) ? -1.56 / TMath::Power(mom, 0.94) + 0.18 : 10.81 / TMath::Power(-1.0 * mom, 0.05) - 9.53;
441}
442
443Float_t BmnTof700Qa::GetSigxMom(FairTrackParam* par) {
444 // SIGMA DX DEUTERONS
445 // 0.66/x^1.21+1.00
446 // SIGMA DX PI-MINUS
447 //-1.01/x^1.00+0.86
448
449 Float_t mom = 1.0 / par->GetQp();
450
451 return (mom > 0) ? 0.66 / TMath::Power(mom, 1.21) + 1.00 : -1.01 / TMath::Power(mom, 1.00) + 0.86;
452}
453
454Float_t BmnTof700Qa::GetSigyMom(FairTrackParam* par) {
455 // SIGMA DY DEUTERONS
456 // 0.98/x^1.20+0.57
457 // SIGMA DY PI-MINUS
458 // -0.87/x^1.00+0.58
459
460 Float_t mom = 1.0 / par->GetQp();
461
462 return (mom > 0) ? 0.98 / TMath::Power(mom, 1.20) + 0.57 : -0.87 / TMath::Power(mom, 1.00) + 0.58;
463}
int i
Definition P4_F32vec4.h:22
__m128 m
Definition P4_F32vec4.h:27
@ kBMNERROR
Definition BmnEnums.h:26
TH2 * H2(const TString &name) const
Return pointer to TH2 histogram.
void WriteToFile()
Write all histograms to current opened file.
void Create2(const TString &name, const TString &title, Int_t nofBinsX, Double_t minBinX, Double_t maxBinX, Int_t nofBinsY, Double_t minBinY, Double_t maxBinY)
Helper function for creation of 2-dimensional histograms and profiles. Template argument is a real ob...
void Create1(const TString &name, const TString &title, Int_t nofBins, Double_t minBin, Double_t maxBin)
Helper function for creation of 1-dimensional histograms and profiles. Template argument is a real ob...
TH1 * H1(const TString &name) const
Return pointer to TH1 histogram.
BmnStatus TGeoTrackPropagate(FairTrackParam *par, Double_t zOut, Int_t pdg, vector< Double_t > *F, Double_t *length, Bool_t isField)
Bool_t PropagateToZ(FairTrackParam *par, Double_t zDst, Int_t pdg, Int_t kalmanType, Double_t *len)
virtual InitStatus Init()
virtual void Exec(Option_t *opt)
virtual void Finish()
Int_t Extrapolate(Double_t z, Double_t *QP0=0, Bool_t line=false)
Access to i-th hit.
void GetTrackParam(FairTrackParam &track)
void SetPID(Int_t pidHypo)
static CbmStsDigiScheme * Instance(int version=1)
FairTrackParam * GetParamLast()
Definition CbmStsTrack.h:70
Int_t GetNStsHits() const
Definition CbmStsTrack.h:60
Int_t GetStsHitIndex(Int_t iHit) const
Definition CbmStsTrack.h:62
FairTrackParam * GetParamFirst()
Definition CbmStsTrack.h:69
Double_t GetZ() const
Definition CbmVertex.h:60
Double_t GetX() const
Definition CbmVertex.h:58
Int_t GetNTracks() const
Definition CbmVertex.h:63
Double_t GetY() const
Definition CbmVertex.h:59