#include <SimpleVertexTree.h>
Public Member Functions | |
void | fill (const TransientVertex &recv, const TrackingVertex *simv=0, reco::RecoToSimCollection *recSimColl=0, const float &time=0.) |
void | fill () |
void | fill (const TransientVertex &recv, const TrackingVertex *simv=0, const float &time=0.) |
void | fill (const TransientVertex &recv, const float &time=0.) |
void | fill (const TrackingVertex *simv) |
SimpleVertexTree (const char *fitterName="VertexFitter", TrackAssociatorByChi2 *associator=0) | |
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 | |
TrackAssociatorByChi2 * | associatorForParamAtPca |
float | chiProb |
float | chiTot |
int | maxTrack |
float | ndf |
int | numberOfVertices |
TString * | parameterNames [5] |
float | recErr [3] |
float | recPos [3] |
VertexFitterResult * | result |
float | simPos [3] |
TString | theFitterName |
bool | trackTest |
TTree * | vertexTree |
Definition at line 41 of file SimpleVertexTree.h.
SimpleVertexTree::SimpleVertexTree | ( | const char * | fitterName = "VertexFitter" , |
TrackAssociatorByChi2 * | associator = 0 |
||
) |
The constructor
fitterName | The name of the TTree, and of the associated histograms. |
Definition at line 10 of file SimpleVertexTree.cc.
References VertexFitterResult::chi2Information(), defineTrackBranch(), 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.
: theFitterName(filterName), associatorForParamAtPca(associator) { vertexTree = new TTree(filterName, "Vertex fit results"); // trackTest // = SimpleConfigurable<bool> (false, "SimpleVertexTree:trackTest").value(); // if (trackTest) { // maxTrack = SimpleConfigurable<int> (100, "SimpleVertexTree:maximumTracksToStore").value(); // } else maxTrack = 0; result = new VertexFitterResult(maxTrack, associator); vertexTree->Branch("vertex",(void *)result->vertexPresent(),"vertex/I"); vertexTree->Branch("simPos",(void *)result->simVertexPos(),"X/F:Y/F:Z/F"); vertexTree->Branch("recPos",(void *)result->recVertexPos(),"X/F:Y/F:Z/F"); vertexTree->Branch("recErr",(void *)result->recVertexErr(),"X/F:Y/F:Z/F"); vertexTree->Branch("nbrTrk",(void *)result->trackInformation(),"Sim/I:Rec/I:Shared/I"); vertexTree->Branch("chiTot",(void *)result->chi2Information (),"chiTot/F"); vertexTree->Branch("ndf",(void *)(result->chi2Information ()+1),"ndf/F"); vertexTree->Branch("chiProb",(void *)(result->chi2Information ()+2),"chiProb/F"); vertexTree->Branch("time",(void *)result->time(),"time/F"); parameterNames[0] = new TString("ptinv"); parameterNames[1] = new TString("theta"); parameterNames[2] = new TString("phi"); parameterNames[3] = new TString("timp"); parameterNames[4] = new TString("limp"); vertexTree->Branch("simTracks",(void *)result->numberSimTracks(),"simTracks/I"); vertexTree->Branch("simTrack_recIndex",(void *)result->simTrack_recIndex(),"simTrack_recIndex[simTracks]/I"); defineTrackBranch("sim", "Par", &VertexFitterResult::simParameters, "simTracks"); vertexTree->Branch("recTracks",(void *)result->numberRecTracks(),"recTracks/I"); vertexTree->Branch("recTrack_simIndex",(void *)result->recTrack_simIndex(),"recTrack_simIndex[recTracks]/I"); vertexTree->Branch("recTrack_weight",(void *)result->recTrackWeight(),"recTrack_weight[recTracks]/F"); defineTrackBranch("rec", "Par", &VertexFitterResult::recParameters, "recTracks"); defineTrackBranch("ref", "Par", &VertexFitterResult::refParameters, "recTracks"); defineTrackBranch("rec", "Err", &VertexFitterResult::recErrors, "recTracks"); defineTrackBranch("ref", "Err", &VertexFitterResult::refErrors, "recTracks"); numberOfVertices = 0; }
SimpleVertexTree::~SimpleVertexTree | ( | ) | [virtual] |
Definition at line 70 of file SimpleVertexTree.cc.
References chiProb, gather_cfg::cout, numberOfVertices, theFitterName, and vertexTree.
{ std::cout << std::endl<< "End of SimpleVertexTree for "<< theFitterName << std::endl; std::cout << std::endl<< "Number of vertices fit: "<< numberOfVertices<<std::endl; // // save current root directory // TDirectory* rootDir = gDirectory; // // close files // vertexTree->GetDirectory()->cd(); vertexTree->Write(); if (numberOfVertices>0) { TH1F *resX = new TH1F(theFitterName + "_ResX","Residual x coordinate: "+theFitterName, 100, -0.03, 0.03); TH1F *resY = new TH1F(theFitterName + "_ResY","Residual y coordinate: "+theFitterName, 100, -0.03, 0.03); TH1F *resZ = new TH1F(theFitterName + "_ResZ","Residual z coordinate: "+theFitterName, 100, -0.03, 0.03); TH1F *pullX = new TH1F(theFitterName + "_PullX","Pull x coordinate: "+theFitterName, 100, -10., 10.); TH1F *pullY = new TH1F(theFitterName + "_PullY","Pull y coordinate: "+theFitterName, 100, -10., 10.); TH1F *pullZ = new TH1F(theFitterName + "_PullZ","Pull z coordinate: "+theFitterName, 100, -10., 10.); TH1F *chiNorm = new TH1F(theFitterName + "_ChiNorm","Normalized chi-square: " +theFitterName, 100, 0., 10.); TH1F *chiProb = new TH1F(theFitterName + "_ChiProb","Chi-square probability: "+theFitterName, 100, 0., 1.); vertexTree->Project(theFitterName + "_ResX", "(simPos.X-recPos.X)"); vertexTree->Project(theFitterName + "_ResY", "(simPos.Y-recPos.Y)"); vertexTree->Project(theFitterName + "_ResZ", "(simPos.Z-recPos.Z)"); vertexTree->Project(theFitterName + "_PullX", "(simPos.X-recPos.X)/recErr.X"); vertexTree->Project(theFitterName + "_PullY", "(simPos.Y-recPos.Y)/recErr.Y"); vertexTree->Project(theFitterName + "_PullZ", "(simPos.Z-recPos.Z)/recErr.Z"); vertexTree->Project(theFitterName + "_ChiNorm", "chiTot/ndf"); vertexTree->Project(theFitterName + "_ChiProb", "chiProb"); std::cout << "Mean of Residual distribution X: "<< resX->GetMean()<<std::endl; std::cout << "Mean of Residual distribution Y: "<< resY->GetMean()<<std::endl; std::cout << "Mean of Residual distribution Z: "<< resZ->GetMean()<<std::endl; std::cout << "RMS of Residual distribution X: "<< resX->GetRMS()<<std::endl; std::cout << "RMS of Residual distribution Y: "<< resY->GetRMS()<<std::endl; std::cout << "RMS of Residual distribution Z: "<< resZ->GetRMS()<<std::endl; std::cout << "Mean of Pull distribution X: "<< pullX->GetMean()<<std::endl; std::cout << "Mean of Pull distribution Y: "<< pullY->GetMean()<<std::endl; std::cout << "Mean of Pull distribution Z: "<< pullZ->GetMean()<<std::endl; std::cout << "RMS of Pull distribution X: "<< pullX->GetRMS()<<std::endl; std::cout << "RMS of Pull distribution Y: "<< pullY->GetRMS()<<std::endl; std::cout << "RMS of Pull distribution Z: "<< pullZ->GetRMS()<<std::endl; std::cout << "Average chi-square probability: "<< chiProb->GetMean()<<std::endl; std::cout << "Average normalized chi-square : "<< chiNorm->GetMean()<<std::endl; resX->Write(); resY->Write(); resZ->Write(); pullX->Write(); pullY->Write(); pullZ->Write(); chiNorm->Write(); chiProb->Write(); } delete vertexTree; // // restore directory // rootDir->cd(); std::cout << std::endl; }
void SimpleVertexTree::defineTrackBranch | ( | const TString & | prefix, |
const TString & | type, | ||
const float *(VertexFitterResult::*)(const int) const | pfunc, | ||
const TString & | index | ||
) | [private] |
Definition at line 57 of file SimpleVertexTree.cc.
References i, parameterNames, result, and vertexTree.
Referenced by SimpleVertexTree().
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 140 of file SimpleVertexTree.cc.
References VertexFitterResult::fill(), fill(), and result.
void SimpleVertexTree::fill | ( | void | ) |
To be used if one wants to record "Failed Fits", e.g. to synchronise two Trees
Definition at line 167 of file SimpleVertexTree.cc.
References numberOfVertices, VertexFitterResult::reset(), result, and vertexTree.
Referenced by fill().
{ ++numberOfVertices; TDirectory* rootDir = gDirectory; // // fill entry // static int nFill(0); vertexTree->GetDirectory()->cd(); vertexTree->Fill(); if ( (++nFill)%1000==0 ) vertexTree->AutoSave(); // // restore directory // rootDir->cd(); result->reset(); }
void SimpleVertexTree::fill | ( | const TrackingVertex * | simv | ) |
Entry for a TkSimVertex, without RecVertex.
Definition at line 146 of file SimpleVertexTree.cc.
References VertexFitterResult::fill(), fill(), and result.
{ result->fill(TransientVertex(), simv, 0); fill(); }
void SimpleVertexTree::fill | ( | const TransientVertex & | recv, |
const TrackingVertex * | simv = 0 , |
||
reco::RecoToSimCollection * | recSimColl = 0 , |
||
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 133 of file SimpleVertexTree.cc.
References VertexFitterResult::fill(), fill(), and result.
Referenced by KVFTest::analyze().
void SimpleVertexTree::fill | ( | const TransientVertex & | recv, |
const TrackingVertex * | simv = 0 , |
||
const float & | time = 0. |
||
) |
Definition at line 108 of file SimpleVertexTree.h.
float SimpleVertexTree::chiProb [private] |
Definition at line 99 of file SimpleVertexTree.h.
Referenced by ~SimpleVertexTree().
float SimpleVertexTree::chiTot [private] |
Definition at line 99 of file SimpleVertexTree.h.
int SimpleVertexTree::maxTrack [private] |
Definition at line 106 of file SimpleVertexTree.h.
Referenced by SimpleVertexTree().
float SimpleVertexTree::ndf [private] |
Definition at line 99 of file SimpleVertexTree.h.
int SimpleVertexTree::numberOfVertices [private] |
Definition at line 100 of file SimpleVertexTree.h.
Referenced by fill(), SimpleVertexTree(), and ~SimpleVertexTree().
TString* SimpleVertexTree::parameterNames[5] [private] |
Definition at line 107 of file SimpleVertexTree.h.
Referenced by defineTrackBranch(), and SimpleVertexTree().
float SimpleVertexTree::recErr[3] [private] |
Definition at line 98 of file SimpleVertexTree.h.
float SimpleVertexTree::recPos[3] [private] |
Definition at line 97 of file SimpleVertexTree.h.
VertexFitterResult* SimpleVertexTree::result [private] |
Definition at line 102 of file SimpleVertexTree.h.
Referenced by defineTrackBranch(), fill(), and SimpleVertexTree().
float SimpleVertexTree::simPos[3] [private] |
Definition at line 96 of file SimpleVertexTree.h.
TString SimpleVertexTree::theFitterName [private] |
Definition at line 103 of file SimpleVertexTree.h.
Referenced by ~SimpleVertexTree().
bool SimpleVertexTree::trackTest [private] |
Definition at line 105 of file SimpleVertexTree.h.
TTree* SimpleVertexTree::vertexTree [private] |
Definition at line 101 of file SimpleVertexTree.h.
Referenced by defineTrackBranch(), fill(), SimpleVertexTree(), and ~SimpleVertexTree().