BmnRoot
Loading...
Searching...
No Matches
BmnParticle.h
Go to the documentation of this file.
1// -------------------------------------------------------------------------
2// ----- BmnParticle header file -----
3// ----- MpdParticle created 21/01/2013 by A. Zinchenko -----
4// ----- BmnParticle created 02/07/2025 by A. Zinchenko -----
5// -------------------------------------------------------------------------
6
14#ifndef BMNPARTICLE_H
15#define BMNPARTICLE_H 1
16
17// #include "CbmKFTrack.h"
18
19// #include "TArrayI.h"
20#include "TMatrixD.h"
21// #include "TMatrixFSym.h"
22#include "TMath.h"
23#include "TVector3.h"
24
25#include <vector>
26// class MpdKalmanTrack;
27class CbmStsTrack;
28// class MpdVertex;
29class CbmVertex;
30// class MpdKalmanFilter;
31class CbmKF;
32class CbmKFTrack;
33
34using namespace std;
35
36class BmnParticle : public TObject
37{
38
39 public:
42 BmnParticle(const BmnParticle& part);
44 Int_t indx = -1,
45 Double_t mass = 0.1396,
46 Double_t* orig = NULL);
48 Int_t indx = -1,
49 Double_t mass = 0.1396,
50 Double_t* orig = NULL);
51
53
55 virtual ~BmnParticle();
56
58 void Print();
59
61 int GetNDF() const { return fNDF; }
62 Int_t GetIndx() const { return fIndx; }
63 Int_t GetPdg() const { return fPdg; }
64 Int_t GetCharge() const { return fCharge; }
65 Double_t GetMass() const { return fMass; }
66 Double_t GetMeas(Int_t i) const { return fMeas(i, 0); }
67 Double_t GetXY(Int_t i) const { return fXY0[i]; }
68 Int_t GetFlag() const { return fFlag; }
69 Double_t GetZ0() const { return fZ; }
70 // Double_t Phi() const { return GetMeas(2); }
71 // Double_t Pt() const { return TMath::Min (TMath::Abs(1./GetMeas(4)*fieldConst), 100.); }
72 // Double_t Theta() const { return GetMeas(3); }
73 // Double_t Momentum() const { return Pt() / TMath::Sin(Theta()); }
74 // TVector3 Momentum3() const { return TVector3(Pt()*TMath::Cos(Phi()), Pt()*TMath::Sin(Phi()),
75 // Momentum()*TMath::Cos(Theta())); }
76 Double_t Phi() const { return TMath::ATan2(fq(1, 0), fq(0, 0)); } // smoothed value
77 Double_t Pt() const { return TMath::Min(Momentum() * TMath::Sin(Theta()), 100.); }
78 Double_t GetPt() const { return Pt(); }
79 // Double_t Theta() const { return fq(1,0); }
80 Double_t Theta() const;
81 Double_t GetTheta() const { return Theta(); }
82 // Double_t Momentum() const { return fCharge == 0 ? fq(2, 0) : Pt() / TMath::Sin(Theta()); }
83 Double_t Momentum() const { return fCharge == 0 ? fq(2, 0) : TMath::Abs(fCharge / fq(2, 0)); }
84 Double_t GetMomentum() const { return Momentum(); }
85 TVector3 Momentum3() const
86 {
87 return TVector3(Pt() * TMath::Cos(Phi()), Pt() * TMath::Sin(Phi()), Momentum() * TMath::Cos(Theta()));
88 }
89 Double_t Energy() const;
90 Double_t Rapidity() const;
91 Double_t Dca() const { return fMeas(0, 0); }
92 Int_t Ndaughters() const { return fDaughtersInds.size(); }
93 const vector<Int_t>& DaughterInds() const { return fDaughtersInds; }
94 void Track2Part(CbmKFTrack& track, Bool_t setWeight, Double_t* orig); // conversion from track to particle
95 const Double_t Chi2Vertex() { return fChi2ver; }
96 Double_t Chi2Vertex(CbmVertex* vtx);
97 Double_t Chi2() const { return fChi2; }
98 Bool_t Point00() const { return fPoint00; }
99 void FillJ(); // fill Jacobian matrix fJ
100 void FillJinv(TVector3& mom3); // fill Jacobian matrix fJinv
101
102 TMatrixD& GetMeas() { return fMeas; }
103 TMatrixD& GetJ() { return fJ; }
104 TMatrixD& GetJinv() { return fJinv; }
105 TMatrixD& GetD() { return fD; }
106 TMatrixD& GetE() { return fE; }
107 TMatrixD& GetA() { return fA; }
108 TMatrixD& GetB() { return fB; }
109 TMatrixD& GetC() { return fC; }
110 TMatrixD& GetG() { return fG; }
111 TMatrixD& GetW() { return fW; }
112 TMatrixD& Getq() { return fq; }
113 TMatrixD& Getx() { return fx; }
114 // TMatrixD& GetW() const { return fW; }
115 Double_t GetX() const { return fx(0, 0); }
116 Double_t GetY() const { return fx(1, 0); }
117 Double_t GetZ() const { return fx(2, 0); }
118 Double_t GetPx() const { return Pt() * TMath::Cos(Phi()); }
119 Double_t GetPy() const { return Pt() * TMath::Sin(Phi()); }
120 Double_t GetPz() const { return Momentum() * TMath::Cos(Theta()); }
121
122 void SetIndx(Int_t indx) { fIndx = indx; }
123 void SetPdg(Int_t pdg)
124 {
125 fPdg = pdg;
126 SetMass();
127 }
128 void SetCharge(Int_t charge) { fCharge = charge; }
129 void SetMass(Double_t mass = -2.0);
130 void AddDaughter(Int_t indx) { fDaughtersInds.push_back(indx); }
131 Double_t BuildMother(vector<BmnParticle*>& vDaught);
132 Double_t BuildMother(vector<CbmStsTrack*>& vTracks, vector<BmnParticle*>& vDaught);
133 void SetChi2(Double_t chi2) { fChi2 = chi2; }
134 void SetMeas(TMatrixD& matr) { fMeas = matr; }
135 void SetCovD(TMatrixD& matr) { fD = matr; }
136 void SetCovE(TMatrixD& matr) { fE = matr; }
137 // void SetCovQ(TMatrixD &matr) { fQ = matr; }
138 void SetA(TMatrixD& matr) { fA = matr; }
139 void SetB(TMatrixD& matr) { fB = matr; }
140 void SetC(TMatrixD& matr) { fC = matr; }
141 void SetG(TMatrixD& matr) { fG = matr; }
142 void SetW(TMatrixD& matr) { fW = matr; }
143 void Setq(TMatrixD& matr) { fq = matr; }
144 void Setx(TMatrixD& matr) { fx = matr; }
145 void SetXY(Double_t x, Double_t y)
146 {
147 fXY0[0] = x;
148 fXY0[1] = y;
149 }
150 void SetFlag(Int_t flag) { fFlag = flag; }
151 void SetNDF(int ndf) { fNDF = ndf; }
152 void SetZ(Double_t z) { fZ = z; }
155
156 private:
157 // FindPca(CbmKFTrackInterface& tr, Double_t* vert);
158 void FindPca(CbmKFTrack& tr, Double_t* vert);
159 TVector3 Parabola(TVector3& pos0, TVector3& pos1, TVector3& pos2);
160 // void WeightAtDca(CbmStsTrack &track, Double_t *vert); // obtain BmnParticle weight at DCA
161 void WeightAtDca(CbmKFTrack& track, Double_t* vert); // obtain BmnParticle weight at DCA
162
163 Int_t fIndx; // index of particle
164 Int_t fPdg; // PDG hypothesis
165 Int_t fCharge; // charge
166 Double_t fMass; // particle mass (GeV)
167 Double_t fieldConst;
168 vector<Int_t> fDaughtersInds; // indices of particles it is created from
169 TMatrixD fMeas; // vector of measurements (params)
170 TMatrixD fq; // geometrical momentum
171 TMatrixD fx; // particle origin (production vertex)
172 Double_t fXY0[2]; // X and Y at DCA
173 TMatrixD fJ;
174 TMatrixD fJinv;
175 TMatrixD fD;
176 TMatrixD fE;
177 // TMatrixD fQ; //! covariance cov(qk,qj)
178 TMatrixD fA;
179 TMatrixD fB;
180 TMatrixD fC;
181 TMatrixD fG; // covariance of params
182 TMatrixD fW;
183 Double_t fChi2; // Chi2 of mother particle
184 Double_t fChi2ver; // Chi2 of particle w.r.t. vertex
185 Int_t fFlag; // status flag
186 Bool_t fPoint00; // flag for tracks extrapolated to (0,0)
187 Double_t fZ; // Z-position
188 int fNDF;
189
190 // MpdKalmanFilter *fKF;
191 CbmKF* fKF;
192
193 ClassDef(BmnParticle, 1);
194};
195
196//__________________________________________________________________________
197
198inline Double_t BmnParticle::Theta() const
199{
200 // Theta angle
201
202 return TMath::ATan(TMath::Sqrt(fq(1, 0) * fq(1, 0) + fq(0, 0) * fq(0, 0)));
203}
204
205//__________________________________________________________________________
206
207#endif
int i
Definition P4_F32vec4.h:22
Double_t GetMass() const
Definition BmnParticle.h:65
Int_t GetFlag() const
Definition BmnParticle.h:68
TMatrixD & GetW()
Double_t GetTheta() const
Definition BmnParticle.h:81
TMatrixD & Getx()
Double_t Phi() const
Definition BmnParticle.h:76
TMatrixD & GetJinv()
Double_t Energy() const
void SetIndx(Int_t indx)
void SetCovD(TMatrixD &matr)
Double_t GetMeas(Int_t i) const
Definition BmnParticle.h:66
void SetG(TMatrixD &matr)
void SetPdg(Int_t pdg)
Double_t BuildMother(vector< BmnParticle * > &vDaught)
void SetXY(Double_t x, Double_t y)
void SetW(TMatrixD &matr)
Double_t GetPz() const
Double_t GetX() const
CbmKFTrack GetKFTrack() const
Double_t Rapidity() const
Bool_t Point00() const
flag for tracks extrapolated to (0,0)
Definition BmnParticle.h:98
Double_t Theta() const
Double_t GetZ() const
BmnParticle()
Default ctor.
const Double_t Chi2Vertex()
return Chi2 w.r.t. vertex
Definition BmnParticle.h:95
BmnParticle(CbmStsTrack &track, Int_t indx=-1, Double_t mass=0.1396, Double_t *orig=NULL)
ctor from STS track
TMatrixD & GetD()
Double_t GetMomentum() const
Definition BmnParticle.h:84
Double_t Pt() const
Definition BmnParticle.h:77
BmnParticle(CbmKFTrack &track, Int_t indx=-1, Double_t mass=0.1396, Double_t *orig=NULL)
ctor from Kalman track
Double_t BuildMother(vector< CbmStsTrack * > &vTracks, vector< BmnParticle * > &vDaught)
void SetCharge(Int_t charge)
Double_t Chi2Vertex(CbmVertex *vtx)
compute Chi2 w.r.t. vertex
Int_t GetCharge() const
Definition BmnParticle.h:64
virtual ~BmnParticle()
void Track2Part(CbmKFTrack &track, Bool_t setWeight, Double_t *orig)
TMatrixD & GetA()
void FillJ()
void SetChi2(Double_t chi2)
void SetA(TMatrixD &matr)
Double_t GetZ0() const
Definition BmnParticle.h:69
TMatrixD & GetMeas()
Int_t GetPdg() const
Definition BmnParticle.h:63
void SetFlag(Int_t flag)
int GetNDF() const
Definition BmnParticle.h:61
TVector3 Momentum3() const
Definition BmnParticle.h:85
void SetZ(Double_t z)
void SetMass(Double_t mass=-2.0)
void SetCovE(TMatrixD &matr)
Double_t Dca() const
signed DCA
Definition BmnParticle.h:91
void SetB(TMatrixD &matr)
Double_t GetY() const
TMatrixD & Getq()
Double_t Chi2() const
Chi2 of mother particle.
Definition BmnParticle.h:97
TMatrixD & GetJ()
Double_t GetPy() const
void SetMeas(TMatrixD &matr)
const vector< Int_t > & DaughterInds() const
Definition BmnParticle.h:93
void SetNDF(int ndf)
void Setq(TMatrixD &matr)
void FillJinv(TVector3 &mom3)
void AddDaughter(Int_t indx)
TMatrixD & GetE()
Double_t GetPt() const
Definition BmnParticle.h:78
Double_t GetPx() const
TMatrixD & GetG()
Int_t Ndaughters() const
Definition BmnParticle.h:92
Double_t GetXY(Int_t i) const
Definition BmnParticle.h:67
void Setx(TMatrixD &matr)
TMatrixD & GetC()
void ParamsAtDca()
void Print()
BmnParticle(const BmnParticle &part)
copy constructor
Int_t GetIndx() const
Definition BmnParticle.h:62
TMatrixD & GetB()
void SetC(TMatrixD &matr)
Double_t Momentum() const
Definition BmnParticle.h:83
BmnParticle & operator=(const BmnParticle &part)
assignment operator
Definition CbmKF.h:28
STL namespace.