CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
SimpleVertexTree Class Reference

#include <SimpleVertexTree.h>

Public Member Functions

void fill (const TransientVertex &recv, const TrackingVertex *simv=nullptr, reco::RecoToSimCollection *recSimColl=nullptr, const float &time=0.)
 
void fill (const TransientVertex &recv, const TrackingVertex *simv=nullptr, const float &time=0.)
 
void fill (const TransientVertex &recv, const float &time=0.)
 
void fill (const TrackingVertex *simv)
 
void fill ()
 
 SimpleVertexTree (const char *fitterName="VertexFitter", const MagneticField *magField=nullptr)
 
virtual ~SimpleVertexTree ()
 

Private Member Functions

void defineTrackBranch (const TString &prefix, const TString &type, const float *(VertexFitterResult::*pfunc)(const int) const, const TString &index)
 

Private Attributes

float chiProb
 
float chiTot
 
int maxTrack
 
float ndf
 
int numberOfVertices
 
TString * parameterNames [5]
 
float recErr [3]
 
float recPos [3]
 
VertexFitterResultresult
 
float simPos [3]
 
TString theFitterName
 
bool trackTest
 
TTree * vertexTree
 

Detailed Description

Definition at line 41 of file SimpleVertexTree.h.

Constructor & Destructor Documentation

◆ SimpleVertexTree()

SimpleVertexTree::SimpleVertexTree ( const char *  fitterName = "VertexFitter",
const MagneticField magField = nullptr 
)

The constructor

Parameters
fitterNameThe name of the TTree, and of the associated histograms.

Definition at line 10 of file SimpleVertexTree.cc.

References VertexFitterResult::chi2Information(), defineTrackBranch(), pdwgDoubleElectron_cfi::filterName, maxTrack, numberOfVertices, VertexFitterResult::numberRecTracks(), VertexFitterResult::numberSimTracks(), parameterNames, VertexFitterResult::recErrors(), VertexFitterResult::recParameters(), VertexFitterResult::recTrack_simIndex(), VertexFitterResult::recTrackWeight(), VertexFitterResult::recVertexErr(), VertexFitterResult::recVertexPos(), VertexFitterResult::refErrors(), VertexFitterResult::refParameters(), result, VertexFitterResult::simParameters(), VertexFitterResult::simTrack_recIndex(), VertexFitterResult::simVertexPos(), VertexFitterResult::time(), VertexFitterResult::trackInformation(), VertexFitterResult::vertexPresent(), and vertexTree.

11  vertexTree = new TTree(filterName, "Vertex fit results");
12  // trackTest
13  // = SimpleConfigurable<bool> (false, "SimpleVertexTree:trackTest").value();
14  // if (trackTest) {
15  // maxTrack = SimpleConfigurable<int> (100, "SimpleVertexTree:maximumTracksToStore").value();
16  // } else
17  maxTrack = 0;
18  result = new VertexFitterResult(maxTrack, magField);
19 
20  vertexTree->Branch("vertex", (void *)result->vertexPresent(), "vertex/I");
21  vertexTree->Branch("simPos", (void *)result->simVertexPos(), "X/F:Y/F:Z/F");
22  vertexTree->Branch("recPos", (void *)result->recVertexPos(), "X/F:Y/F:Z/F");
23  vertexTree->Branch("recErr", (void *)result->recVertexErr(), "X/F:Y/F:Z/F");
24  vertexTree->Branch("nbrTrk", (void *)result->trackInformation(), "Sim/I:Rec/I:Shared/I");
25  vertexTree->Branch("chiTot", (void *)result->chi2Information(), "chiTot/F");
26  vertexTree->Branch("ndf", (void *)(result->chi2Information() + 1), "ndf/F");
27  vertexTree->Branch("chiProb", (void *)(result->chi2Information() + 2), "chiProb/F");
28  vertexTree->Branch("time", (void *)result->time(), "time/F");
29 
30  parameterNames[0] = new TString("ptinv");
31  parameterNames[1] = new TString("theta");
32  parameterNames[2] = new TString("phi");
33  parameterNames[3] = new TString("timp");
34  parameterNames[4] = new TString("limp");
35 
36  vertexTree->Branch("simTracks", (void *)result->numberSimTracks(), "simTracks/I");
37  vertexTree->Branch("simTrack_recIndex", (void *)result->simTrack_recIndex(), "simTrack_recIndex[simTracks]/I");
38  defineTrackBranch("sim", "Par", &VertexFitterResult::simParameters, "simTracks");
39 
40  vertexTree->Branch("recTracks", (void *)result->numberRecTracks(), "recTracks/I");
41  vertexTree->Branch("recTrack_simIndex", (void *)result->recTrack_simIndex(), "recTrack_simIndex[recTracks]/I");
42  vertexTree->Branch("recTrack_weight", (void *)result->recTrackWeight(), "recTrack_weight[recTracks]/F");
43  defineTrackBranch("rec", "Par", &VertexFitterResult::recParameters, "recTracks");
44  defineTrackBranch("ref", "Par", &VertexFitterResult::refParameters, "recTracks");
45  defineTrackBranch("rec", "Err", &VertexFitterResult::recErrors, "recTracks");
46  defineTrackBranch("ref", "Err", &VertexFitterResult::refErrors, "recTracks");
47 
48  numberOfVertices = 0;
49 }
void defineTrackBranch(const TString &prefix, const TString &type, const float *(VertexFitterResult::*pfunc)(const int) const, const TString &index)
const int * simTrack_recIndex()
const int * recTrack_simIndex()
const float * simVertexPos() const
const float * refErrors(const int i) const
const float * recParameters(const int i) const
const float * recTrackWeight()
const float * recErrors(const int i) const
const int * numberSimTracks()
const float * chi2Information() const
const float * simParameters(const int i) const
const float * recVertexErr() const
const int * numberRecTracks()
const float * time() const
const float * refParameters(const int i) const
TString * parameterNames[5]
const float * recVertexPos() const
const int * trackInformation() const
const int * vertexPresent() const
VertexFitterResult * result

◆ ~SimpleVertexTree()

SimpleVertexTree::~SimpleVertexTree ( )
virtual

Definition at line 63 of file SimpleVertexTree.cc.

References chiProb, gather_cfg::cout, numberOfVertices, theFitterName, and vertexTree.

63  {
64  std::cout << std::endl << "End of SimpleVertexTree for " << theFitterName << std::endl;
65  std::cout << std::endl << "Number of vertices fit: " << numberOfVertices << std::endl;
66 
67  //
68  // save current root directory
69  //
70  TDirectory *rootDir = gDirectory;
71  //
72  // close files
73  //
74  vertexTree->GetDirectory()->cd();
75  vertexTree->Write();
76  if (numberOfVertices > 0) {
77  TH1F *resX = new TH1F(theFitterName + "_ResX", "Residual x coordinate: " + theFitterName, 100, -0.03, 0.03);
78  TH1F *resY = new TH1F(theFitterName + "_ResY", "Residual y coordinate: " + theFitterName, 100, -0.03, 0.03);
79  TH1F *resZ = new TH1F(theFitterName + "_ResZ", "Residual z coordinate: " + theFitterName, 100, -0.03, 0.03);
80  TH1F *pullX = new TH1F(theFitterName + "_PullX", "Pull x coordinate: " + theFitterName, 100, -10., 10.);
81  TH1F *pullY = new TH1F(theFitterName + "_PullY", "Pull y coordinate: " + theFitterName, 100, -10., 10.);
82  TH1F *pullZ = new TH1F(theFitterName + "_PullZ", "Pull z coordinate: " + theFitterName, 100, -10., 10.);
83  TH1F *chiNorm = new TH1F(theFitterName + "_ChiNorm", "Normalized chi-square: " + theFitterName, 100, 0., 10.);
84  TH1F *chiProb = new TH1F(theFitterName + "_ChiProb", "Chi-square probability: " + theFitterName, 100, 0., 1.);
85  vertexTree->Project(theFitterName + "_ResX", "(simPos.X-recPos.X)");
86  vertexTree->Project(theFitterName + "_ResY", "(simPos.Y-recPos.Y)");
87  vertexTree->Project(theFitterName + "_ResZ", "(simPos.Z-recPos.Z)");
88  vertexTree->Project(theFitterName + "_PullX", "(simPos.X-recPos.X)/recErr.X");
89  vertexTree->Project(theFitterName + "_PullY", "(simPos.Y-recPos.Y)/recErr.Y");
90  vertexTree->Project(theFitterName + "_PullZ", "(simPos.Z-recPos.Z)/recErr.Z");
91  vertexTree->Project(theFitterName + "_ChiNorm", "chiTot/ndf");
92  vertexTree->Project(theFitterName + "_ChiProb", "chiProb");
93  std::cout << "Mean of Residual distribution X: " << resX->GetMean() << std::endl;
94  std::cout << "Mean of Residual distribution Y: " << resY->GetMean() << std::endl;
95  std::cout << "Mean of Residual distribution Z: " << resZ->GetMean() << std::endl;
96  std::cout << "RMS of Residual distribution X: " << resX->GetRMS() << std::endl;
97  std::cout << "RMS of Residual distribution Y: " << resY->GetRMS() << std::endl;
98  std::cout << "RMS of Residual distribution Z: " << resZ->GetRMS() << std::endl;
99  std::cout << "Mean of Pull distribution X: " << pullX->GetMean() << std::endl;
100  std::cout << "Mean of Pull distribution Y: " << pullY->GetMean() << std::endl;
101  std::cout << "Mean of Pull distribution Z: " << pullZ->GetMean() << std::endl;
102  std::cout << "RMS of Pull distribution X: " << pullX->GetRMS() << std::endl;
103  std::cout << "RMS of Pull distribution Y: " << pullY->GetRMS() << std::endl;
104  std::cout << "RMS of Pull distribution Z: " << pullZ->GetRMS() << std::endl;
105  std::cout << "Average chi-square probability: " << chiProb->GetMean() << std::endl;
106  std::cout << "Average normalized chi-square : " << chiNorm->GetMean() << std::endl;
107  resX->Write();
108  resY->Write();
109  resZ->Write();
110  pullX->Write();
111  pullY->Write();
112  pullZ->Write();
113  chiNorm->Write();
114  chiProb->Write();
115  }
116  delete vertexTree;
117  //
118  // restore directory
119  //
120  rootDir->cd();
121  std::cout << std::endl;
122 }

Member Function Documentation

◆ defineTrackBranch()

void SimpleVertexTree::defineTrackBranch ( const TString &  prefix,
const TString &  type,
const float *(VertexFitterResult::*)(const int) const  pfunc,
const TString &  index 
)
private

Definition at line 51 of file SimpleVertexTree.cc.

References electrons_cff::branchName, mps_fire::i, parameterNames, hcallasereventfilter2012_cfi::prefix, result, and vertexTree.

Referenced by SimpleVertexTree().

54  {
55  TString branchName, branchVariables;
56  for (int i = 0; i < 5; i++) {
57  branchName = prefix + type + '_' + *parameterNames[i];
58  branchVariables = branchName + '[' + index + "]/F";
59  vertexTree->Branch(branchName, (void *)(result->*pfunc)(i), branchVariables);
60  }
61 }
TString * parameterNames[5]
VertexFitterResult * result

◆ fill() [1/5]

void SimpleVertexTree::fill ( const TransientVertex recv,
const TrackingVertex simv = nullptr,
reco::RecoToSimCollection recSimColl = nullptr,
const float &  time = 0. 
)

Entry for a RecVertex. If the vertex was not associated to a TkSimVertex, an empty pointer can be given (would be identical to the next method). Timing information for the fit can also be provided.

Definition at line 124 of file SimpleVertexTree.cc.

References VertexFitterResult::fill(), fill(), result, and hcalRecHitTable_cff::time.

127  {
128  result->fill(recv, simv, recSimColl, time);
129  fill();
130 }
void fill(const TransientVertex &recv, const TrackingVertex *simv=nullptr, reco::RecoToSimCollection *recSimColl=nullptr, const float &time=0)
VertexFitterResult * result

◆ fill() [2/5]

void SimpleVertexTree::fill ( const TransientVertex recv,
const TrackingVertex simv = nullptr,
const float &  time = 0. 
)

◆ fill() [3/5]

void SimpleVertexTree::fill ( const TransientVertex recv,
const float &  time = 0. 
)

Entry for a RecVertex, without associated vertex. Timing information for the fit can also be provided.

Definition at line 132 of file SimpleVertexTree.cc.

References VertexFitterResult::fill(), fill(), result, and hcalRecHitTable_cff::time.

132  {
133  result->fill(recv, nullptr, nullptr, time);
134  fill();
135 }
void fill(const TransientVertex &recv, const TrackingVertex *simv=nullptr, reco::RecoToSimCollection *recSimColl=nullptr, const float &time=0)
VertexFitterResult * result

◆ fill() [4/5]

void SimpleVertexTree::fill ( const TrackingVertex simv)

Entry for a TkSimVertex, without RecVertex.

Definition at line 137 of file SimpleVertexTree.cc.

References VertexFitterResult::fill(), fill(), and result.

137  {
138  result->fill(TransientVertex(), simv, nullptr);
139  fill();
140 }
void fill(const TransientVertex &recv, const TrackingVertex *simv=nullptr, reco::RecoToSimCollection *recSimColl=nullptr, const float &time=0)
VertexFitterResult * result

◆ fill() [5/5]

void SimpleVertexTree::fill ( void  )

To be used if one wants to record "Failed Fits", e.g. to synchronise two Trees

Definition at line 157 of file SimpleVertexTree.cc.

References numberOfVertices, VertexFitterResult::reset(), result, and vertexTree.

Referenced by fill().

157  {
159  static std::atomic<int> nFill{0};
160  TDirectory *rootDir = gDirectory;
161  //
162  // fill entry
163  //
164  vertexTree->GetDirectory()->cd();
165  vertexTree->Fill();
166  if ((++nFill) % 1000 == 0)
167  vertexTree->AutoSave();
168  //
169  // restore directory
170  //
171  rootDir->cd();
172  result->reset();
173 }
VertexFitterResult * result

Member Data Documentation

◆ chiProb

float SimpleVertexTree::chiProb
private

Definition at line 97 of file SimpleVertexTree.h.

Referenced by ~SimpleVertexTree().

◆ chiTot

float SimpleVertexTree::chiTot
private

Definition at line 97 of file SimpleVertexTree.h.

◆ maxTrack

int SimpleVertexTree::maxTrack
private

Definition at line 104 of file SimpleVertexTree.h.

Referenced by SimpleVertexTree().

◆ ndf

float SimpleVertexTree::ndf
private

Definition at line 97 of file SimpleVertexTree.h.

◆ numberOfVertices

int SimpleVertexTree::numberOfVertices
private

Definition at line 98 of file SimpleVertexTree.h.

Referenced by fill(), SimpleVertexTree(), and ~SimpleVertexTree().

◆ parameterNames

TString* SimpleVertexTree::parameterNames[5]
private

Definition at line 105 of file SimpleVertexTree.h.

Referenced by defineTrackBranch(), and SimpleVertexTree().

◆ recErr

float SimpleVertexTree::recErr[3]
private

Definition at line 96 of file SimpleVertexTree.h.

◆ recPos

float SimpleVertexTree::recPos[3]
private

Definition at line 95 of file SimpleVertexTree.h.

◆ result

VertexFitterResult* SimpleVertexTree::result
private

Definition at line 100 of file SimpleVertexTree.h.

Referenced by defineTrackBranch(), fill(), and SimpleVertexTree().

◆ simPos

float SimpleVertexTree::simPos[3]
private

Definition at line 94 of file SimpleVertexTree.h.

◆ theFitterName

TString SimpleVertexTree::theFitterName
private

Definition at line 101 of file SimpleVertexTree.h.

Referenced by ~SimpleVertexTree().

◆ trackTest

bool SimpleVertexTree::trackTest
private

Definition at line 103 of file SimpleVertexTree.h.

◆ vertexTree

TTree* SimpleVertexTree::vertexTree
private

Definition at line 99 of file SimpleVertexTree.h.

Referenced by defineTrackBranch(), fill(), SimpleVertexTree(), and ~SimpleVertexTree().