Public Member Functions | Private Member Functions | Private Attributes | List of all members
garana::StructuredRecoTree Class Reference

#include <StructuredRecoTree.h>

Inheritance diagram for garana::StructuredRecoTree:
garana::RecoTree garana::TreeReader

Public Member Functions

 StructuredRecoTree (TTree *tree=0)
 
const size_t NTrack () const override
 number of tracks in this event More...
 
const size_t NVertex () const override
 number of vertices in this event More...
 
const size_t NVee () const override
 number of vees in this event More...
 
const size_t NCalCluster () const override
 number of ECal clusters in this event More...
 
const TLorentzVector * TrackVertex (const size_t &itrack) const override
 4-position of track's assumed start point More...
 
const TLorentzVector * TrackEnd (const size_t &itrack) const override
 4-position of track's assumed end point More...
 
const size_t NTrackHit (const size_t &itrack) const override
 number of reconstructed hits used in track fit More...
 
const TVector3 * TrackMomBeg (const size_t &itrack) const override
 momentum 3-vector as track's assumed start point More...
 
const TVector3 * TrackMomEnd (const size_t &itrack) const override
 momentum 3-vector as track's assumed end point More...
 
const float TrackVtxDirectionX (const size_t &itrack) const override
 x-direction cosine at track vertex More...
 
const float TrackVtxDirectionY (const size_t &itrack) const override
 y-direction cosine at track vertex More...
 
const float TrackVtxDirectionZ (const size_t &itrack) const override
 z-direction cosine at track vertex More...
 
const float TrackEndDirectionX (const size_t &itrack) const override
 x-direction cosine at track end More...
 
const float TrackEndDirectionY (const size_t &itrack) const override
 y-direction cosine at track end More...
 
const float TrackEndDirectionZ (const size_t &itrack) const override
 z-direction cosine at track end More...
 
const float TrackLenFwd (const size_t &itrack) const override
 track length from forward fit More...
 
const float TrackLenBkd (const size_t &itrack) const override
 track length from backward fit More...
 
const float TrackIonizFwd (const size_t &itrack) const override
 track average ionization rate from forward fit More...
 
const float TrackIonizBkd (const size_t &itrack) const override
 track average ionization rate from backward fit More...
 
const int TrackChiSqrFwd (const size_t &itrack) const override
 chi-squared of track fit in forward direction More...
 
const int TrackChiSqrBkd (const size_t &itrack) const override
 chi-squared of track fit in backward direction More...
 
const int TrackChgFwd (const size_t &itrack) const override
 charge sign of track if current hypothesis correct More...
 
const int TrackChgBkd (const size_t &itrack) const override
 charge sign of track if current hypothesis wrong More...
 
void TrackParBeg (const size_t &itrack, float pars[5]) const override
 track fit parameters at the track's assumed start More...
 
void TrackParEnd (const size_t &itrack, float pars[5]) const override
 track fit parameters at the track's assumed end More...
 
void TrackCovarBeg (const size_t &itrack, float pars[15]) const override
 track fit covariance matrix (assume symmetry) at track's assumed start More...
 
void TrackCovarEnd (const size_t &itrack, float pars[15]) const override
 track fit covariance matrix (assume symmetry) at track's assumed end More...
 
const TLorentzVector * TrackTruePosBeg (const size_t &itrack) const override
 true 4-position at track vertex [cm,ns] More...
 
const TLorentzVector * TrackTruePosEnd (const size_t &itrack) const override
 true 4-position at track end [cm,ns] More...
 
const TLorentzVector * TrackTrueMomBeg (const size_t &itrack) const override
 true 4-momentum at track vertex [GeV/c,GeV] More...
 
const TLorentzVector * TrackTrueMomEnd (const size_t &itrack) const override
 true 4-momentum at track end [GeV/c,GeV] More...
 
const float TrackTrueEnergy (const size_t &icluster) const override
 total associated true energy deposited with ith cluster More...
 
const size_t TrackNTrueTrack (const size_t &icluster) const override
 number of MCParticles associated with ith cluster More...
 
const int TrackTrkIdMaxDeposit (const size_t &icluster) const override
 trackID of the MCParticle depositing the most energy in ith cluster More...
 
const float TrackMaxDeposit (const size_t &icluster) const override
 maximum true deposited energy from a single MCParticle More...
 
const pair< int, float > * TrackTrueDeposit (const size_t &itrack, size_t &itrue) const override
 
const vector< pair< int, float > > * TrackTrueDeposits (const size_t &itrack) const override
 
const TLorentzVector * GetVertex (const size_t &ivertex) const override
 vertex 4-position for vertex with index ivertex More...
 
void VertexCovariance (const size_t &ivertex, float covar[][3]) const override
 given a vertex index, fill given position covariance matrix More...
 
const TLorentzVector * VeeVertex (const size_t &ivee) const override
 =============== Vee ======================= More...
 
void VeeCovariance (const size_t &ivee, float covar[][3]) const override
 given a vee index, fill given position covariance matrix More...
 
const vector< TLorentzVector > * VeeMomentumPerHypothesis (const size_t &ivee) const override
 4-momentum for vee if it matches {Kshort, Lambda1, Lambda2} More...
 
const float VeeChiSquared (const size_t &ivee) const override
 
const CaloClusterGetCalCluster (const size_t &icluster) const
 ================ ECal cluster ====================== More...
 
const TLorentzVector * CalClustPosition (const size_t &icluster) const override
 
const float CalClustEnergy (const size_t &icluster) const override
 reconstructed ECal cluster energy More...
 
const float CalClustEnergyError (const size_t &icluster) const override
 
const float CalClustTrueEnergy (const size_t &icluster) const override
 total associated true energy deposited with ith cluster More...
 
const size_t CalClustNTrueTrack (const size_t &icluster) const override
 number of MCParticles associated with ith cluster More...
 
const int CalClustTrkIdMaxDeposit (const size_t &icluster) const override
 trackID of the MCParticle depositing the most energy in ith cluster More...
 
const float CalClustMaxDeposit (const size_t &icluster) const override
 maximum true deposited energy from a single MCParticle More...
 
const std::pair< int, float > * CalClustTrueDeposit (const size_t &icluster, const size_t &itrack) const override
 
const float CalClustTimeDifference (const size_t &icluster) const override
 
const float * CalClustShape (const size_t &icluster) const override
 
const float CalClustTheta (const size_t &icluster) const override
 
const float CalClustPhi (const size_t &icluster) const override
 
const vector< TVector3 > * CalClustEigenVecs (const size_t &icluster) const override
 
- Public Member Functions inherited from garana::RecoTree
virtual ~RecoTree ()
 
const float TrackMaxDepositFrac (const size_t &itrack) const
 largest fraction of total energy contributed by single MCParticle More...
 
const float CalClustMaxDepositFrac (const size_t &icluster) const
 largest fraction of total energy contributed by single MCParticle More...
 
void GetTrackG4PIndices (const size_t &itrk, vector< UInt_t > &ig4ps) const
 given a track index, fill a given vector with matched G4 particle indices More...
 
void GetVertexTrackIndices (const size_t &ivtx, vector< UInt_t > &itracks) const
 given a vertex index, fill a given vector with matched track indices More...
 
void GetVeeTrackIndices (const size_t &ivee, vector< UInt_t > &itracks) const
 given a vee index, fill a given vector with matched G4 particle indices More...
 
void GetCalClusterTrackIndices (const size_t &iclust, vector< UInt_t > &itracks) const
 given a calocluster index, fill a given vector with matched track indices More...
 
void GetCalClusterG4Indices (const size_t &iclust, vector< UInt_t > &ig4ps) const
 given a calocluster index, fill a given vector with matched G4Particle indices More...
 
- Public Member Functions inherited from garana::TreeReader
virtual ~TreeReader ()
 
void SetupRead (TTree *tree)
 
TTree * GetInputTree ()
 
size_t NEntries () const
 
virtual void GetEntry (const UInt_t &ientry)
 
Int_t Event () const
 
const TObjArray * GetBranchList () const
 
void Fill ()
 
void Write ()
 
void CheckOpt (char opt)
 
bool BlockWrite () const
 

Private Member Functions

bool SetBranchAddresses () override
 

Private Attributes

vector< Track > * fTracks = nullptr
 
vector< Vee > * fVees = nullptr
 
vector< Vertex > * fVertices = nullptr
 
vector< CaloCluster > * fCalClusters = nullptr
 
vector< vector< Int_t > > * fVertTrackEnds = nullptr
 
vector< vector< Int_t > > * fVeeTrackEnds = nullptr
 
TBranch * b_Tracks = nullptr
 
TBranch * b_Vees = nullptr
 
TBranch * b_Vertices = nullptr
 
TBranch * b_CalClusters = nullptr
 
TBranch * b_VertTrackEnds = nullptr
 
TBranch * b_VeeTrackEnds = nullptr
 

Additional Inherited Members

- Protected Attributes inherited from garana::RecoTree
vector< vector< UInt_t > > * fTrackG4PIndices = nullptr
 
vector< vector< UInt_t > > * fVertTrackIndices = nullptr
 
vector< vector< UInt_t > > * fVeeTrackIndices = nullptr
 
vector< vector< UInt_t > > * fCalClusterTrackIndices = nullptr
 
vector< vector< UInt_t > > * fCalClusterG4Indices = nullptr
 
TBranch * b_TrackG4PIndices = nullptr
 
TBranch * b_VertTrackIndices = nullptr
 
TBranch * b_VeeTrackIndices = nullptr
 
TBranch * b_CalClusterTrackIndices = nullptr
 
TBranch * b_CalClusterG4Indices = nullptr
 
- Protected Attributes inherited from garana::TreeReader
const std::string treename
 
char fOpt = 'r'
 
TTree * fTreeIn = nullptr
 pointer to the analyzed TTree or TChain More...
 
UInt_t fCurrentEntry = UINT_MAX
 
TBranch * b_Event = nullptr
 
Int_t fEvent = -1
 event number for tree entry More...
 

Detailed Description

Definition at line 18 of file StructuredRecoTree.h.

Constructor & Destructor Documentation

StructuredRecoTree::StructuredRecoTree ( TTree *  tree = 0)

Definition at line 13 of file StructuredRecoTree.cxx.

13  {
14  SetupRead(tree); //initialize tree pointer in TreeReader instance and set branch addresses
15 }
void SetupRead(TTree *tree)
Definition: TreeReader.cxx:6

Member Function Documentation

const vector< TVector3 > * StructuredRecoTree::CalClustEigenVecs ( const size_t &  icluster) const
overridevirtual

Implements garana::RecoTree.

Definition at line 373 of file StructuredRecoTree.cxx.

373  {
374  return fCalClusters->at(icluster).EigenVecs();
375 }
vector< CaloCluster > * fCalClusters
const float StructuredRecoTree::CalClustEnergy ( const size_t &  icluster) const
overridevirtual

reconstructed ECal cluster energy

Implements garana::RecoTree.

Definition at line 331 of file StructuredRecoTree.cxx.

331  {
332  return GetCalCluster(icluster)->Energy();
333 }
float const & Energy() const
Definition: CaloCluster.cxx:58
const CaloCluster * GetCalCluster(const size_t &icluster) const
================ ECal cluster ======================
const float StructuredRecoTree::CalClustEnergyError ( const size_t &  icluster) const
overridevirtual

Implements garana::RecoTree.

Definition at line 335 of file StructuredRecoTree.cxx.

335  {
336  return GetCalCluster(icluster)->EnergyError();
337 }
float const & EnergyError() const
Definition: CaloCluster.cxx:62
const CaloCluster * GetCalCluster(const size_t &icluster) const
================ ECal cluster ======================
const float StructuredRecoTree::CalClustMaxDeposit ( const size_t &  icluster) const
overridevirtual

maximum true deposited energy from a single MCParticle

Implements garana::RecoTree.

Definition at line 348 of file StructuredRecoTree.cxx.

348  {
349  if(GetCalCluster(icluster))
350  return GetCalCluster(icluster)->MaxDeposit();
351  else
352  return 0.;
353 }
const CaloCluster * GetCalCluster(const size_t &icluster) const
================ ECal cluster ======================
float const & MaxDeposit() const
const size_t StructuredRecoTree::CalClustNTrueTrack ( const size_t &  icluster) const
overridevirtual

number of MCParticles associated with ith cluster

Implements garana::RecoTree.

Definition at line 342 of file StructuredRecoTree.cxx.

342  {
343  return GetCalCluster(icluster)->NIdes();
344 }
const CaloCluster * GetCalCluster(const size_t &icluster) const
================ ECal cluster ======================
size_t NIdes() const
Definition: CaloCluster.cxx:86
const float StructuredRecoTree::CalClustPhi ( const size_t &  icluster) const
overridevirtual

Implements garana::RecoTree.

Definition at line 369 of file StructuredRecoTree.cxx.

369  {
370  return fCalClusters->at(icluster).Phi();
371 }
vector< CaloCluster > * fCalClusters
const TLorentzVector * StructuredRecoTree::CalClustPosition ( const size_t &  icluster) const
overridevirtual

Implements garana::RecoTree.

Definition at line 327 of file StructuredRecoTree.cxx.

327  {
328  return GetCalCluster(icluster)->Position();
329 }
const CaloCluster * GetCalCluster(const size_t &icluster) const
================ ECal cluster ======================
const TLorentzVector * Position() const
Definition: CaloCluster.cxx:50
const float * StructuredRecoTree::CalClustShape ( const size_t &  icluster) const
overridevirtual

Implements garana::RecoTree.

Definition at line 361 of file StructuredRecoTree.cxx.

361  {
362  return fCalClusters->at(icluster).Shape();
363 }
vector< CaloCluster > * fCalClusters
const float StructuredRecoTree::CalClustTheta ( const size_t &  icluster) const
overridevirtual

Implements garana::RecoTree.

Definition at line 365 of file StructuredRecoTree.cxx.

365  {
366  return fCalClusters->at(icluster).Theta();
367 }
vector< CaloCluster > * fCalClusters
const float StructuredRecoTree::CalClustTimeDifference ( const size_t &  icluster) const
overridevirtual

Implements garana::RecoTree.

Definition at line 357 of file StructuredRecoTree.cxx.

357  {
358  return fCalClusters->at(icluster).TimeDifference();
359 }
vector< CaloCluster > * fCalClusters
const int StructuredRecoTree::CalClustTrkIdMaxDeposit ( const size_t &  icluster) const
overridevirtual

trackID of the MCParticle depositing the most energy in ith cluster

Implements garana::RecoTree.

Definition at line 345 of file StructuredRecoTree.cxx.

345  {
346  return GetCalCluster(icluster)->TrackIdMaxDep();
347 }
const CaloCluster * GetCalCluster(const size_t &icluster) const
================ ECal cluster ======================
int const & TrackIdMaxDep() const
const std::pair< int, float > * StructuredRecoTree::CalClustTrueDeposit ( const size_t &  icluster,
const size_t &  itrack 
) const
overridevirtual

Implements garana::RecoTree.

Definition at line 354 of file StructuredRecoTree.cxx.

354  {
355  return GetCalCluster(icluster)->GetTrackIdEdep(itrack);
356 }
const CaloCluster * GetCalCluster(const size_t &icluster) const
================ ECal cluster ======================
const std::pair< int, float > * GetTrackIdEdep(const size_t &iide) const
Definition: CaloCluster.cxx:90
const float StructuredRecoTree::CalClustTrueEnergy ( const size_t &  icluster) const
overridevirtual

total associated true energy deposited with ith cluster

Implements garana::RecoTree.

Definition at line 339 of file StructuredRecoTree.cxx.

339  {
340  return GetCalCluster(icluster)->TotalTrueEnergy();
341 }
const CaloCluster * GetCalCluster(const size_t &icluster) const
================ ECal cluster ======================
float TotalTrueEnergy() const
Definition: CaloCluster.cxx:94
const CaloCluster * StructuredRecoTree::GetCalCluster ( const size_t &  icluster) const

================ ECal cluster ======================

Definition at line 321 of file StructuredRecoTree.cxx.

321  {
322  if(fCalClusters->size()==0)
323  return nullptr;
324 
325  return &(fCalClusters->at(icluster));
326 }
vector< CaloCluster > * fCalClusters
const TLorentzVector * StructuredRecoTree::GetVertex ( const size_t &  ivertex) const
overridevirtual

vertex 4-position for vertex with index ivertex

Implements garana::RecoTree.

Definition at line 294 of file StructuredRecoTree.cxx.

294  {
295  return fVertices->at(ivertex).GetVertex();
296 }
const size_t StructuredRecoTree::NCalCluster ( ) const
overridevirtual

number of ECal clusters in this event

Implements garana::RecoTree.

Definition at line 69 of file StructuredRecoTree.cxx.

69  {
70  return fCalClusters->size();
71 }
vector< CaloCluster > * fCalClusters
const size_t StructuredRecoTree::NTrack ( ) const
overridevirtual

number of tracks in this event

Implements garana::RecoTree.

Definition at line 60 of file StructuredRecoTree.cxx.

60  {
61  return fTracks->size();
62 }
const size_t StructuredRecoTree::NTrackHit ( const size_t &  itrack) const
overridevirtual

number of reconstructed hits used in track fit

Implements garana::RecoTree.

Definition at line 81 of file StructuredRecoTree.cxx.

81  {
82  return fTracks->at(itrack).fNHits;
83 }
const size_t StructuredRecoTree::NVee ( ) const
overridevirtual

number of vees in this event

Implements garana::RecoTree.

Definition at line 66 of file StructuredRecoTree.cxx.

66  {
67  return fVees->size();
68 }
const size_t StructuredRecoTree::NVertex ( ) const
overridevirtual

number of vertices in this event

Implements garana::RecoTree.

Definition at line 63 of file StructuredRecoTree.cxx.

63  {
64  return fVertices->size();
65 }
bool StructuredRecoTree::SetBranchAddresses ( )
overrideprivatevirtual

Implements garana::TreeReader.

Definition at line 17 of file StructuredRecoTree.cxx.

17  {
18 
19  cout << "setting StructuredRecoTree branch addresses" << endl;
20 
21  fTreeIn->SetBranchAddress("Event", &fEvent, &b_Event );
22  cout << "set Event" << endl;
23  fTreeIn->SetBranchAddress("Tracks", &fTracks , &b_Tracks );
24  cout << "set Track" << endl;
25  fTreeIn->SetBranchAddress("Vees", &fVees , &b_Vees );
26  cout << "set Vees" << endl;
27  fTreeIn->SetBranchAddress("Vertices", &fVertices , &b_Vertices );
28  cout << "set Vertices" << endl;
29  fTreeIn->SetBranchAddress("CalClusters", &fCalClusters , &b_CalClusters );
30  cout << "set CalClusters" << endl;
31  fTreeIn->SetBranchAddress("TrackG4Indices", &fTrackG4PIndices , &b_TrackG4PIndices );
32  cout << "set TrackG4Indices" << endl;
33  fTreeIn->SetBranchAddress("VertTrackIndices", &fVertTrackIndices, &b_VertTrackIndices );
34  fTreeIn->SetBranchAddress("VertTrackEnds", &fVertTrackEnds , &b_VertTrackEnds );
35  fTreeIn->SetBranchAddress("VeeTrackIndices", &fVeeTrackIndices , &b_VeeTrackIndices );
36  fTreeIn->SetBranchAddress("VeeTrackEnds", &fVeeTrackEnds , &b_VeeTrackEnds );
37  fTreeIn->SetBranchAddress("CalTrackIndices", &fCalClusterTrackIndices , &b_CalClusterTrackIndices );
38  //fTreeIn->SetBranchAddress("CalTrackEnds", &fCalTrackEnds , &b_CalTrackEnds );
39  fTreeIn->SetBranchAddress("CalG4Indices", &fCalClusterG4Indices , &b_CalClusterG4Indices );
40  /*if(fGeo->HasMuonDetector()){
41  fTreeIn->SetBranchAddress("MuIDClusters", "vector<gar::rec::Cluster>", &fMuClusters);
42  }*/
43 
44  //TODO: add Assns for reco r&d products
45  /*if(fAnaMode == "reco"){
46  fRecoTree->Branch("TPCHits", "vector<gar::rec::Hit>", &fTPCHits);
47  fRecoTree->Branch("TPCClusters", "vector<gar::rec::TPCCluster>", &fTPCClusters);
48  fRecoTree->Branch("CalHits", "vector<gar::rec::CaloHit>", &fCalHits);
49  if(fGeo->HasMuonDetector()){
50  fRecoTree->Branch("MuIDHits", "vector<gar::rec::CaloHit>", &fMuHits);
51  }
52  }*/
53 
54  cout << "done." << endl;
55 
56  return true;
57 }// end const'or
vector< vector< UInt_t > > * fTrackG4PIndices
Definition: RecoTree.h:108
vector< vector< Int_t > > * fVertTrackEnds
vector< vector< UInt_t > > * fCalClusterTrackIndices
Definition: RecoTree.h:111
vector< vector< UInt_t > > * fVertTrackIndices
Definition: RecoTree.h:109
TTree * fTreeIn
pointer to the analyzed TTree or TChain
Definition: TreeReader.h:51
vector< CaloCluster > * fCalClusters
TBranch * b_VeeTrackIndices
Definition: RecoTree.h:116
TBranch * b_VertTrackIndices
Definition: RecoTree.h:115
TBranch * b_TrackG4PIndices
Definition: RecoTree.h:114
vector< vector< Int_t > > * fVeeTrackEnds
vector< vector< UInt_t > > * fCalClusterG4Indices
Definition: RecoTree.h:112
TBranch * b_Event
Definition: TreeReader.h:54
TBranch * b_CalClusterTrackIndices
Definition: RecoTree.h:117
vector< vector< UInt_t > > * fVeeTrackIndices
Definition: RecoTree.h:110
Int_t fEvent
event number for tree entry
Definition: TreeReader.h:55
QTextStream & endl(QTextStream &s)
TBranch * b_CalClusterG4Indices
Definition: RecoTree.h:118
const int StructuredRecoTree::TrackChgBkd ( const size_t &  itrack) const
overridevirtual

charge sign of track if current hypothesis wrong

Implements garana::RecoTree.

Definition at line 148 of file StructuredRecoTree.cxx.

148  {
149  return fTracks->at(itrack).fChgBac;
150 }
const int StructuredRecoTree::TrackChgFwd ( const size_t &  itrack) const
overridevirtual

charge sign of track if current hypothesis correct

Implements garana::RecoTree.

Definition at line 144 of file StructuredRecoTree.cxx.

144  {
145  return fTracks->at(itrack).fChgFwd;
146 }
const int StructuredRecoTree::TrackChiSqrBkd ( const size_t &  itrack) const
overridevirtual

chi-squared of track fit in backward direction

Implements garana::RecoTree.

Definition at line 140 of file StructuredRecoTree.cxx.

140  {
141  return fTracks->at(itrack).fChiBac;
142 }
const int StructuredRecoTree::TrackChiSqrFwd ( const size_t &  itrack) const
overridevirtual

chi-squared of track fit in forward direction

Implements garana::RecoTree.

Definition at line 136 of file StructuredRecoTree.cxx.

136  {
137  return fTracks->at(itrack).fChiFwd;
138 }
void StructuredRecoTree::TrackCovarBeg ( const size_t &  itrack,
float  pars[15] 
) const
overridevirtual

track fit covariance matrix (assume symmetry) at track's assumed start

Implements garana::RecoTree.

Definition at line 162 of file StructuredRecoTree.cxx.

162  {
163  for(size_t i=0; i<15; i++)
164  covar[i] = fTracks->at(itrack).fCovMatBeg[i];
165 }
void StructuredRecoTree::TrackCovarEnd ( const size_t &  itrack,
float  pars[15] 
) const
overridevirtual

track fit covariance matrix (assume symmetry) at track's assumed end

Implements garana::RecoTree.

Definition at line 167 of file StructuredRecoTree.cxx.

167  {
168  for(size_t i=0; i<15; i++)
169  covar[i] = fTracks->at(itrack).fCovMatEnd[i];
170 }
const TLorentzVector * StructuredRecoTree::TrackEnd ( const size_t &  itrack) const
overridevirtual

4-position of track's assumed end point

Implements garana::RecoTree.

Definition at line 78 of file StructuredRecoTree.cxx.

78  {
79  return &(fTracks->at(itrack).fEnd);
80 }
const float StructuredRecoTree::TrackEndDirectionX ( const size_t &  itrack) const
overridevirtual

x-direction cosine at track end

Implements garana::RecoTree.

Definition at line 108 of file StructuredRecoTree.cxx.

108  {
109  return fTracks->at(itrack).fEndDir.X();
110 }
const float StructuredRecoTree::TrackEndDirectionY ( const size_t &  itrack) const
overridevirtual

y-direction cosine at track end

Implements garana::RecoTree.

Definition at line 112 of file StructuredRecoTree.cxx.

112  {
113  return fTracks->at(itrack).fEndDir.Y();
114 }
const float StructuredRecoTree::TrackEndDirectionZ ( const size_t &  itrack) const
overridevirtual

z-direction cosine at track end

Implements garana::RecoTree.

Definition at line 116 of file StructuredRecoTree.cxx.

116  {
117  return fTracks->at(itrack).fEndDir.Z();
118 }
const float StructuredRecoTree::TrackIonizBkd ( const size_t &  itrack) const
overridevirtual

track average ionization rate from backward fit

Implements garana::RecoTree.

Definition at line 132 of file StructuredRecoTree.cxx.

132  {
133  return fTracks->at(itrack).fLenBac;
134 }
const float StructuredRecoTree::TrackIonizFwd ( const size_t &  itrack) const
overridevirtual

track average ionization rate from forward fit

Implements garana::RecoTree.

Definition at line 128 of file StructuredRecoTree.cxx.

128  {
129  return fTracks->at(itrack).fIonFwd;
130 }
const float StructuredRecoTree::TrackLenBkd ( const size_t &  itrack) const
overridevirtual

track length from backward fit

Implements garana::RecoTree.

Definition at line 124 of file StructuredRecoTree.cxx.

124  {
125  return fTracks->at(itrack).fLenBac;
126 }
const float StructuredRecoTree::TrackLenFwd ( const size_t &  itrack) const
overridevirtual

track length from forward fit

Implements garana::RecoTree.

Definition at line 120 of file StructuredRecoTree.cxx.

120  {
121  return fTracks->at(itrack).fLenFwd;
122 }
const float StructuredRecoTree::TrackMaxDeposit ( const size_t &  itrack) const
overridevirtual

maximum true deposited energy from a single MCParticle

Implements garana::RecoTree.

Definition at line 268 of file StructuredRecoTree.cxx.

268  {
269  //TODO check if we actually want total EDep rather than just from the leading contributor
270  auto it = std::max_element(fTracks->at(itrack).fTrueEnergy.begin(),fTracks->at(itrack).fTrueEnergy.end(),
271  [](const std::pair<int,float>& lhs,const std::pair<int,float>& rhs) -> bool {return lhs.second < rhs.second; }
272  );
273  if(it==fTracks->at(itrack).fTrueEnergy.end()) {
274  cout << "max element not found for itrack = " << itrack << "! (this is just a warning; is it all noise hits?)" << endl;
275  //cout << "fTracks->at(itrack).fTrueEnergy.size() = " << fTracks->at(itrack).fTrueEnergy.size() << endl;
276  return 0.;
277  }
278  return (*it).second;
279  //return fTracks->at(itrack).fTrueEnergy.at(TrackTrkIdMaxDeposit(itrack)).second;
280 }
QTextStream & endl(QTextStream &s)
const TVector3 * StructuredRecoTree::TrackMomBeg ( const size_t &  itrack) const
overridevirtual

momentum 3-vector as track's assumed start point

Implements garana::RecoTree.

Definition at line 84 of file StructuredRecoTree.cxx.

84  {
85  TVector3* v = new TVector3(fTracks->at(itrack).fVtxDir);
86  (*v) *= fTracks->at(itrack).fMomBeg;
87  return v;
88 }
const TVector3 * StructuredRecoTree::TrackMomEnd ( const size_t &  itrack) const
overridevirtual

momentum 3-vector as track's assumed end point

Implements garana::RecoTree.

Definition at line 90 of file StructuredRecoTree.cxx.

90  {
91  TVector3* v = new TVector3(fTracks->at(itrack).fEndDir);
92  (*v) *= fTracks->at(itrack).fMomEnd;
93  return v;
94 }
const size_t StructuredRecoTree::TrackNTrueTrack ( const size_t &  itrack) const
overridevirtual

number of MCParticles associated with ith cluster

Implements garana::RecoTree.

Definition at line 239 of file StructuredRecoTree.cxx.

239  {
240  return fTracks->at(itrack).fTrueEnergy.size();
241 }
void StructuredRecoTree::TrackParBeg ( const size_t &  itrack,
float  pars[5] 
) const
overridevirtual

track fit parameters at the track's assumed start

Implements garana::RecoTree.

Definition at line 152 of file StructuredRecoTree.cxx.

152  {
153  for(size_t i=0; i<5; i++)
154  pars[i] = fTracks->at(itrack).fTrackParBeg[i];
155 }
void StructuredRecoTree::TrackParEnd ( const size_t &  itrack,
float  pars[5] 
) const
overridevirtual

track fit parameters at the track's assumed end

Implements garana::RecoTree.

Definition at line 157 of file StructuredRecoTree.cxx.

157  {
158  for(size_t i=0; i<5; i++)
159  pars[i] = fTracks->at(itrack).fTrackParEnd[i];
160 }
const int StructuredRecoTree::TrackTrkIdMaxDeposit ( const size_t &  itrack) const
overridevirtual

trackID of the MCParticle depositing the most energy in ith cluster

Implements garana::RecoTree.

Definition at line 243 of file StructuredRecoTree.cxx.

243  {
244  /*int imax=-1;
245  float emax = -1.;
246  for(size_t i=0; i<fTracks->at(itrack).fTrueEnergy.size(); i++) {
247 
248  if(fTracks->at(itrack).fTrueEnergy.at(i).second > emax){
249  emax = fTracks->at(itrack).fTrueEnergy.at(i).second;
250  imax = i;
251  }
252  }
253 
254  return fTracks->at(itrack).fTrueEnergy.at(imax).first;*/
255 
256  auto it = std::max_element(fTracks->at(itrack).fTrueEnergy.begin(),fTracks->at(itrack).fTrueEnergy.end(),
257  [](const std::pair<int,float>& lhs,const std::pair<int,float>& rhs) -> bool {return lhs.second < rhs.second; }
258  );
259  if(it==fTracks->at(itrack).fTrueEnergy.end()) {
260  cout << "max element not found for itrack = " << itrack << "! (this is just a warning; is it all noise hits?)" << endl;
261  //cout << "fTracks->at(itrack).fTrueEnergy.size() = " << fTracks->at(itrack).fTrueEnergy.size() << endl;
262  return -INT_MAX;
263  }
264  return (*it).first;
265 
266 }
QTextStream & endl(QTextStream &s)
const pair< int, float > * StructuredRecoTree::TrackTrueDeposit ( const size_t &  itrack,
size_t &  itrue 
) const
overridevirtual

Implements garana::RecoTree.

Definition at line 282 of file StructuredRecoTree.cxx.

282  {
283 
284  return &(fTracks->at(itrack).fTrueEnergy.at(itrue));
285 
286 }
const vector< pair< int, float > > * StructuredRecoTree::TrackTrueDeposits ( const size_t &  itrack) const
overridevirtual

Implements garana::RecoTree.

Definition at line 288 of file StructuredRecoTree.cxx.

288  {
289 
290  return &(fTracks->at(itrack).fTrueEnergy);
291 }
const float StructuredRecoTree::TrackTrueEnergy ( const size_t &  itrack) const
overridevirtual

total associated true energy deposited with ith cluster

Implements garana::RecoTree.

Definition at line 230 of file StructuredRecoTree.cxx.

230  {
231  float etrue = 0.;
232  for(auto const& e : fTracks->at(itrack).fTrueEnergy)
233  etrue += e.second;
234 
235  //return FLT_MAX-itrack*0;
236  return etrue;
237 }
const double e
const TLorentzVector * StructuredRecoTree::TrackTrueMomBeg ( const size_t &  itrack) const
overridevirtual

true 4-momentum at track vertex [GeV/c,GeV]

Implements garana::RecoTree.

Definition at line 202 of file StructuredRecoTree.cxx.

202  {
203  TLorentzVector outvec(0,0,0,0);
204  float ptot = 0.;
205  for(size_t imctrk=0; imctrk<fTracks->at(itrack).fTrueMomVtx.size(); imctrk++) {
206  TLorentzVector tmpvec = fTracks->at(itrack).fTrueMomVtx[imctrk].second;
207  tmpvec *= fTracks->at(itrack).fTrueMomVtx[imctrk].second.P();
208  outvec += tmpvec;
209  ptot += fTracks->at(itrack).fTrueMomVtx[imctrk].second.P();
210  }
211 
212  outvec *= 1.0/ptot;
213  return new TLorentzVector(outvec);
214 }
const TLorentzVector * StructuredRecoTree::TrackTrueMomEnd ( const size_t &  itrack) const
overridevirtual

true 4-momentum at track end [GeV/c,GeV]

Implements garana::RecoTree.

Definition at line 216 of file StructuredRecoTree.cxx.

216  {
217  TLorentzVector outvec(0,0,0,0);
218  float ptot = 0.;
219  for(size_t imctrk=0; imctrk<fTracks->at(itrack).fTrueMomEnd.size(); imctrk++) {
220  TLorentzVector tmpvec = fTracks->at(itrack).fTrueMomEnd[imctrk].second;
221  tmpvec *= fTracks->at(itrack).fTrueMomEnd[imctrk].second.P();
222  outvec += tmpvec;
223  ptot += fTracks->at(itrack).fTrueMomEnd[imctrk].second.P();
224  }
225 
226  outvec *= 1.0/ptot;
227  return new TLorentzVector(outvec);
228 }
const TLorentzVector * StructuredRecoTree::TrackTruePosBeg ( const size_t &  itrack) const
overridevirtual

true 4-position at track vertex [cm,ns]

Implements garana::RecoTree.

Definition at line 172 of file StructuredRecoTree.cxx.

172  {
173  TLorentzVector outvec(0,0,0,0);
174  float ptot = 0.;
175  for(size_t imctrk=0; imctrk<fTracks->at(itrack).fTruePosVtx.size(); imctrk++ ) {
176  TLorentzVector tmpvec = fTracks->at(itrack).fTruePosVtx[imctrk].second;
177  tmpvec *= fTracks->at(itrack).fTrueMomVtx[imctrk].second.P();
178  outvec += tmpvec;
179  ptot += fTracks->at(itrack).fTrueMomVtx[imctrk].second.P();
180  }
181 
182  outvec *= 1.0/ptot;
183  return new TLorentzVector(outvec);
184 }
const TLorentzVector * StructuredRecoTree::TrackTruePosEnd ( const size_t &  itrack) const
overridevirtual

true 4-position at track end [cm,ns]

Implements garana::RecoTree.

Definition at line 186 of file StructuredRecoTree.cxx.

186  {
187  TLorentzVector outvec(0,0,0,0);
188  float ptot = 0.;
189  std::cout << "loop over -- " << fTracks->at(itrack).fTruePosEnd.size() << " -- true positions" << std::endl;
190  for(size_t imctrk=0; imctrk<fTracks->at(itrack).fTruePosEnd.size(); imctrk++) {
191  TLorentzVector tmpvec = fTracks->at(itrack).fTruePosEnd[imctrk].second;
192  std::cout << " submomentum: " << fTracks->at(itrack).fTrueMomEnd[imctrk].second.P() << std::endl;
193  tmpvec *= fTracks->at(itrack).fTrueMomEnd[imctrk].second.P();
194  outvec += tmpvec;
195  ptot += fTracks->at(itrack).fTrueMomEnd[imctrk].second.P();
196  }
197  std::cout << "total true momentum: " << std::to_string(ptot) << std::endl;
198  outvec *= 1.0/ptot;
199  return new TLorentzVector(outvec);
200 }
std::string to_string(ModuleType const mt)
Definition: ModuleType.h:34
QTextStream & endl(QTextStream &s)
const TLorentzVector * StructuredRecoTree::TrackVertex ( const size_t &  itrack) const
overridevirtual

4-position of track's assumed start point

Implements garana::RecoTree.

Definition at line 74 of file StructuredRecoTree.cxx.

74  {
75  return &(fTracks->at(itrack).fVtx);
76 }
const float StructuredRecoTree::TrackVtxDirectionX ( const size_t &  itrack) const
overridevirtual

x-direction cosine at track vertex

Implements garana::RecoTree.

Definition at line 96 of file StructuredRecoTree.cxx.

96  {
97  return fTracks->at(itrack).fVtxDir.X();
98 }
const float StructuredRecoTree::TrackVtxDirectionY ( const size_t &  itrack) const
overridevirtual

y-direction cosine at track vertex

Implements garana::RecoTree.

Definition at line 100 of file StructuredRecoTree.cxx.

100  {
101  return fTracks->at(itrack).fVtxDir.Y();
102 }
const float StructuredRecoTree::TrackVtxDirectionZ ( const size_t &  itrack) const
overridevirtual

z-direction cosine at track vertex

Implements garana::RecoTree.

Definition at line 104 of file StructuredRecoTree.cxx.

104  {
105  return fTracks->at(itrack).fVtxDir.Z();
106 }
const float StructuredRecoTree::VeeChiSquared ( const size_t &  ivee) const
overridevirtual

Implements garana::RecoTree.

Definition at line 316 of file StructuredRecoTree.cxx.

316  {
317  return fVees->at(ivee).GetChiSqr();
318 }
void StructuredRecoTree::VeeCovariance ( const size_t &  ivee,
float  covar[][3] 
) const
overridevirtual

given a vee index, fill given position covariance matrix

Implements garana::RecoTree.

Definition at line 308 of file StructuredRecoTree.cxx.

308  {
309  return fVees->at(ivee).GetCovar(covar);
310 }
const vector< TLorentzVector > * StructuredRecoTree::VeeMomentumPerHypothesis ( const size_t &  ivee) const
overridevirtual

4-momentum for vee if it matches {Kshort, Lambda1, Lambda2}

Implements garana::RecoTree.

Definition at line 312 of file StructuredRecoTree.cxx.

312  {
313  return fVees->at(ivee).GetMomentaPerHypothesis();
314 }
const TLorentzVector * StructuredRecoTree::VeeVertex ( const size_t &  ivee) const
overridevirtual

=============== Vee =======================

Implements garana::RecoTree.

Definition at line 304 of file StructuredRecoTree.cxx.

304  {
305  return fVees->at(ivee).GetVertex();
306 }
void StructuredRecoTree::VertexCovariance ( const size_t &  ivertex,
float  covar[][3] 
) const
overridevirtual

given a vertex index, fill given position covariance matrix

Implements garana::RecoTree.

Definition at line 298 of file StructuredRecoTree.cxx.

298  {
299  fVertices->at(ivertex).GetCovar(covar);
300 }

Member Data Documentation

TBranch* garana::StructuredRecoTree::b_CalClusters = nullptr
private

Definition at line 109 of file StructuredRecoTree.h.

TBranch* garana::StructuredRecoTree::b_Tracks = nullptr
private

Definition at line 106 of file StructuredRecoTree.h.

TBranch* garana::StructuredRecoTree::b_Vees = nullptr
private

Definition at line 107 of file StructuredRecoTree.h.

TBranch* garana::StructuredRecoTree::b_VeeTrackEnds = nullptr
private

Definition at line 112 of file StructuredRecoTree.h.

TBranch* garana::StructuredRecoTree::b_Vertices = nullptr
private

Definition at line 108 of file StructuredRecoTree.h.

TBranch* garana::StructuredRecoTree::b_VertTrackEnds = nullptr
private

Definition at line 111 of file StructuredRecoTree.h.

vector<CaloCluster>* garana::StructuredRecoTree::fCalClusters = nullptr
private

Definition at line 101 of file StructuredRecoTree.h.

vector<Track>* garana::StructuredRecoTree::fTracks = nullptr
private

Definition at line 98 of file StructuredRecoTree.h.

vector<Vee>* garana::StructuredRecoTree::fVees = nullptr
private

Definition at line 99 of file StructuredRecoTree.h.

vector<vector<Int_t> >* garana::StructuredRecoTree::fVeeTrackEnds = nullptr
private

Definition at line 104 of file StructuredRecoTree.h.

vector<Vertex>* garana::StructuredRecoTree::fVertices = nullptr
private

Definition at line 100 of file StructuredRecoTree.h.

vector<vector<Int_t> >* garana::StructuredRecoTree::fVertTrackEnds = nullptr
private

Definition at line 103 of file StructuredRecoTree.h.


The documentation for this class was generated from the following files: