PandaRoot
RhoCandidate.h
Go to the documentation of this file.
1 //****************************************************************************
2 //* This file is part of PandaRoot. *
3 //* *
4 //* PandaRoot is distributed under the terms of the *
5 //* GNU General Public License (GPL) version 3, *
6 //* copied verbatim in the file "LICENSE". *
7 //* *
8 //* Copyright (C) 2006 - 2024 FAIR GmbH and copyright holders of PandaRoot *
9 //* The copyright holders are listed in the file "COPYRIGHTHOLDERS". *
10 //* The authors are listed in the file "AUTHORS". *
11 //****************************************************************************
12 
13 #ifndef RHOCANDIDATE_H
14 #define RHOCANDIDATE_H
15 // //
17 // RhoCandidate //
18 // //
19 // Contains (polymorphic) contents for RhoCandidate objects //
20 // Candidate "Tracks" or "Particles" for analysis use //
21 // //
22 // Author List: //
23 // Sascha Berger, RUB, Feb.99 //
24 // Marcel Kunze, RUB, Aug.99 //
25 // Copyright (C) 1999-2001, Ruhr-University Bochum. //
26 // Ralf Kliemt, HIM/GSI Feb.2013 (Cleanup & Restructuring) //
27 // //
29 
30 #include <iostream>
31 #include <vector>
32 
33 #include "TNtuple.h"
34 #include "TVector3.h"
35 #include "RhoVector3Err.h"
36 #include "TLorentzVector.h"
37 #include "TObject.h"
38 #include "RhoCandList.h"
39 #include "PndPidCandidate.h"
40 #include "RhoLorentzVectorErr.h"
41 
42 // class VAbsVertex;
43 // class VAbsTruth;
44 class TParticlePDG;
45 class RhoError;
47 
48 #define MATRIXSIZE 28
49 #define MAXDAUG 10
50 
51 // ---------------------
52 // -- Class Interface --
53 // ---------------------
54 class RhoCandidate : public TObject {
55 
56  protected:
57  //--------------------
58  // Instance Members --
59  //--------------------
60 
61  // Fast mode : Do not use error matrix
62  Bool_t fFastMode;
63 
64  // Lock : if true, modifications are a fatal error
65  Bool_t fLocked;
66 
67  // The mother
69  // Counted reference to the vertex
71 
72  // Identity
73  const TParticlePDG *fPdtEntry;
74  int fPdgCode;
75 
76  // is this a resonance ?
77  Bool_t fIsAResonance;
78 
79  // Monte-Carlo truth
80  // VAbsTruth* fTruth; //! Pointer to MCTruth info
81 
82  // Interface to objects storable in micro database
83  PndPidCandidate *fMicroCand; // !Pointer to micro data
84 
85  UInt_t fTrackNumber;
86  UInt_t fUid;
87 
88  // Daughters
89  // std::vector<RhoCandidate*> fDaughters; //! List of Daughters
91  Int_t fNDaug;
92  // Constraints
93  // TConstraint* fConstraints[5]; //! Array of constraints
94  Short_t fNCons;
95 
96  UInt_t fMarker[4];
97 
98  // added by K Goetzen
99  double fPidLH[30];
100 
101  private:
102  Double_t fChi2;
103 
104  //[ralfk:may2013] changed mc truth access to direct pointers
105  // int fMcIdx;
106  RhoCandidate *fMcTruth;
107 
108  RhoCandidate *fFit;
109  // the params
110  Char_t fCharge; // The electrical charge
111  Float_t fXposition; // The origin in x
112  Float_t fYposition; // The origin in y
113  Float_t fZposition; // The origin in z
114  Double_t fXmomentum; // The momentum in x
115  Double_t fYmomentum; // The momentum in y
116  Double_t fZmomentum; // The momentum in z
117  Double_t fEnergy; // The total energy
118  Float_t fErrP7[MATRIXSIZE]; // The symmetric 7*7 error matrix
119 
120  public:
121  //--------------------
122  // Public interface --
123  //--------------------
124 
125  //
126  // Constructors
127  //
128 
130  RhoCandidate();
131 
139  RhoCandidate(const TLorentzVector &v, Double_t charge = 0, RhoVector3Err *vp = nullptr);
140 
147  RhoCandidate(const TVector3 &v, const TParticlePDG *pdt, RhoVector3Err *vp = nullptr);
148 
150  RhoCandidate(const RhoCandidate &);
151  // RhoCandidate ( const RhoCandidate* );
152 
153  // Special constructor from MicroCandidate
154  RhoCandidate(PndPidCandidate &a, Int_t n);
155  RhoCandidate(PndPidCandidate &a, Int_t n, RhoVector3Err &vp, Bool_t fast = kFALSE);
156 
157  // RhoCandidate ( TLorentzVector p4,
158  // RhoError& p4Err,
159  // RhoCandListIterator& iterDau,
160  // RhoVector3Err& theVertex,
161  // const TParticlePDG* hypo = 0 );
162 
163  //
164  // Destructor
165  //
166  virtual ~RhoCandidate();
167 
168  //
169  // Accessors to contained information
170  //
171  Double_t Mass() const;
172  Double_t GetMass() const { return Mass(); }
173 
174  const TParticlePDG *PdtEntry() const;
175  int PdgCode() { return fPdgCode; }
176 
177  //
178  // By convention, the 4-momentum is given :
179  // - at the point of closest approach from the origin for non-composites
180  // - at the decay vertex for composite candidates
181  //
182 
183  // The accessors and modifiers have been noved from TFitParams
184  // accessors
185  Double_t GetCharge() const { return fCharge; }
186  Double_t Charge() const { return fCharge; }
187  TVector3 GetPosition() const { return TVector3(fXposition, fYposition, fZposition); }
188  TVector3 Pos() const { return GetPosition(); } // Position where the fourmomentum is defined
189  TVector3 GetDecayPos() const { return TVector3(fDecayVtx); } // position of decay (set by vertexing, if applicable)
190  TVector3 GetMomentum() const { return TVector3(fXmomentum, fYmomentum, fZmomentum); }
191  TVector3 GetVect() const { return TVector3(fXmomentum, fYmomentum, fZmomentum); }
192  Double_t Px() const { return fXmomentum; }
193  Double_t Py() const { return fYmomentum; }
194  Double_t Pz() const { return fZmomentum; }
195  RhoVector3Err PosWCov() const;
196  RhoError PosCov() const;
197  TLorentzVector P4() const { return TLorentzVector(fXmomentum, fYmomentum, fZmomentum, fEnergy); }
198  RhoError P4Cov() const;
199  RhoError P4Err() const { return P4Cov(); }
200  RhoLorentzVectorErr P4WErr() const;
201  TVector3 P3() const { return GetMomentum(); }
202  RhoError P3Cov() const;
203  RhoVector3Err P3WErr() const;
204  Double_t E() const { return fEnergy; }
205  Double_t Energy() const { return fEnergy; }
206  Double_t GetEnergy() const { return fEnergy; }
207  Float_t *GetErrP7() { return fErrP7; }
208  Double_t EVar() const;
209  Double_t M() const;
210  Double_t P() const;
211  TMatrixD Cov7() const;
212  TMatrixD XPCov() const;
213  Double_t Pt() const { return TMath::Sqrt(fXmomentum * fXmomentum + fYmomentum * fYmomentum); }
214  TMatrixD GetDecayPosCov() const { return TMatrixD(fDecayVtx.CovMatrix()); }
215 
216  // Modifiers
217  //
218  // Set type of candidate from a particle data table entry.
219  // The setType function sets the type and therefore the mass
220  // for any RhoCandidate which has no daughter.
221  // Otherwise, i.e. for composite RhoCandidates,
222  // the setType function only sets the type, _the mass is not set_.
223  void SetType(const TParticlePDG *pdt);
224  void SetType(const char *name);
225  void SetType(int pdgcode);
226 
227  // modifiers
228 
229  void SetCharge(Double_t charge) { fCharge = (Char_t)charge; }
230  void SetMass(Double_t mass);
231  void SetEnergy(Double_t newE); // Set the energy (at constant mass, by default)
232  void SetE(Double_t energy) { fEnergy = energy; }
233  void SetMassAndEnergy(Double_t mass, Double_t energy);
234  void SetPosition(const TVector3 &pos);
235  void SetPos(const TVector3 &pos) { SetPosition(pos); }
236  void SetMomentum(Double_t newP); // Set the momentum (at constant mass, by default)
237  void SetP3(const TVector3 &p3);
238  void SetVect(const TVector3 &p3) { SetP3(p3); }
239  void SetP4(Double_t mass, const TVector3 &p3);
240  void SetP4(const TLorentzVector &p4);
241  void SetP7(const TVector3 &pos, const TLorentzVector &p4);
242  void SetCovPos(const TMatrixD &covPos);
243  void SetCovP4(const TMatrixD &covP4);
244  void SetCov7(const TMatrixD &cov7);
245  void SetCov7(const TMatrixD &covPos, const TMatrixD &covP4);
246  void SetCov7(const TMatrixD &covPos, const TMatrixD &covP4, const TMatrixD &covPosP4);
247  void Set(const TVector3 &pos, const TLorentzVector &p4, const TMatrixD &cov7);
248  void Set(Double_t mass, const RhoVector3Err &posErr, const RhoVector3Err &p3Err, const TMatrixD &xpErr);
249  void SetErr(Float_t *err)
250  {
251  if (err != 0)
252  for (int i = 0; i < MATRIXSIZE; i++) {
253  fErrP7[i] = err[i];
254  }
255  }
256 
257  void SetRecoCandidate(PndPidCandidate &micro) { fMicroCand = &micro; }
258  // Allow the candidate to fly or not to fly
259  // This overrides the default which is to consider the
260  // candidate a resonance (a non-flying state) if it
261  // has a proper decay length c*tau of less than a nanometer
262  // according to the pdt table.
263  // This function MUST be called before any call to setVertex,
264  // or inconsistencies might occur. A inconsistent call to
265  // one of these functions when the decay vertex is already set
266  // will abort.
267  void SetFly();
268  void SetNoFly();
269 
270  void Boost(Double_t bx, Double_t by, Double_t bz);
271  void Boost(const TVector3 &p) { Boost(p.X(), p.Y(), p.Z()); }
272 
273  //
274  // Genealogy
275  //
276 
277  // access to the mother
278  const RhoCandidate *TheMother() const { return fTheMother; }
280  // access to daughters
281  Int_t NDaughters() const;
282  // RhoCandListIterator DaughterIterator() const;
283  // const RhoCandidate* Daughter ( Int_t n ) const;
284  RhoCandidate *Daughter(Int_t n);
285  void RemoveAssociations();
286 
287  Bool_t IsComposite() const; // true when nDaugthers > 0
288  Bool_t IsAResonance() const; // true when c*tau is very small ( < 1nm )
289 
290  //
291  // Access to vertex information
292  //
293 
295 
296  Double_t Chi2() const { return fChi2; }
297  void SetChi2(Double_t chi2) { fChi2 = chi2; }
298  void SetFit(RhoCandidate *b) { fFit = b; }
299  RhoCandidate *GetFit() const { return fFit; }
300  // recursive function that invalidates the vertex fit.
301  // it sets all the vertices to the UnFitted status.
302  // warning : the function may trig lots of copies and clones.
303  void InvalidateFit(); // defunct a.t.m.
304 
305  //
306  // Check for overlap
307  //
308 
309  // the function isCloneOf returns true if the two
310  // Candidates have shared in the past a common Base.
311  // For instance, an identified muon candidate is the
312  // clone of the original Candidate it originates from
313  // (note that "clone" must be taken in an enlarged sense
314  // here. In that example for instance, the two clones have
315  // different masses).
316  Bool_t IsCloneOf(const RhoCandidate &, Bool_t checkType = kFALSE) const;
317 
318  // two Candidates are equal if they share the same Base
319  Bool_t operator==(const RhoCandidate *) const;
320 
321  // and different if they don not
322  Bool_t operator!=(RhoCandidate *) const;
323 
324  // this function returns the pointer of the first clone
325  // of a given Candidate found in the decay tree of the
326  // present Candidate, if present, else the null pointer.
327  const RhoCandidate *CloneInTree(const RhoCandidate &) const;
328 
329  //
330  // Locks
331  //
332  // set the flag that prevents a cand from being changed
333 
334  void Lock() { fLocked = kTRUE; }
335  void UnLock() { fLocked = kFALSE; }
336  Bool_t IsLocked() { return fLocked; }
337  //
338  // Constraints
339  //
340 
341  // modifiers
342  // TConstraint& AddConstraint( TConstraint& );
343  // TConstraint& AddConstraint( TConstraint::Type );
344  // void RemoveConstraint( const TConstraint& );
345  // void RemoveConstraint( TConstraint::Type );
346 
347  // access
348  // Int_t NConstraints() const;
349  // Int_t NumberOfConstraints() const { return NConstraints(); }
350  // const TConstraint* Constraint( Int_t i ) ;
351  // const TConstraint* Constraint( TConstraint::Type ) const;
352 
353  //
354  // Origin Point
355  // - The origin point is the position of the production vertex
356  // when present, else the "origin" (0.,0.,0.)
357  // - The origin point is not necessarily the point
358  // where the 4-momentum P4() is defined. For composite
359  // candidates for instance, P4() is given at the decay point,
360  // not the production point.
361  // - To get the momentum and covariance at the origin point :
362  //
363  // TLorentzVector p4origin = cand->p4( cand->origin() );
364  //
365 
366  TVector3 Origin() const;
367 
368  //
369  // Operations
370  //
371 
373 
374  //
375  // Prints
376  //
377  void PrintOn(std::ostream &o = std::cout) const;
378 
379  void SetFast(Bool_t yesno) { fFastMode = yesno; }
380  Bool_t IsFast() const { return fFastMode; }
381 
383 
385  //************** added Combine for more candidates K.Goetzen, 05/2008
393  RhoCandidate *
396  RhoCandidate *c9, RhoCandidate *c10);
397 
398  // two Candidates overlap if they are identical
399  // (same pointers), equal (same Base),
400  // clones (same Uid), representing a same
401  // reconstructed object, or having daughters
402  // that are overlapping
403 
404  Bool_t Overlaps(const RhoCandidate *c) const
405  {
406  return ((fMarker[0] & c->fMarker[0]) != 0 || (fMarker[1] & c->fMarker[1]) != 0 || (fMarker[2] & c->fMarker[2]) != 0 || (fMarker[3] & c->fMarker[3]) != 0);
407  }
408 
409  Bool_t Equals(const RhoCandidate *c) const
410  {
411  return ((fMarker[0] == c->fMarker[0]) && (fMarker[1] == c->fMarker[1]) && (fMarker[2] == c->fMarker[2]) && (fMarker[3] == c->fMarker[3]));
412  }
413 
414  UInt_t GetMarker(UInt_t m = 0) const
415  {
416  if (m < 4) {
417  return fMarker[m];
418  } else {
419  return 0;
420  }
421  }
422 
423  void SetMarker(UInt_t l, UInt_t m);
424  void SetMarker(UInt_t n);
425 
426  Int_t GetTrackNumber() const { return fTrackNumber; }
427  void SetTrackNumber(Int_t trnum = -1) { fTrackNumber = trnum; };
428  Int_t Uid() const { return fUid; }
429  void SetUid(UInt_t uid = 0);
430 
431  // void SetTrajectory ( const TLorentzVector& p4, const RhoError& p4Err,
432  // Int_t charge,const TParticlePDG* hypo,
433  // RhoVector3Err dVtx );
434 
435  void SetPidInfo(double *pidinfo = nullptr);
436  void SetPidInfo(int hypo, double value);
437  double GetPidInfo(int hypo);
438  const double *GetPidInfo() const;
439 
440  //[ralfk:may2013] changed mc truth access to direct pointers
441  // void SetMcIdx ( int idx ) {fMcIdx=idx;}
442  // int GetMcIdx() {return fMcIdx;}
443  void SetMcTruth(RhoCandidate *mct) { fMcTruth = mct; }
444  RhoCandidate *GetMcTruth() const { return fMcTruth; }
445 
446  Bool_t IsLocal() const { return kTRUE; }
447 
448  // Set the decay vertex - operators can do that
449  void SetDecayVtx(RhoVector3Err theVtx);
450 
451  // Sets the mother link and adds a daughter link in the mother
452  void SetMotherLink(RhoCandidate *m, bool verbose = true);
453 
454  // Drop the mother link
455  void DropMotherLink();
456 
457  // Add a daughter link and set the daughters mother link
458  // void AddDaughterLink ( const RhoCandidate* );
459 
460  // Add a daughter link without touching the daughters
461  void AddDaughterLinkSimple(const RhoCandidate *, bool verbose = true);
462 
463  // Remove a daughter
465 
466  Double_t Correlation(Int_t x1, Int_t x2, const TMatrixD &m, const TMatrixD &cov) const;
467 
468  ClassDef(RhoCandidate, 3) // Candidate base class
469 
470  // friend class TBooster;
471  // friend class TOperatorBase;
472 };
473 
474 // standalone print
475 std::ostream &operator<<(std::ostream &o, const RhoCandidate &);
476 
477 #endif
UInt_t GetMarker(UInt_t m=0) const
Definition: RhoCandidate.h:414
int fPdgCode
Pointer to particle database.
Definition: RhoCandidate.h:74
void AddDaughterLinkSimple(const RhoCandidate *, bool verbose=true)
UInt_t fTrackNumber
Definition: RhoCandidate.h:85
void SetP7(const TVector3 &pos, const TLorentzVector &p4)
TLorentzVector P4() const
Definition: RhoCandidate.h:197
Int_t Uid() const
Definition: RhoCandidate.h:428
Double_t P() const
Double_t Correlation(Int_t x1, Int_t x2, const TMatrixD &m, const TMatrixD &cov) const
void SetPos(const TVector3 &pos)
Definition: RhoCandidate.h:235
__m128 m
Definition: P4_F32vec4.h:38
Bool_t IsAResonance() const
RhoVector3Err fDecayVtx
Do not stream.
Definition: RhoCandidate.h:70
void SetMotherLink(RhoCandidate *m, bool verbose=true)
Double_t GetEnergy() const
Definition: RhoCandidate.h:206
void SetUid(UInt_t uid=0)
Bool_t operator==(const RhoCandidate *) const
void SetMassAndEnergy(Double_t mass, Double_t energy)
RhoCandidate * GetMcTruth() const
Definition: RhoCandidate.h:444
Int_t GetTrackNumber() const
Definition: RhoCandidate.h:426
static T Sqrt(const T &x)
Definition: PndCAMath.h:57
TVector3 Origin() const
Double_t GetCharge() const
Definition: RhoCandidate.h:185
Bool_t fLocked
Do not stream.
Definition: RhoCandidate.h:65
RhoCandidate & operator=(const RhoCandidate &)
void SetVect(const TVector3 &p3)
Definition: RhoCandidate.h:238
RhoVector3Err PosWCov() const
Float_t * GetErrP7()
Definition: RhoCandidate.h:207
void RemoveDaughter(RhoCandidate *)
void Boost(Double_t bx, Double_t by, Double_t bz)
TMatrixD GetDecayPosCov() const
Definition: RhoCandidate.h:214
TVector3 GetVect() const
Definition: RhoCandidate.h:191
void SetNoFly()
void DropMotherLink()
RhoCandidate * fTheMother
Do not stream.
Definition: RhoCandidate.h:68
void SetPosition(const TVector3 &pos)
__m128 v
Definition: P4_F32vec4.h:15
unsigned int i
Definition: P4_F32vec4.h:33
void SetE(Double_t energy)
Definition: RhoCandidate.h:232
Bool_t Overlaps(const RhoCandidate *c) const
Definition: RhoCandidate.h:404
const RhoError & CovMatrix() const
Definition: RhoVector3Err.h:68
void SetRecoCandidate(PndPidCandidate &micro)
Definition: RhoCandidate.h:257
Bool_t IsLocal() const
Definition: RhoCandidate.h:446
Int_t fNDaug
List of Daughters.
Definition: RhoCandidate.h:91
void SetFast(Bool_t yesno)
Definition: RhoCandidate.h:379
Double_t E() const
Definition: RhoCandidate.h:204
void SetType(const TParticlePDG *pdt)
Bool_t IsComposite() const
void SetP4(Double_t mass, const TVector3 &p3)
std::ostream & operator<<(std::ostream &o, const RhoCandidate &)
Bool_t fIsAResonance
Definition: RhoCandidate.h:77
Double_t GetMass() const
Definition: RhoCandidate.h:172
void SetMomentum(Double_t newP)
RhoCandidate * Combine(RhoCandidate *c)
#define MAXDAUG
Definition: RhoCandidate.h:49
RhoError P3Cov() const
PndPidCandidate * fMicroCand
Rsonance flag.
Definition: RhoCandidate.h:83
void SetP3(const TVector3 &p3)
void InvalidateFit()
double fPidLH[30]
Overlap.
Definition: RhoCandidate.h:99
RhoVector3Err P3WErr() const
Double_t Energy() const
Definition: RhoCandidate.h:205
RhoCandidate * Daughter(Int_t n)
Bool_t IsCloneOf(const RhoCandidate &, Bool_t checkType=kFALSE) const
Bool_t operator!=(RhoCandidate *) const
Bool_t IsLocked()
Definition: RhoCandidate.h:336
void RemoveAssociations()
RhoCandidate * fDaughters[MAXDAUG]
unique number
Definition: RhoCandidate.h:90
RhoCandidate * TheMother()
Definition: RhoCandidate.h:279
TMatrixD XPCov() const
Double_t Pt() const
Definition: RhoCandidate.h:213
void SetErr(Float_t *err)
Definition: RhoCandidate.h:249
TVector3 Pos() const
Definition: RhoCandidate.h:188
void Set(const TVector3 &pos, const TLorentzVector &p4, const TMatrixD &cov7)
void SetChi2(Double_t chi2)
Definition: RhoCandidate.h:297
const RhoCandidate * CloneInTree(const RhoCandidate &) const
void SetMass(Double_t mass)
RhoLorentzVectorErr P4WErr() const
virtual ~RhoCandidate()
RhoError P4Cov() const
Double_t M() const
void Boost(const TVector3 &p)
Definition: RhoCandidate.h:271
Double_t EVar() const
Double_t Chi2() const
Definition: RhoCandidate.h:296
void SetFit(RhoCandidate *b)
Definition: RhoCandidate.h:298
void SetPidInfo(double *pidinfo=nullptr)
const TParticlePDG * PdtEntry() const
void SetMcTruth(RhoCandidate *mct)
Definition: RhoCandidate.h:443
void SetCovP4(const TMatrixD &covP4)
Double_t Mass() const
void SetTrackNumber(Int_t trnum=-1)
Definition: RhoCandidate.h:427
#define MATRIXSIZE
Definition: RhoCandidate.h:48
void SetEnergy(Double_t newE)
Double_t Pz() const
Definition: RhoCandidate.h:194
RhoError PosCov() const
TVector3 GetDecayPos() const
Definition: RhoCandidate.h:189
void PrintOn(std::ostream &o=std::cout) const
Int_t NDaughters() const
const TParticlePDG * fPdtEntry
Vertex.
Definition: RhoCandidate.h:73
void SetCharge(Double_t charge)
Definition: RhoCandidate.h:229
Double_t Px() const
Definition: RhoCandidate.h:192
TMatrixD Cov7() const
void SetMarker(UInt_t l, UInt_t m)
Double_t Py() const
Definition: RhoCandidate.h:193
void SetDecayVtx(RhoVector3Err theVtx)
Bool_t Equals(const RhoCandidate *c) const
Definition: RhoCandidate.h:409
Bool_t IsFast() const
Definition: RhoCandidate.h:380
Bool_t fFastMode
Definition: RhoCandidate.h:62
UInt_t fUid
Micro association.
Definition: RhoCandidate.h:86
const double * GetPidInfo() const
void SetCovPos(const TMatrixD &covPos)
RhoError P4Err() const
Definition: RhoCandidate.h:199
RhoVector3Err DecayVtx()
Definition: RhoCandidate.h:294
void SetCov7(const TMatrixD &cov7)
Double_t Charge() const
Definition: RhoCandidate.h:186
TVector3 GetPosition() const
Definition: RhoCandidate.h:187
TVector3 GetMomentum() const
Definition: RhoCandidate.h:190
const RhoCandidate * TheMother() const
Definition: RhoCandidate.h:278
TMatrixT< double > TMatrixD
Definition: PndLmdDim.h:64
UInt_t fMarker[4]
Number of constraints.
Definition: RhoCandidate.h:96
Short_t fNCons
Definition: RhoCandidate.h:94
TVector3 P3() const
Definition: RhoCandidate.h:201
PndPidCandidate * GetRecoCandidate() const
Definition: RhoCandidate.h:382
RhoCandidate * GetFit() const
Definition: RhoCandidate.h:299