#include <RecoVertex/KalmanVertexFit/interface/SimpleVertexTree.h>
Public Member Functions | |
void | fill () |
To be used if one wants to record "Failed Fits", e.g. | |
void | fill (const TrackingVertex *simv) |
Entry for a TkSimVertex, without RecVertex. | |
void | fill (const TransientVertex &recv, const float &time=0.) |
Entry for a RecVertex, without associated vertex. | |
void | fill (const TransientVertex &recv, const TrackingVertex *simv=0, const float &time=0.) |
void | fill (const TransientVertex &recv, const TrackingVertex *simv=0, reco::RecoToSimCollection *recSimColl=0, const float &time=0.) |
Entry for a RecVertex. | |
SimpleVertexTree (const char *fitterName="VertexFitter", TrackAssociatorByChi2 *associator=0) | |
The constructor . | |
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.
00011 : 00012 theFitterName(filterName), associatorForParamAtPca(associator) 00013 { 00014 00015 vertexTree = new TTree(filterName, "Vertex fit results"); 00016 // trackTest 00017 // = SimpleConfigurable<bool> (false, "SimpleVertexTree:trackTest").value(); 00018 // if (trackTest) { 00019 // maxTrack = SimpleConfigurable<int> (100, "SimpleVertexTree:maximumTracksToStore").value(); 00020 // } else 00021 maxTrack = 0; 00022 result = new VertexFitterResult(maxTrack, associator); 00023 00024 vertexTree->Branch("vertex",(void *)result->vertexPresent(),"vertex/I"); 00025 vertexTree->Branch("simPos",(void *)result->simVertexPos(),"X/F:Y/F:Z/F"); 00026 vertexTree->Branch("recPos",(void *)result->recVertexPos(),"X/F:Y/F:Z/F"); 00027 vertexTree->Branch("recErr",(void *)result->recVertexErr(),"X/F:Y/F:Z/F"); 00028 vertexTree->Branch("nbrTrk",(void *)result->trackInformation(),"Sim/I:Rec/I:Shared/I"); 00029 vertexTree->Branch("chiTot",(void *)result->chi2Information (),"chiTot/F"); 00030 vertexTree->Branch("ndf",(void *)(result->chi2Information ()+1),"ndf/F"); 00031 vertexTree->Branch("chiProb",(void *)(result->chi2Information ()+2),"chiProb/F"); 00032 vertexTree->Branch("time",(void *)result->time(),"time/F"); 00033 00034 00035 parameterNames[0] = new TString("ptinv"); 00036 parameterNames[1] = new TString("theta"); 00037 parameterNames[2] = new TString("phi"); 00038 parameterNames[3] = new TString("timp"); 00039 parameterNames[4] = new TString("limp"); 00040 00041 vertexTree->Branch("simTracks",(void *)result->numberSimTracks(),"simTracks/I"); 00042 vertexTree->Branch("simTrack_recIndex",(void *)result->simTrack_recIndex(),"simTrack_recIndex[simTracks]/I"); 00043 defineTrackBranch("sim", "Par", &VertexFitterResult::simParameters, "simTracks"); 00044 00045 vertexTree->Branch("recTracks",(void *)result->numberRecTracks(),"recTracks/I"); 00046 vertexTree->Branch("recTrack_simIndex",(void *)result->recTrack_simIndex(),"recTrack_simIndex[recTracks]/I"); 00047 vertexTree->Branch("recTrack_weight",(void *)result->recTrackWeight(),"recTrack_weight[recTracks]/F"); 00048 defineTrackBranch("rec", "Par", &VertexFitterResult::recParameters, "recTracks"); 00049 defineTrackBranch("ref", "Par", &VertexFitterResult::refParameters, "recTracks"); 00050 defineTrackBranch("rec", "Err", &VertexFitterResult::recErrors, "recTracks"); 00051 defineTrackBranch("ref", "Err", &VertexFitterResult::refErrors, "recTracks"); 00052 00053 numberOfVertices = 0; 00054 }
SimpleVertexTree::~SimpleVertexTree | ( | ) | [virtual] |
Definition at line 70 of file SimpleVertexTree.cc.
References chiProb, GenMuonPlsPt100GeV_cfg::cout, lat::endl(), numberOfVertices, theFitterName, and vertexTree.
00071 { 00072 cout << endl<< "End of SimpleVertexTree for "<< theFitterName << endl; 00073 cout << endl<< "Number of vertices fit: "<< numberOfVertices<<endl; 00074 00075 // 00076 // save current root directory 00077 // 00078 TDirectory* rootDir = gDirectory; 00079 // 00080 // close files 00081 // 00082 vertexTree->GetDirectory()->cd(); 00083 vertexTree->Write(); 00084 if (numberOfVertices>0) { 00085 TH1F *resX = new TH1F(theFitterName + "_ResX","Residual x coordinate: "+theFitterName, 100, -0.03, 0.03); 00086 TH1F *resY = new TH1F(theFitterName + "_ResY","Residual y coordinate: "+theFitterName, 100, -0.03, 0.03); 00087 TH1F *resZ = new TH1F(theFitterName + "_ResZ","Residual z coordinate: "+theFitterName, 100, -0.03, 0.03); 00088 TH1F *pullX = new TH1F(theFitterName + "_PullX","Pull x coordinate: "+theFitterName, 100, -10., 10.); 00089 TH1F *pullY = new TH1F(theFitterName + "_PullY","Pull y coordinate: "+theFitterName, 100, -10., 10.); 00090 TH1F *pullZ = new TH1F(theFitterName + "_PullZ","Pull z coordinate: "+theFitterName, 100, -10., 10.); 00091 TH1F *chiNorm = new TH1F(theFitterName + "_ChiNorm","Normalized chi-square: " +theFitterName, 100, 0., 10.); 00092 TH1F *chiProb = new TH1F(theFitterName + "_ChiProb","Chi-square probability: "+theFitterName, 100, 0., 1.); 00093 vertexTree->Project(theFitterName + "_ResX", "(simPos.X-recPos.X)"); 00094 vertexTree->Project(theFitterName + "_ResY", "(simPos.Y-recPos.Y)"); 00095 vertexTree->Project(theFitterName + "_ResZ", "(simPos.Z-recPos.Z)"); 00096 vertexTree->Project(theFitterName + "_PullX", "(simPos.X-recPos.X)/recErr.X"); 00097 vertexTree->Project(theFitterName + "_PullY", "(simPos.Y-recPos.Y)/recErr.Y"); 00098 vertexTree->Project(theFitterName + "_PullZ", "(simPos.Z-recPos.Z)/recErr.Z"); 00099 vertexTree->Project(theFitterName + "_ChiNorm", "chiTot/ndf"); 00100 vertexTree->Project(theFitterName + "_ChiProb", "chiProb"); 00101 cout << "Mean of Residual distribution X: "<< resX->GetMean()<<endl; 00102 cout << "Mean of Residual distribution Y: "<< resY->GetMean()<<endl; 00103 cout << "Mean of Residual distribution Z: "<< resZ->GetMean()<<endl; 00104 cout << "RMS of Residual distribution X: "<< resX->GetRMS()<<endl; 00105 cout << "RMS of Residual distribution Y: "<< resY->GetRMS()<<endl; 00106 cout << "RMS of Residual distribution Z: "<< resZ->GetRMS()<<endl; 00107 cout << "Mean of Pull distribution X: "<< pullX->GetMean()<<endl; 00108 cout << "Mean of Pull distribution Y: "<< pullY->GetMean()<<endl; 00109 cout << "Mean of Pull distribution Z: "<< pullZ->GetMean()<<endl; 00110 cout << "RMS of Pull distribution X: "<< pullX->GetRMS()<<endl; 00111 cout << "RMS of Pull distribution Y: "<< pullY->GetRMS()<<endl; 00112 cout << "RMS of Pull distribution Z: "<< pullZ->GetRMS()<<endl; 00113 cout << "Average chi-square probability: "<< chiProb->GetMean()<<endl; 00114 cout << "Average normalized chi-square : "<< chiNorm->GetMean()<<endl; 00115 resX->Write(); 00116 resY->Write(); 00117 resZ->Write(); 00118 pullX->Write(); 00119 pullY->Write(); 00120 pullZ->Write(); 00121 chiNorm->Write(); 00122 chiProb->Write(); 00123 } 00124 delete vertexTree; 00125 // 00126 // restore directory 00127 // 00128 rootDir->cd(); 00129 cout << endl; 00130 00131 }
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().
00060 { 00061 TString branchName, branchVariables; 00062 for ( int i=0; i<5; i++ ) { 00063 branchName = prefix + type + '_' + *parameterNames[i]; 00064 branchVariables = branchName + '[' + index + "]/F"; 00065 vertexTree->Branch(branchName,(void *)(result->*pfunc)(i),branchVariables); 00066 } 00067 }
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().
00168 { 00169 ++numberOfVertices; 00170 00171 TDirectory* rootDir = gDirectory; 00172 // 00173 // fill entry 00174 // 00175 static int nFill(0); 00176 vertexTree->GetDirectory()->cd(); 00177 vertexTree->Fill(); 00178 if ( (++nFill)%1000==0 ) vertexTree->AutoSave(); 00179 // 00180 // restore directory 00181 // 00182 rootDir->cd(); 00183 result->reset(); 00184 00185 00186 }
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.
00147 { 00148 result->fill(TransientVertex(), simv, 0); 00149 fill(); 00150 }
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 | ( | const TransientVertex & | recv, | |
const TrackingVertex * | simv = 0 , |
|||
const float & | time = 0. | |||
) |
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().
Definition at line 108 of file SimpleVertexTree.h.
float SimpleVertexTree::chiProb [private] |
float SimpleVertexTree::chiTot [private] |
Definition at line 99 of file SimpleVertexTree.h.
int SimpleVertexTree::maxTrack [private] |
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] |
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().