00001
00002
00003
00004
00005
00016
00017
00018
00019
00020
00021
00022
00023
00024 #include <memory>
00025
00026
00027
00028 #include "Alignment/HIPAlignmentAlgorithm/interface/LhcTrackAnalyzer.h"
00029 #include "FWCore/Framework/interface/Frameworkfwd.h"
00030 #include "FWCore/Framework/interface/EDAnalyzer.h"
00031
00032 #include "SimDataFormats/EncodedEventId/interface/EncodedEventId.h"
00033 #include "FWCore/Framework/interface/Event.h"
00034 #include <FWCore/Framework/interface/ESHandle.h>
00035 #include "FWCore/Framework/interface/MakerMacros.h"
00036
00037 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00038
00039
00040 #include "TH1F.h"
00041 #include "TH2F.h"
00042 #include "TFile.h"
00043 #include "TROOT.h"
00044 #include "TChain.h"
00045 #include "TNtuple.h"
00046 #include "FWCore/ServiceRegistry/interface/Service.h"
00047 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00048 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00049 #include <Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h>
00050 #include <DataFormats/GeometrySurface/interface/Surface.h>
00051 #include <DataFormats/GeometrySurface/interface/GloballyPositioned.h>
00052 #include <Geometry/CommonDetUnit/interface/GeomDet.h>
00053 #include "MagneticField/Engine/interface/MagneticField.h"
00054 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00055
00056 #include "DataFormats/TrackReco/interface/Track.h"
00057 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00058 #include "DataFormats/VertexReco/interface/Vertex.h"
00059 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00060
00061
00062
00063 LhcTrackAnalyzer::LhcTrackAnalyzer(const edm::ParameterSet& iConfig)
00064
00065 {
00066
00067 debug_ = iConfig.getParameter<bool> ("Debug");
00068 TrackCollectionTag_ = iConfig.getParameter<edm::InputTag>("TrackCollectionTag");
00069 PVtxCollectionTag_ = iConfig.getParameter<edm::InputTag>("PVtxCollectionTag");
00070 filename_ = iConfig.getParameter<std::string>("OutputFileName");
00071 }
00072
00073
00074 LhcTrackAnalyzer::~LhcTrackAnalyzer()
00075 {
00076
00077
00078
00079
00080 }
00081
00082
00083
00084
00085
00086
00087
00088 void
00089 LhcTrackAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00090 {
00091 using namespace edm;
00092 using namespace reco;
00093 using namespace std;
00094
00095
00096
00097
00098
00099 SetVarToZero();
00100
00101
00102
00103
00104
00105 Handle< TrackCollection> trackCollectionHandle;
00106 iEvent.getByLabel(TrackCollectionTag_, trackCollectionHandle);
00107 Handle<VertexCollection> vertexCollectionHandle;
00108 iEvent.getByLabel(PVtxCollectionTag_, vertexCollectionHandle);
00109 for(VertexCollection::const_iterator vtx = vertexCollectionHandle->begin();vtx!=vertexCollectionHandle->end(); ++vtx)
00110 {
00111 if(vtx==vertexCollectionHandle->begin()){
00112 if(vtx->isFake())goodvtx_=false;
00113 else goodvtx_=true;
00114 }
00115 else break;
00116 }
00117
00118
00119
00120 goodbx_=true;
00121
00122
00123
00124
00125
00126 run_=iEvent.id().run();
00127 event_=iEvent.id().event();
00128
00129 if(debug_)
00130 cout<<"LhcTrackAnalyzer::analyze() looping over "<< trackCollectionHandle->size()<< "tracks." << endl;
00131
00132
00133 bool toomanytracks=false;
00134 for(TrackCollection::const_iterator track = trackCollectionHandle->begin(); track!= trackCollectionHandle->end(); ++track)
00135 {
00136
00137 if ( nTracks_ >= nMaxtracks_ ) {
00138 std::cout << " LhcTrackAnalyzer::analyze() : Warning - Run "<< run_<<" Event "<< event_<<"\tNumber of tracks: " << trackCollectionHandle->size()<< " , greater than " << nMaxtracks_ << std::endl;
00139 toomanytracks=true;
00140 continue;
00141
00142 }
00143
00144 pt_[nTracks_] = track->pt();
00145 eta_[nTracks_] = track->eta();
00146 phi_[nTracks_] = track->phi();
00147 chi2_[nTracks_] = track->chi2();
00148 chi2ndof_[nTracks_] = track->normalizedChi2();
00149 charge_[nTracks_] = track->charge();
00150 qoverp_[nTracks_] = track->qoverp();
00151 dz_[nTracks_] = track->dz();
00152 dxy_[nTracks_] = track->dxy();
00153 xPCA_[nTracks_] = track->vertex().x();
00154 yPCA_[nTracks_] = track->vertex().y();
00155 zPCA_[nTracks_] = track->vertex().z();
00156 validhits_[nTracks_][0]=track->numberOfValidHits();
00157 validhits_[nTracks_][1]=track->hitPattern().numberOfValidPixelBarrelHits();
00158 validhits_[nTracks_][2]=track->hitPattern().numberOfValidPixelEndcapHits();
00159 validhits_[nTracks_][3]=track->hitPattern().numberOfValidStripTIBHits();
00160 validhits_[nTracks_][4]=track->hitPattern().numberOfValidStripTIDHits();
00161 validhits_[nTracks_][5]=track->hitPattern().numberOfValidStripTOBHits();
00162 validhits_[nTracks_][6]=track->hitPattern().numberOfValidStripTECHits();
00163
00164
00165
00166 int myalgo=-88;
00167 if(track->algo()==reco::TrackBase::undefAlgorithm)myalgo=0;
00168 if(track->algo()==reco::TrackBase::ctf)myalgo=1;
00169 if(track->algo()==reco::TrackBase::iter0)myalgo=4;
00170 if(track->algo()==reco::TrackBase::iter1)myalgo=5;
00171 if(track->algo()==reco::TrackBase::iter2)myalgo=6;
00172 if(track->algo()==reco::TrackBase::iter3)myalgo=7;
00173 if(track->algo()==reco::TrackBase::iter4)myalgo=8;
00174 if(track->algo()==reco::TrackBase::iter5)myalgo=9;
00175 if(track->algo()==reco::TrackBase::iter6)myalgo=10;
00176 if(track->algo()==reco::TrackBase::iter7)myalgo=11;
00177 trkAlgo_[nTracks_] = myalgo;
00178
00179 int myquality=-99;
00180 if(track->quality(reco::TrackBase::undefQuality))myquality=-1;
00181 if(track->quality(reco::TrackBase::loose))myquality=0;
00182 if(track->quality(reco::TrackBase::tight))myquality=1;
00183 if(track->quality(reco::TrackBase::highPurity))myquality=2;
00184
00185
00186
00187 trkQuality_[nTracks_]= myquality;
00188
00189 if(track->quality(reco::TrackBase::highPurity))isHighPurity_[nTracks_]=1;
00190 else isHighPurity_[nTracks_]=0;
00191 nTracks_++;
00192
00193
00194
00195
00196
00197
00198
00199
00200 }
00201
00202 for(int d=0;d<nTracks_;++d){
00203 if(abs(trkQuality_[d])>5)cout<<"MYQUALITY!!! " <<trkQuality_[d] <<" at track # "<<d<<"/"<< nTracks_<<endl;
00204 }
00205
00206
00207
00208
00209 rootTree_->Fill();
00210 }
00211
00212
00213
00214 void LhcTrackAnalyzer::beginJob()
00215 {
00216 edm::LogInfo("beginJob") << "Begin Job" << std::endl;
00217
00218 rootFile_ = new TFile(filename_.c_str(),"recreate");
00219 rootTree_ = new TTree("tree","Lhc Track tree");
00220
00221
00222 rootTree_->Branch("run",&run_,"run/I");
00223 rootTree_->Branch("event",&event_,"event/I");
00224 rootTree_->Branch("goodbx",&goodbx_,"goodbx/O");
00225 rootTree_->Branch("goodvtx",&goodvtx_,"goodvtx/O");
00226 rootTree_->Branch("nTracks",&nTracks_,"nTracks/I");
00227 rootTree_->Branch("pt",&pt_,"pt[nTracks]/D");
00228 rootTree_->Branch("eta",&eta_,"eta[nTracks]/D");
00229 rootTree_->Branch("phi",&phi_,"phi[nTracks]/D");
00230 rootTree_->Branch("chi2",&chi2_,"chi2[nTracks]/D");
00231 rootTree_->Branch("chi2ndof",&chi2ndof_,"chi2ndof[nTracks]/D");
00232 rootTree_->Branch("charge",&charge_,"charge[nTracks]/I");
00233 rootTree_->Branch("qoverp",&qoverp_,"qoverp[nTracks]/D");
00234 rootTree_->Branch("dz",&dz_,"dz[nTracks]/D");
00235 rootTree_->Branch("dxy",&dxy_,"dxy[nTracks]/D");
00236 rootTree_->Branch("xPCA",&xPCA_,"xPCA[nTracks]/D");
00237 rootTree_->Branch("yPCA",&yPCA_,"yPCA[nTracks]/D");
00238 rootTree_->Branch("zPCA",&zPCA_,"zPCA[nTracks]/D");
00239 rootTree_->Branch("isHighPurity",&isHighPurity_,"isHighPurity[nTracks]/I");
00240 rootTree_->Branch("trkQuality",&trkQuality_,"trkQuality[nTracks]/I");
00241 rootTree_->Branch("trkAlgo",&trkAlgo_,"trkAlgo[nTracks]/I");
00242 rootTree_->Branch("nValidHits",&validhits_,"nValidHits[nTracks][7]/I");
00243
00244
00245 }
00246
00247
00248 void LhcTrackAnalyzer::endJob()
00249 {
00250 if ( rootFile_ ) {
00251 rootFile_->Write();
00252 rootFile_->Close();
00253 }
00254
00255
00256
00257 }
00258
00259 void LhcTrackAnalyzer::SetVarToZero() {
00260 run_=-1;
00261 event_=-99;
00262 nTracks_ = 0;
00263 for ( int i=0; i<nMaxtracks_; ++i ) {
00264 pt_[i] = 0;
00265 eta_[i] = 0;
00266 phi_[i] = 0;
00267 chi2_[i] = 0;
00268 chi2ndof_[i] = 0;
00269 charge_[i] = 0;
00270 qoverp_[i] = 0;
00271 dz_[i] = 0;
00272 dxy_[i] = 0;
00273 xPCA_[i] = 0;
00274 yPCA_[i] = 0;
00275 zPCA_[i] = 0;
00276 trkQuality_[i] = 0;
00277 trkAlgo_[i] = -1;
00278 isHighPurity_[i]=-3;
00279 for(int j=0;j<7;j++){
00280 validhits_[nTracks_][j]=-1*j;
00281 }
00282 }
00283
00284
00285 }
00286
00287
00288 DEFINE_FWK_MODULE(LhcTrackAnalyzer);