BmnRoot
Loading...
Searching...
No Matches
BmnPairFinder.cxx
Go to the documentation of this file.
1
2#include <TCanvas.h>
3#include "TFile.h"
4#include <TGeoManager.h>
5#include "BmnPairFinder.h"
6#include "CbmMCTrack.h"
7#include "CbmStsTrack.h"
8#include "CbmStsPoint.h"
9#include "BmnVertex.h"
10#include "CbmVertex.h"
11
12const Float_t pathCutMin = 1.0;
13const Float_t pathCutMax = 100.0;
14const Float_t DCA12CutMin = 0.0;
15const Float_t DCA12CutMax = 2.0;
16const Float_t DCA1CutMin = 0.0;
17const Float_t DCA1CutMax = 10.0;
18const Float_t DCA2CutMin = 0.0;
19const Float_t DCA2CutMax = 10.0;
20const Float_t DCA0CutMin = 0.0;
21const Float_t DCA0CutMax = 1.0;
22const Float_t PVXCutMin = -1.0;
23const Float_t PVXCutMax = 1.5;
24const Float_t PVYCutMin = -1.0;
25const Float_t PVYCutMax = 1.0;
26const Float_t PVZCutMin = -0.4;
27const Float_t PVZCutMax = 0.4;
28const Int_t NTrCutMin = 0;
29const Int_t NTrCutMax = 100;
30
31BmnPairFinder::BmnPairFinder(TParticlePDG * particle)
32{
33 cout << "\nBmnPairFinder::Constructor()" << endl;
34 fDstHeader = nullptr;
35 fKalman = nullptr;
36 fParticlePair = nullptr;
37 fVertex = nullptr;
38 evNumber = 0;
39 fRunId = -10;
40 fParticle = unique_ptr<TParticlePDG>(particle);
41}
42
44{
45 cout << "\nBmnPairFinder::Destructor()" << endl;
46}
47
49{
50 cout << "\nBmnPairFinder::Init()" << endl;
51 switch (fParticle->PdgCode()){
52 case 310:
53 fDecayFirstParticle = unique_ptr<TParticlePDG>(TDatabasePDG::Instance()->GetParticle(211));
54 fDecaySecondParticle = unique_ptr<TParticlePDG>(TDatabasePDG::Instance()->GetParticle(-211));
55 lowV0Z = -10.0;
56 break;
57 case 3122:
58 fDecayFirstParticle = unique_ptr<TParticlePDG>(TDatabasePDG::Instance()->GetParticle(2212));
59 fDecaySecondParticle = unique_ptr<TParticlePDG>(TDatabasePDG::Instance()->GetParticle(-211));
60 lowV0Z = 0.0;
61 break;
62 default:
63 LOGF(error, "Unsupported particle %d %s !", fParticle->PdgCode(), fParticle->GetName());
64 return kERROR;
65 break;
66 }
67 if (fParticle)
68 LOGF(info, "Search decays: %s -> %s %s",
69 fParticle->GetName(), fDecayFirstParticle->GetName(), fDecaySecondParticle->GetName());
70
71 FairRootManager *ioman = FairRootManager::Instance();
72 fDstHeader = (DstEventHeader *)ioman->GetObject("DstEventHeader.");
73 fStsTracks = (TClonesArray *)ioman->GetObject("StsVector");
74 fVertex = (CbmVertex *)ioman->GetObject("PrimaryVertex.");
75
76 fParticlePair = new TClonesArray("BmnParticlePair");
77 ioman->Register("ParticlePair", "Lambda", fParticlePair, kTRUE);
78
79 fKalman = new BmnKalmanFilter();
80
81 return kSUCCESS;
82}
83
84// ------------------------------------------------------------------
85
86void BmnPairFinder::Exec(Option_t *option)
87{
88 // cout << "\nBmnPairFinder::Exec()" << endl;
89 // cout << "Event = " << evNumber++ << endl;
90
91 evNumber++;
92 if (evNumber % 1000 == 0) printf("Event number: %d\n", evNumber);
93
94 fParticlePair->Delete();
95
96 if (fVertex->GetNTracks() < 2)
97 return;
98
99
100 if (fVertex->GetZ() > PVZCutMax || fVertex->GetZ() < PVZCutMin) return;
101 if (fVertex->GetX() > PVXCutMax || fVertex->GetX() < PVXCutMin) return;
102 if (fVertex->GetY() > PVYCutMax || fVertex->GetY() < PVYCutMin) return;
103
104 if (fStsTracks->GetEntriesFast() > NTrCutMax || fStsTracks->GetEntriesFast() < NTrCutMin) return;
105
106 for (Int_t iTrack = 0; iTrack < fStsTracks->GetEntriesFast(); iTrack++)
107 {
108 CbmStsTrack *cbmststrack1 = (CbmStsTrack *)fStsTracks->At(iTrack);
109
110 BmnGlobalTrack globTr1;
111
112 //globTr1.SetGemTrackIndex(i);
113 globTr1.SetParamFirst(*(cbmststrack1->GetParamFirst()));
114 globTr1.SetParamLast(*(cbmststrack1->GetParamLast()));
115 Int_t nHitsPos = cbmststrack1->GetNStsHits();
116 //if (nHitsPos < 5) continue;
117 globTr1.SetNHits(nHitsPos);
118 globTr1.SetNDF(cbmststrack1->GetNDF());
119 globTr1.SetChi2(cbmststrack1->GetChi2());
120
121 FairTrackParam par_prot_Vp = *(globTr1.GetParamFirst());
122
123 if (par_prot_Vp.GetQp() < 0) continue;
124 if (fKalman->TGeoTrackPropagate(&par_prot_Vp, fVertex->GetZ(), fDecayFirstParticle->PdgCode(), nullptr, nullptr, kTRUE) == kBMNERROR) continue;
125
126 for (Int_t jTrack = 0; jTrack < fStsTracks->GetEntriesFast(); jTrack++)
127 {
128 CbmStsTrack *cbmststrack2 = (CbmStsTrack *)fStsTracks->At(jTrack);
129
130 BmnGlobalTrack globTr2;
131
132 globTr2.SetParamFirst(*(cbmststrack2->GetParamFirst()));
133 globTr2.SetParamLast(*(cbmststrack2->GetParamLast()));
134 Int_t nHitsNeg = cbmststrack2->GetNStsHits();
135 //if (nHitsNeg < 4) continue;
136 globTr2.SetNHits(cbmststrack2->GetNStsHits());
137 globTr2.SetNDF(cbmststrack2->GetNDF());
138 globTr2.SetChi2(cbmststrack2->GetChi2());
139
140 FairTrackParam par_pion_Vp = *(globTr2.GetParamFirst());
141 if (par_pion_Vp.GetQp() > 0) continue; //suppose negative particle is pi-
142 if (fKalman->TGeoTrackPropagate(&par_pion_Vp, fVertex->GetZ(), fDecaySecondParticle->PdgCode(), nullptr, nullptr, kTRUE) == kBMNERROR) continue;
143 //Float_t vDist = TMath::Sqrt(Sq(par_pion_Vp.GetX() - fVertex->GetX()) + Sq(par_pion_Vp.GetY() - fVertex->GetY()));
144 //if (vDist > 2.0) continue;
145
147
148 Double_t V0Z = FindV0ByVirtualPlanes(&globTr1, &globTr2, 20, 50);
149 if (V0Z < lowV0Z || V0Z > 50) continue;
150
151 FairTrackParam par_prot_V0 = *(globTr1.GetParamFirst());
152 FairTrackParam par_pion_V0 = *(globTr2.GetParamFirst());
153 if (fKalman->TGeoTrackPropagate(&par_prot_V0, V0Z, fDecayFirstParticle->PdgCode(), nullptr, nullptr, kTRUE) == kBMNERROR) continue;
154 if (fKalman->TGeoTrackPropagate(&par_pion_V0, V0Z, fDecaySecondParticle->PdgCode(), nullptr, nullptr, kTRUE) == kBMNERROR) continue;
155
156 if (fKalman->TGeoTrackPropagate(&par_prot_Vp, fVertex->GetZ(), fDecayFirstParticle->PdgCode(), nullptr, nullptr, kTRUE) == kBMNERROR) continue;
157 if (fKalman->TGeoTrackPropagate(&par_pion_Vp, fVertex->GetZ(), fDecaySecondParticle->PdgCode(), nullptr, nullptr, kTRUE) == kBMNERROR) continue;
158
159 Double_t V0X = .5 * (par_prot_V0.GetX() + par_pion_V0.GetX());
160 if (V0X < -2.0 || V0X > 2.0) continue;
161 Double_t V0Y = .5 * (par_prot_V0.GetY() + par_pion_V0.GetY());
162 if (V0Y < -2.0 || V0Y > 2.0) continue;
163
164 Float_t Path = GetPath(fVertex->GetX(), fVertex->GetY(), fVertex->GetZ(), V0X, V0Y, V0Z);
165 if (Path < pathCutMin || Path > pathCutMax) continue;
166
167 Float_t DCA0 = FindDCA0(par_prot_V0, par_pion_V0, fVertex->GetX(), fVertex->GetY(), fVertex->GetZ(), V0X, V0Y, V0Z);
168 if (DCA0 < DCA0CutMin || DCA0 > DCA0CutMax) continue;
169
170 Float_t DCA1 = TMath::Sqrt(Sq(par_prot_Vp.GetX() - fVertex->GetX()) + Sq(par_prot_Vp.GetY() - fVertex->GetY()));
171 if (DCA1 < DCA1CutMin || DCA1 > DCA1CutMax) continue;
172
173 Float_t DCA2 = TMath::Sqrt(Sq(par_pion_Vp.GetX() - fVertex->GetX()) + Sq(par_pion_Vp.GetY() - fVertex->GetY()));
174 if (DCA2 < DCA2CutMin || DCA2 > DCA2CutMax) continue;
175
176 Float_t DCA12 = TMath::Sqrt(Sq(par_prot_V0.GetX() - par_pion_V0.GetX()) + Sq(par_prot_V0.GetY() - par_pion_V0.GetY()));
177 if (DCA12 < DCA12CutMin || DCA12 > DCA12CutMax) continue;
178
179 pair.SetPath(Path);
180 pair.SetDCA0(DCA0);
181 pair.SetDCA12(DCA12);
182 pair.SetDCA1(DCA1);
183 pair.SetDCA2(DCA2);
184
185 pair.SetV0X(V0X);
186 pair.SetV0Y(V0Y);
187 pair.SetV0Z(V0Z);
188
189 pair.SetTxPair(par_prot_V0.GetTx(), par_pion_V0.GetTx());
190 pair.SetTyPair(par_prot_V0.GetTy(), par_pion_V0.GetTy());
191
192
193 pair.SetMomPair(globTr1.GetP(), globTr2.GetP());
194 pair.SetNHitsPair(nHitsPos, nHitsNeg);
195 pair.SetChi2Pair(globTr1.GetChi2(), globTr2.GetChi2());
196 pair.SetNDFPair(globTr1.GetNDF(), globTr2.GetNDF());
197 pair.SetInvMass(FindInvMass(par_prot_V0, par_pion_V0));
198 pair.SetAngleDecayProducts(FindAngleDecayProducts(par_prot_V0, par_pion_V0));
199 pair.SetEtaPair(FindEta(par_prot_V0), FindEta(par_pion_V0));
200 new ((*fParticlePair)[fParticlePair->GetEntriesFast()]) BmnParticlePair(pair);
201 }
202 }
203}
204// -------------------------------------------------------------------
205
206Double_t BmnPairFinder::FindV0ByVirtualPlanes(BmnGlobalTrack *track1, BmnGlobalTrack *track2, Double_t z_0, Double_t range)
207
208{
209
210 const Int_t nPlanes = 5; // FIXME
211
212 // if (V0Z < Vpz - 10. || V0Z > fDetector->GetGemStation(5)->GetZPosition()) // FIXME!
213 // continue;
214
215 while (range >= 0.1)
216 {
217 Double_t zMax = z_0 + range;
218 Double_t zMin = z_0 - range;
219 Double_t zStep = (zMax - zMin) / nPlanes;
220
221 Double_t zPlane[nPlanes];
222 Double_t Dist[nPlanes];
223
224 FairTrackParam par1 = *(track1->GetParamFirst());
225 FairTrackParam par2 = *(track2->GetParamFirst());
226
227 for (Int_t iPlane = 0; iPlane < nPlanes; ++iPlane)
228 {
229 zPlane[iPlane] = zMax - iPlane * zStep;
230 BmnStatus stat1 = fKalman->TGeoTrackPropagate(&par1, zPlane[iPlane], fDecayFirstParticle->PdgCode(), nullptr, nullptr, kTRUE);
231 BmnStatus stat2 = fKalman->TGeoTrackPropagate(&par2, zPlane[iPlane], fDecaySecondParticle->PdgCode(), nullptr, nullptr, kTRUE);
232
233 if (stat1 == kBMNERROR || stat2 == kBMNERROR)
234 {
235 return DBL_MAX;
236 }
237
238 Dist[iPlane] = TMath::Sqrt(Sq(par1.GetX() - par2.GetX()) + Sq(par1.GetY() - par2.GetY()));
239 }
240
241 map<Double_t, Double_t> dist_zPlanes; // distance --> z-position
242
243 for (Int_t iPlane = 0; iPlane < nPlanes; iPlane++)
244 dist_zPlanes[Dist[iPlane]] = zPlane[iPlane];
245
246 z_0 = dist_zPlanes.begin()->second;
247
248 range /= 2;
249 }
250
251 return z_0;
252}
253
254Double_t BmnPairFinder::FindInvMass(FairTrackParam pi_plus_meson_V0, FairTrackParam pi_minus_meson_V0)
255{
256
257 Double_t Tx1, Ty1, Tx2, Ty2, p1, p2;
258 Double_t A1, A2;
259 TLorentzVector lPos, lNeg;
260
261 Tx1 = pi_plus_meson_V0.GetTx();
262 Ty1 = pi_plus_meson_V0.GetTy();
263 Tx2 = pi_minus_meson_V0.GetTx();
264 Ty2 = pi_minus_meson_V0.GetTy();
265 p1 = 1. / pi_plus_meson_V0.GetQp();
266 //cout<<"small_test_1:"<<""<<p1<<endl;
267 p2 = 1. / pi_minus_meson_V0.GetQp();
268 //cout<<"small_test_2:"<<""<<p2<<endl;
269 // armenPodol = ArmenterosPodol(proton_V0, pion_V0);
270
271 A1 = 1. / TMath::Sqrt(Tx1 * Tx1 + Ty1 * Ty1 + 1);
272 lPos.SetXYZM(Tx1 * A1 * p1, Ty1 * A1 * p1, p1 * A1, fDecayFirstParticle->Mass());
273
274 p2 *= -1.; // Since in the calculations pos. mom. values should be used
275
276 A2 = 1. / TMath::Sqrt(Tx2 * Tx2 + Ty2 * Ty2 + 1);
277 lNeg.SetXYZM(Tx2 * A2 * p2, Ty2 * A2 * p2, p2 * A2, fDecaySecondParticle->Mass());
278
279 // partPair.SetAlpha(armenPodol.X());
280 // partPair.SetPtPodol(armenPodol.Y());
281 return (TLorentzVector((lPos + lNeg)).Mag());
282}
283
284Double_t BmnPairFinder::GetPath(Double_t vPx, Double_t vPy, Double_t vPz, Double_t v0x, Double_t v0y, Double_t v0z)
285{
286 TVector3 V0(v0x, v0y, v0z);
287 TVector3 Vp(vPx, vPy, vPz);
288 return TVector3(V0 - Vp).Mag();
289}
290
291Double_t BmnPairFinder::FindDCA0(FairTrackParam pi_plus_meson_V0, FairTrackParam pi_minus_meson_V0, Double_t vPx, Double_t vPy, Double_t vPz, Double_t v0x, Double_t v0y, Double_t v0z)
292{
293 Double_t Tx1, Ty1, Tx2, Ty2, p1, p2;
294 Double_t A1, A2;
295 TVector3 V0(v0x, v0y, v0z);
296 TVector3 Vp(vPx, vPy, vPz);
297
298 Tx1 = pi_plus_meson_V0.GetTx();
299 Ty1 = pi_plus_meson_V0.GetTy();
300 Tx2 = pi_minus_meson_V0.GetTx();
301 Ty2 = pi_minus_meson_V0.GetTy();
302 p1 = 1. / pi_plus_meson_V0.GetQp();
303 p2 = 1. / pi_minus_meson_V0.GetQp();
304
305 A1 = 1. / TMath::Sqrt(Tx1 * Tx1 + Ty1 * Ty1 + 1);
306 p2 *= -1.;
307 A2 = 1. / TMath::Sqrt(Tx2 * Tx2 + Ty2 * Ty2 + 1);
308
309 Double_t PzPart1 = Abs(p1) * A1;
310 Double_t PzPart2 = Abs(p2) * A2;
311
312 Double_t PzPart0 = PzPart1 + PzPart2;
313 Double_t PxPart0 = PzPart1 * Tx1 + PzPart2 * Tx2;
314 Double_t PyPart0 = PzPart1 * Ty1 + PzPart2 * Ty2;
315
316 Double_t TxPart0 = PxPart0 / PzPart0;
317 Double_t TyPart0 = PyPart0 / PzPart0;
318
319 Double_t xPart0 = TxPart0 * (Vp.Z() - V0.Z()) + V0.X();
320 Double_t yPart0 = TyPart0 * (Vp.Z() - V0.Z()) + V0.Y();
321
322 Double_t dca0 = TMath::Sqrt(Sq(xPart0 - Vp.X()) + Sq(yPart0 - Vp.Y()));
323
324 return dca0;
325}
326
327Double_t BmnPairFinder::FindEta(FairTrackParam pi_plus_meson_V0)
328{
329 Double_t Tx = pi_plus_meson_V0.GetTx();
330 Double_t Ty = pi_plus_meson_V0.GetTy();
331 Double_t p = 1. / pi_plus_meson_V0.GetQp();
332
333 Double_t Pz = Abs(p) / TMath::Sqrt(1 + Tx * Tx + Ty * Ty);
334
335 Double_t eta = 0.5 * Log((Abs(p) + Pz) / (Abs(p) - Pz));
336
337 return eta;
338}
339
340Double_t BmnPairFinder::FindAngleDecayProducts(FairTrackParam proton_V0, FairTrackParam pi_minus_meson_V0)
341{
342 Double_t Tx1, Ty1, Tx2, Ty2, p1, p2;
343 Double_t A1, A2;
344 TVector3 Pproton;
345 TVector3 Pnegpimeson;
346
347 Tx1 = proton_V0.GetTx();
348 Ty1 = proton_V0.GetTy();
349 Tx2 = pi_minus_meson_V0.GetTx();
350 Ty2 = pi_minus_meson_V0.GetTy();
351
352 p1 = 1. / proton_V0.GetQp();
353 p2 = 1. / pi_minus_meson_V0.GetQp();
354
355 A1 = 1. / TMath::Sqrt(Tx1 * Tx1 + Ty1 * Ty1 + 1);
356 p2 *= -1.;
357 A2 = 1. / TMath::Sqrt(Tx2 * Tx2 + Ty2 * Ty2 + 1);
358
359 Pproton.SetX(Tx1 * A1 * p1);
360 Pproton.SetY(Ty1 * A1 * p1);
361 Pproton.SetZ(p1 * A1);
362
363 Pnegpimeson.SetX(Tx2 * A2 * p2);
364 Pnegpimeson.SetY(Ty2 * A2 * p2);
365 Pnegpimeson.SetZ(p2 * A2);
366
367 Double_t dot_product = Pproton.Dot(Pnegpimeson);
368 //cout<<"1:"<<dot_product<<endl;
369 Double_t mag_Pproton = Pproton.Mag();
370 //cout<<"2:"<<mag_Pproton<<endl;
371 Double_t mag_Pnegpimeson = Pnegpimeson.Mag();
372 //cout<<"3:"<<mag_Pnegpimeson<<endl;
373 //cout<<"4:"<<ACos(dot_product/(mag_Pproton*mag_Pnegpimeson))<<endl;
374 //cout<<"5:"<<ACos(dot_product/(mag_Pproton*mag_Pnegpimeson))*180/Pi()<<endl;
375 return ACos(dot_product/(mag_Pproton*mag_Pnegpimeson))*180/Pi();
376}
377
379{
380
381 cout << "\n-I- [BmnPairFinder::Finish] " << endl;
382}
Float_t Dist(Float_t x1, Float_t y1, Float_t x2, Float_t y2)
Definition BmnMath.cxx:432
const Float_t PVZCutMin
const Float_t DCA0CutMax
const Float_t pathCutMin
const Int_t NTrCutMin
const Float_t DCA1CutMin
const Float_t DCA12CutMax
const Float_t DCA12CutMin
const Float_t DCA1CutMax
const Float_t DCA0CutMin
const Float_t PVYCutMin
const Float_t DCA2CutMax
const Int_t NTrCutMax
const Float_t PVXCutMin
const Float_t PVZCutMax
const Float_t pathCutMax
const Float_t PVYCutMax
const Float_t PVXCutMax
const Float_t DCA2CutMin
BmnStatus
Definition BmnEnums.h:24
@ kBMNERROR
Definition BmnEnums.h:26
BmnStatus TGeoTrackPropagate(FairTrackParam *par, Double_t zOut, Int_t pdg, vector< Double_t > *F, Double_t *length, Bool_t isField)
virtual InitStatus Init()
BmnPairFinder(TParticlePDG *particle=TDatabasePDG::Instance() ->GetParticle("Lambda0"))
virtual void Exec(Option_t *option)
virtual void Finish()
virtual ~BmnPairFinder()
void SetDCA2(Double_t val)
void SetTxPair(Double_t val1, Double_t val2)
void SetNHitsPair(vector< Int_t > part1, vector< Int_t > part2)
void SetTyPair(Double_t val1, Double_t val2)
void SetInvMass(Double_t val)
void SetChi2Pair(Double_t val1, Double_t val2)
void SetV0Y(Double_t val)
void SetV0Z(Double_t val)
void SetDCA1(Double_t val)
void SetAngleDecayProducts(Double_t val)
void SetDCA12(Double_t val)
void SetMomPair(Double_t val1, Double_t val2)
void SetDCA0(Double_t val)
void SetEtaPair(Double_t val1, Double_t val2)
void SetPath(Double_t val)
void SetV0X(Double_t val)
void SetNDFPair(Int_t val1, Int_t val2)
Int_t GetNDF() const
Definition BmnTrack.h:60
void SetChi2(Double_t chi2)
Definition BmnTrack.h:97
void SetParamLast(FairTrackParam &par)
Definition BmnTrack.h:89
FairTrackParam * GetParamFirst()
Definition BmnTrack.h:72
Float_t GetChi2() const
Definition BmnTrack.h:56
Double_t GetP()
Definition BmnTrack.h:80
void SetNHits(Int_t n)
Definition BmnTrack.h:105
void SetParamFirst(FairTrackParam &par)
Definition BmnTrack.h:85
void SetNDF(Int_t ndf)
Definition BmnTrack.h:101
FairTrackParam * GetParamLast()
Definition CbmStsTrack.h:70
Int_t GetNStsHits() const
Definition CbmStsTrack.h:60
Double_t GetChi2() const
Definition CbmStsTrack.h:66
FairTrackParam * GetParamFirst()
Definition CbmStsTrack.h:69
Int_t GetNDF() const
Definition CbmStsTrack.h:67
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