00001 #include "RecoVertex/KalmanVertexFit/interface/SimpleVertexTree.h"
00002
00003 #include "TROOT.h"
00004 #include "TTree.h"
00005 #include "TH1F.h"
00006 #include "iostream"
00007
00008 using namespace std;
00009
00010 SimpleVertexTree::SimpleVertexTree(const char * filterName,
00011 TrackAssociatorByChi2 * associator) :
00012 theFitterName(filterName), associatorForParamAtPca(associator)
00013 {
00014
00015 vertexTree = new TTree(filterName, "Vertex fit results");
00016
00017
00018
00019
00020
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 }
00055
00056
00057 void SimpleVertexTree::defineTrackBranch(const TString& prefix, const TString& type,
00058 const float* (VertexFitterResult::*pfunc)(const int) const,
00059 const TString& index)
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 }
00068
00069
00070 SimpleVertexTree::~SimpleVertexTree()
00071 {
00072 std::cout << std::endl<< "End of SimpleVertexTree for "<< theFitterName << std::endl;
00073 std::cout << std::endl<< "Number of vertices fit: "<< numberOfVertices<<std::endl;
00074
00075
00076
00077
00078 TDirectory* rootDir = gDirectory;
00079
00080
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 std::cout << "Mean of Residual distribution X: "<< resX->GetMean()<<std::endl;
00102 std::cout << "Mean of Residual distribution Y: "<< resY->GetMean()<<std::endl;
00103 std::cout << "Mean of Residual distribution Z: "<< resZ->GetMean()<<std::endl;
00104 std::cout << "RMS of Residual distribution X: "<< resX->GetRMS()<<std::endl;
00105 std::cout << "RMS of Residual distribution Y: "<< resY->GetRMS()<<std::endl;
00106 std::cout << "RMS of Residual distribution Z: "<< resZ->GetRMS()<<std::endl;
00107 std::cout << "Mean of Pull distribution X: "<< pullX->GetMean()<<std::endl;
00108 std::cout << "Mean of Pull distribution Y: "<< pullY->GetMean()<<std::endl;
00109 std::cout << "Mean of Pull distribution Z: "<< pullZ->GetMean()<<std::endl;
00110 std::cout << "RMS of Pull distribution X: "<< pullX->GetRMS()<<std::endl;
00111 std::cout << "RMS of Pull distribution Y: "<< pullY->GetRMS()<<std::endl;
00112 std::cout << "RMS of Pull distribution Z: "<< pullZ->GetRMS()<<std::endl;
00113 std::cout << "Average chi-square probability: "<< chiProb->GetMean()<<std::endl;
00114 std::cout << "Average normalized chi-square : "<< chiNorm->GetMean()<<std::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
00127
00128 rootDir->cd();
00129 std::cout << std::endl;
00130
00131 }
00132
00133 void SimpleVertexTree::fill(const TransientVertex & recv, const TrackingVertex * simv,
00134 reco::RecoToSimCollection *recSimColl, const float &time)
00135 {
00136 result->fill(recv, simv, recSimColl, time);
00137 fill();
00138 }
00139
00140 void SimpleVertexTree::fill(const TransientVertex & recv, const float &time)
00141 {
00142 result->fill(recv, 0, 0, time);
00143 fill();
00144 }
00145
00146 void SimpleVertexTree::fill(const TrackingVertex * simv)
00147 {
00148 result->fill(TransientVertex(), simv, 0);
00149 fill();
00150 }
00151
00152
00153
00154
00155
00156
00157
00158
00159
00160
00161
00162
00163
00164
00165
00166
00167 void SimpleVertexTree::fill()
00168 {
00169 ++numberOfVertices;
00170
00171 TDirectory* rootDir = gDirectory;
00172
00173
00174
00175 static int nFill(0);
00176 vertexTree->GetDirectory()->cd();
00177 vertexTree->Fill();
00178 if ( (++nFill)%1000==0 ) vertexTree->AutoSave();
00179
00180
00181
00182 rootDir->cd();
00183 result->reset();
00184
00185
00186 }