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