CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/DQMOffline/Alignment/src/TkAlCaRecoMonitor.cc

Go to the documentation of this file.
00001 /*
00002  *  See header file for a description of this class.
00003  *
00004  */
00005 
00006 #include "FWCore/Framework/interface/ESHandle.h"
00007 
00008 #include "DataFormats/TrackReco/interface/Track.h"
00009 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00010 #include "MagneticField/Engine/interface/MagneticField.h"
00011 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00012 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00013 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00014 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
00015 
00016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00017 #include "FWCore/Utilities/interface/InputTag.h"
00018 #include "DQMServices/Core/interface/DQMStore.h"
00019 #include "DQMOffline/Alignment/interface/TkAlCaRecoMonitor.h"
00020 
00021 #include <DataFormats/JetReco/interface/CaloJet.h>
00022 #include <DataFormats/Math/interface/deltaR.h>
00023 
00024 #include <string>
00025 //#include <sstream>
00026 #include "TLorentzVector.h"
00027 
00028 TkAlCaRecoMonitor::TkAlCaRecoMonitor(const edm::ParameterSet& iConfig) {
00029   dqmStore_ = edm::Service<DQMStore>().operator->();
00030   conf_ = iConfig;
00031 }
00032 
00033 TkAlCaRecoMonitor::~TkAlCaRecoMonitor() { } 
00034 
00035 void TkAlCaRecoMonitor::beginJob() {
00036 
00037   std::string histname;  //for naming the histograms according to algorithm used
00038 
00039   trackProducer_ = conf_.getParameter<edm::InputTag>("TrackProducer");
00040   referenceTrackProducer_ = conf_.getParameter<edm::InputTag>("ReferenceTrackProducer");
00041 
00042   std::string AlgoName     = conf_.getParameter<std::string>("AlgoName");
00043   std::string MEFolderName = conf_.getParameter<std::string>("FolderName"); 
00044 
00045   daughterMass_ = conf_.getParameter<double>("daughterMass");
00046 
00047   maxJetPt_ = conf_.getParameter<double>("maxJetPt");
00048 
00049   dqmStore_->setCurrentFolder(MEFolderName+"/TkAlignmentSpecific");
00050   fillInvariantMass_ = conf_.getParameter<bool>("fillInvariantMass");
00051   runsOnReco_ = conf_.getParameter<bool>("runsOnReco");
00052   useSignedR_ = conf_.getParameter<bool>("useSignedR");
00053   fillRawIdMap_ = conf_.getParameter<bool>("fillRawIdMap");
00054 
00055   //    
00056   unsigned int MassBin = conf_.getParameter<unsigned int>("MassBin");
00057   double MassMin = conf_.getParameter<double>("MassMin");
00058   double MassMax = conf_.getParameter<double>("MassMax");
00059 
00060   if(fillInvariantMass_){
00061     histname = "InvariantMass_";
00062     invariantMass_ = dqmStore_->book1D(histname+AlgoName, histname+AlgoName, MassBin, MassMin, MassMax);
00063     invariantMass_->setAxisTitle("invariant Mass / GeV");
00064   }else{
00065     invariantMass_ = 0;
00066   }
00067 
00068   unsigned int TrackPtPositiveBin = conf_.getParameter<unsigned int>("TrackPtBin");
00069   double TrackPtPositiveMin = conf_.getParameter<double>("TrackPtMin");
00070   double TrackPtPositiveMax = conf_.getParameter<double>("TrackPtMax");
00071 
00072   histname = "TrackPtPositive_";
00073   TrackPtPositive_ = dqmStore_->book1D(histname+AlgoName, histname+AlgoName, TrackPtPositiveBin, TrackPtPositiveMin, TrackPtPositiveMax);
00074   TrackPtPositive_->setAxisTitle("p_{T} of tracks charge > 0");
00075 
00076   unsigned int TrackPtNegativeBin = conf_.getParameter<unsigned int>("TrackPtBin");
00077   double TrackPtNegativeMin = conf_.getParameter<double>("TrackPtMin");
00078   double TrackPtNegativeMax = conf_.getParameter<double>("TrackPtMax");
00079 
00080   histname = "TrackPtNegative_";
00081   TrackPtNegative_ = dqmStore_->book1D(histname+AlgoName, histname+AlgoName, TrackPtNegativeBin, TrackPtNegativeMin, TrackPtNegativeMax);
00082   TrackPtNegative_->setAxisTitle("p_{T} of tracks charge < 0");
00083 
00084   histname = "TrackQuality_";
00085   TrackQuality_ = dqmStore_->book1D(histname+AlgoName, histname+AlgoName, 
00086                                  reco::TrackBase::qualitySize, -0.5, reco::TrackBase::qualitySize-0.5);
00087   TrackQuality_->setAxisTitle("quality");
00088   for ( int i = 0; i<reco::TrackBase::qualitySize; ++i){
00089     TrackQuality_->getTH1()->GetXaxis()->SetBinLabel(i+1,
00090                     reco::TrackBase::qualityName( reco::TrackBase::TrackQuality(i) ).c_str());
00091   } 
00092 
00093   unsigned int SumChargeBin = conf_.getParameter<unsigned int>("SumChargeBin");
00094   double SumChargeMin = conf_.getParameter<double>("SumChargeMin");
00095   double SumChargeMax = conf_.getParameter<double>("SumChargeMax");
00096 
00097   histname = "SumCharge_";
00098   sumCharge_ = dqmStore_->book1D(histname+AlgoName, histname+AlgoName, SumChargeBin, SumChargeMin, SumChargeMax);
00099   sumCharge_->setAxisTitle("#SigmaCharge");
00100 
00101   unsigned int TrackCurvatureBin = conf_.getParameter<unsigned int>("TrackCurvatureBin");
00102   double TrackCurvatureMin = conf_.getParameter<double>("TrackCurvatureMin");
00103   double TrackCurvatureMax = conf_.getParameter<double>("TrackCurvatureMax");
00104   
00105   histname = "TrackCurvature_";
00106   TrackCurvature_ = dqmStore_->book1D(histname+AlgoName, histname+AlgoName, TrackCurvatureBin, TrackCurvatureMin, TrackCurvatureMax);
00107   TrackCurvature_->setAxisTitle("#kappa track");
00108 
00109   if( runsOnReco_ ){
00110 
00111     unsigned int    JetPtBin = conf_.getParameter<unsigned int>("JetPtBin");
00112     double JetPtMin = conf_.getParameter<double>("JetPtMin");
00113     double JetPtMax = conf_.getParameter<double>("JetPtMax");
00114 
00115     histname = "JetPt_";
00116     jetPt_ = dqmStore_->book1D(histname+AlgoName, histname+AlgoName, JetPtBin, JetPtMin, JetPtMax);
00117     jetPt_->setAxisTitle("jet p_{T} / GeV");
00118 
00119     unsigned int    MinJetDeltaRBin = conf_.getParameter<unsigned int>("MinJetDeltaRBin");
00120     double MinJetDeltaRMin = conf_.getParameter<double>("MinJetDeltaRMin");
00121     double MinJetDeltaRMax = conf_.getParameter<double>("MinJetDeltaRMax");
00122 
00123     histname = "MinJetDeltaR_";
00124     minJetDeltaR_ = dqmStore_->book1D(histname+AlgoName, histname+AlgoName, MinJetDeltaRBin, MinJetDeltaRMin, MinJetDeltaRMax);
00125     minJetDeltaR_->setAxisTitle("minimal Jet #DeltaR / rad");
00126   }else{
00127     jetPt_ = NULL;
00128     minJetDeltaR_ = NULL;
00129   }
00130 
00131   unsigned int    MinTrackDeltaRBin = conf_.getParameter<unsigned int>("MinTrackDeltaRBin");
00132   double MinTrackDeltaRMin = conf_.getParameter<double>("MinTrackDeltaRMin");
00133   double MinTrackDeltaRMax = conf_.getParameter<double>("MinTrackDeltaRMax");
00134 
00135   histname = "MinTrackDeltaR_";
00136   minTrackDeltaR_ = dqmStore_->book1D(histname+AlgoName, histname+AlgoName, MinTrackDeltaRBin, MinTrackDeltaRMin, MinTrackDeltaRMax);
00137   minTrackDeltaR_->setAxisTitle("minimal Track #DeltaR / rad");
00138 
00139   unsigned int TrackEfficiencyBin = conf_.getParameter<unsigned int>("TrackEfficiencyBin");
00140   double TrackEfficiencyMin = conf_.getParameter<double>("TrackEfficiencyMin");
00141   double TrackEfficiencyMax = conf_.getParameter<double>("TrackEfficiencyMax");
00142 
00143   histname = "AlCaRecoTrackEfficiency_";
00144   AlCaRecoTrackEfficiency_ = dqmStore_->book1D(histname+AlgoName, histname+AlgoName, TrackEfficiencyBin, TrackEfficiencyMin, TrackEfficiencyMax);
00145   AlCaRecoTrackEfficiency_->setAxisTitle("n("+trackProducer_.label()+") / n("+referenceTrackProducer_.label()+")");
00146 
00147   int zBin =  conf_.getParameter<unsigned int>("HitMapsZBin"); //300
00148   double zMax = conf_.getParameter<double>("HitMapZMax"); //300.0; //cm
00149   
00150   int rBin = conf_.getParameter<unsigned int>("HitMapsRBin");//120;
00151   double rMax = conf_.getParameter<double>("HitMapRMax"); //120.0; //cm
00152 
00153   histname = "Hits_ZvsR_";
00154   double rMin = 0.0;
00155   if( useSignedR_ )
00156     rMin = -rMax;
00157 
00158   Hits_ZvsR_ = dqmStore_->book2D(histname+AlgoName, histname+AlgoName, zBin, -zMax, zMax, rBin, rMin, rMax);
00159 
00160   histname = "Hits_XvsY_";
00161   Hits_XvsY_ = dqmStore_->book2D(histname+AlgoName, histname+AlgoName, rBin, -rMax, rMax, rBin, -rMax, rMax);
00162 
00163   if( fillRawIdMap_){
00164     histname = "Hits_perDetId_";
00165     
00166     //leads to differences in axsis between samples??
00167     //int nModules = binByRawId_.size();
00168     //Hits_perDetId_ = dqmStore_->book1D(histname+AlgoName, histname+AlgoName, nModules, static_cast<double>(nModules) -0.5, static_cast<double>(nModules) -0.5);
00169     Hits_perDetId_ = dqmStore_->book1D(histname+AlgoName, histname+AlgoName, 16601,-0.5, 16600.5 );
00170     Hits_perDetId_->setAxisTitle("rawId Bins");
00171 
00173     //  std::stringstream binLabel;
00174     //  for( std::map<int,int>::iterator it = binByRawId_.begin(); it != binByRawId_.end(); ++it ){
00175     //    binLabel.str() = "";
00176     //    binLabel << (*it).first;
00177     //    Hits_perDetId_->getTH1()->GetXaxis()->SetBinLabel( (*it).second +1, binLabel.str().c_str());
00178     //  }
00179   }
00180 }
00181 //
00182 // -- Analyse
00183 //
00184 void TkAlCaRecoMonitor::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
00185 
00186   edm::Handle<reco::TrackCollection> trackCollection;
00187   iEvent.getByLabel(trackProducer_, trackCollection);
00188   if (!trackCollection.isValid()){
00189     edm::LogError("Alignment")<<"invalid trackcollection encountered!";
00190     return;
00191   }
00192 
00193   edm::Handle<reco::TrackCollection> referenceTrackCollection;
00194   iEvent.getByLabel(referenceTrackProducer_, referenceTrackCollection);
00195   if (!trackCollection.isValid()){
00196     edm::LogError("Alignment")<<"invalid reference track-collection encountered!";
00197     return;
00198   }
00199 
00200   edm::ESHandle<TrackerGeometry> geometry;
00201   iSetup.get<TrackerDigiGeometryRecord>().get(geometry);
00202   if(! geometry.isValid()){
00203     edm::LogError("Alignment")<<"invalid geometry found in event setup!";
00204   }
00205 
00206   edm::ESHandle<MagneticField> magneticField;
00207   iSetup.get<IdealMagneticFieldRecord>().get(magneticField);
00208   if (!magneticField.isValid()){
00209     edm::LogError("Alignment")<<"invalid magnetic field configuration encountered!";
00210     return;
00211   }
00212 
00213   edm::Handle<reco::CaloJetCollection> jets;
00214   if(runsOnReco_){
00215     edm::InputTag jetCollection = conf_.getParameter<edm::InputTag>("CaloJetCollection");
00216     iEvent.getByLabel(jetCollection, jets);
00217     if(! jets.isValid()){
00218       edm::LogError("Alignment")<<"no jet collection found in event!";
00219     }
00220   }
00221   // fill only once - not yet in beginJob since no access to geometry
00222   if (fillRawIdMap_ && binByRawId_.empty()) this->fillRawIdMap(*geometry);
00223 
00224   AlCaRecoTrackEfficiency_->Fill( static_cast<double>((*trackCollection).size()) / (*referenceTrackCollection).size() );
00225   
00226   double sumOfCharges = 0;
00227   for( reco::TrackCollection::const_iterator track = (*trackCollection).begin(); track < (*trackCollection).end(); ++track ){
00228     double dR = 0;  
00229     if(runsOnReco_){
00230       double minJetDeltaR = 10; // some number > 2pi
00231       for(reco::CaloJetCollection::const_iterator itJet = jets->begin(); itJet != jets->end() ; ++itJet) {
00232         jetPt_->Fill( (*itJet).pt() );
00233         dR = deltaR((*track),(*itJet));
00234         if((*itJet).pt() > maxJetPt_ && dR < minJetDeltaR)
00235           minJetDeltaR = dR; 
00236         
00237         //edm::LogInfo("Alignment") <<">  isolated: "<< isolated << " jetPt "<< (*itJet).pt() <<" deltaR: "<< deltaR(*(*it),(*itJet)) ;
00238       }
00239       minJetDeltaR_->Fill( minJetDeltaR );
00240     }
00241     
00242     double minTrackDeltaR = 10; // some number > 2pi
00243     for( reco::TrackCollection::const_iterator track2 = (*trackCollection).begin(); track2 < (*trackCollection).end(); ++track2 ){
00244       dR = deltaR((*track),(*track2));
00245       if(dR < minTrackDeltaR && dR > 1e-6)
00246         minTrackDeltaR = dR;
00247     }
00248 
00249     for ( int i = 0; i<reco::TrackBase::qualitySize; ++i){
00250       if( (*track).quality( reco::TrackBase::TrackQuality(i) ) ){
00251         TrackQuality_->Fill( i );
00252       }
00253     }
00254     
00255     GlobalPoint gPoint((*track).vx(), (*track).vy(), (*track).vz());
00256     double B = magneticField->inTesla(gPoint).z();
00257     double curv = -(*track).charge()*0.002998*B/(*track).pt();
00258     TrackCurvature_->Fill( curv );
00259 
00260     if( (*track).charge() > 0 )
00261       TrackPtPositive_->Fill( (*track).pt() );
00262     if( (*track).charge() < 0 )
00263       TrackPtNegative_->Fill( (*track).pt() );
00264 
00265     minTrackDeltaR_->Fill( minTrackDeltaR );
00266     fillHitmaps( *track, *geometry );
00267     sumOfCharges += (*track).charge();
00268   }
00269 
00270   sumCharge_->Fill( sumOfCharges );
00271 
00272   if(fillInvariantMass_){
00273     if((*trackCollection).size() == 2){
00274       TLorentzVector track0((*trackCollection).at(0).px(),(*trackCollection).at(0).py(),(*trackCollection).at(0).pz(),
00275                             sqrt(((*trackCollection).at(0).p()*(*trackCollection).at(0).p())+daughterMass_*daughterMass_));
00276       TLorentzVector track1((*trackCollection).at(1).px(),(*trackCollection).at(1).py(),(*trackCollection).at(1).pz(),
00277                             sqrt(((*trackCollection).at(1).p()*(*trackCollection).at(1).p())+daughterMass_*daughterMass_));
00278       TLorentzVector mother = track0+track1;
00279       
00280       invariantMass_->Fill( mother.M() );
00281     }else{
00282       edm::LogInfo("Alignment")<<"wrong number of tracks trackcollection encountered: "<<(*trackCollection).size();
00283     }
00284   }
00285 }
00286 
00287 void TkAlCaRecoMonitor::fillHitmaps(const reco::Track &track, const TrackerGeometry& geometry)
00288 {
00289   for (trackingRecHit_iterator iHit = track.recHitsBegin(); iHit != track.recHitsEnd(); ++iHit) {    
00290     if( (*iHit)->isValid() ){
00291       const TrackingRecHit *hit = (*iHit).get();
00292       const DetId geoId(hit->geographicalId());
00293       const GeomDet* gd = geometry.idToDet(geoId);
00294       // since 2_1_X local hit positions are transient. taking center of the hit module for now.
00295       // The alternative would be the coarse estimation or a refit.
00296       //const GlobalPoint globP( gd->toGlobal( hit->localPosition() ) );
00297       const GlobalPoint globP( gd->toGlobal( Local3DPoint(0.,0.,0.) ) );
00298       double r = sqrt( globP.x()*globP.x() + globP.y()*globP.y() );
00299       if( useSignedR_ )
00300         r*= globP.y() / fabs( globP.y() );
00301 
00302       Hits_ZvsR_->Fill( globP.z(), r );
00303       Hits_XvsY_->Fill( globP.x(), globP.y() );
00304 
00305       if( fillRawIdMap_)
00306          Hits_perDetId_->Fill( binByRawId_[ geoId.rawId() ]);  
00307     }
00308   }
00309 }
00310 
00311 void TkAlCaRecoMonitor::fillRawIdMap(const TrackerGeometry &geometry)
00312 {
00313   std::vector<int> sortedRawIds;
00314   for (std::vector<DetId>::const_iterator iDetId = geometry.detUnitIds().begin();
00315        iDetId != geometry.detUnitIds().end(); ++iDetId) {
00316     sortedRawIds.push_back((*iDetId).rawId());
00317   }
00318   std::sort(sortedRawIds.begin(), sortedRawIds.end());
00319   
00320   int i = 0;
00321   for (std::vector<int>::iterator iRawId = sortedRawIds.begin();
00322        iRawId != sortedRawIds.end(); ++iRawId){
00323     binByRawId_[*iRawId] = i;
00324     ++i;
00325   }
00326 }
00327 
00328 
00329 void TkAlCaRecoMonitor::endJob(void) {
00330   bool outputMEsInRootFile = conf_.getParameter<bool>("OutputMEsInRootFile");
00331   std::string outputFileName = conf_.getParameter<std::string>("OutputFileName");
00332   if(outputMEsInRootFile){
00333     //dqmStore_->showDirStructure();
00334     dqmStore_->save(outputFileName);
00335   }
00336 }
00337 
00338 DEFINE_FWK_MODULE(TkAlCaRecoMonitor);
00339