CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/DQM/Physics/src/QcdUeDQM.cc

Go to the documentation of this file.
00001 
00002 /*
00003     This is the DQM code for UE physics plots
00004     11/12/2009 Sunil Bansal
00005 */
00006 #include "DQM/Physics/src/QcdUeDQM.h"
00007 #include "DataFormats/Common/interface/TriggerResults.h"
00008 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00009 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
00010 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
00011 #include "DataFormats/SiPixelDetId/interface/PXFDetId.h"
00012 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHitCollection.h"
00013 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00014 #include "DataFormats/VertexReco/interface/Vertex.h"
00015 #include "DQMServices/Core/interface/DQMStore.h"
00016 #include "DQMServices/Core/interface/MonitorElement.h"
00017 #include "FWCore/Framework/interface/Event.h"
00018 #include "FWCore/Framework/interface/EventSetup.h"
00019 #include "FWCore/Framework/interface/ESHandle.h"
00020 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00022 #include "FWCore/ServiceRegistry/interface/Service.h"
00023 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
00024 #include "CommonTools/RecoAlgos/src/TrackToRefCandidate.h"
00025 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00026 #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetType.h"
00027 #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetUnit.h"
00028 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00029 #include "Geometry/TrackerGeometryBuilder/interface/RectangularPixelTopology.h" 
00030 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
00031 #include "AnalysisDataFormats/TrackInfo/interface/TrackInfo.h"
00032 #include "AnalysisDataFormats/TrackInfo/interface/TrackInfoTrackAssociation.h"
00033 #include "DataFormats/TrackReco/interface/Track.h"
00034 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00035 #include "DataFormats/Candidate/interface/CandMatchMap.h"
00036 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00037 #include "DataFormats/Math/interface/deltaPhi.h"
00038 #include <TString.h>
00039 #include <TMath.h>
00040 #include <TH1F.h>
00041 #include <TH2F.h>
00042 #include <TH3F.h>
00043 #include <TProfile.h>
00044 using namespace std;
00045 using namespace edm;
00046 
00047 #define CP(level) \
00048   if (level>=verbose_)
00049 
00050 struct deleter {
00051   void operator()(TH3F *&h) { delete h; h=0;}
00052 };
00053 
00054 
00055 //--------------------------------------------------------------------------------------------------
00056 QcdUeDQM::QcdUeDQM(const ParameterSet &parameters) :
00057   hltResName_(parameters.getUntrackedParameter<string>("hltTrgResults")),
00058   verbose_(parameters.getUntrackedParameter<int>("verbose",3)),
00059   tgeo_(0),
00060   theDbe_(0),
00061   repSumMap_(0),
00062   repSummary_(0),
00063   h2TrigCorr_(0),
00064   ptMin_(parameters.getParameter<double>("ptMin")),
00065   minRapidity_(parameters.getParameter<double>("minRapidity")),
00066   maxRapidity_(parameters.getParameter<double>("maxRapidity")),
00067   tip_(parameters.getParameter<double>("tip")),
00068   lip_(parameters.getParameter<double>("lip")),
00069   diffvtxbs_(parameters.getParameter<double>("diffvtxbs")),
00070   ptErr_pt_(parameters.getParameter<double>("ptErr_pt")),
00071   vtxntk_(parameters.getParameter<double>("vtxntk")),
00072   minHit_(parameters.getParameter<int>("minHit")),
00073   pxlLayerMinCut_(parameters.getParameter<double>("pxlLayerMinCut")),
00074   requirePIX1_(parameters.getParameter<bool>("requirePIX1")),
00075   min3DHit_(parameters.getParameter<int>("min3DHit")),
00076   maxChi2_(parameters.getParameter<double>("maxChi2")),
00077   bsuse_(parameters.getParameter<bool>("bsuse")),
00078   allowTriplets_(parameters.getParameter<bool>("allowTriplets")),
00079   bsPos_(parameters.getParameter<double>("bsPos")),
00080   caloJetLabel_(parameters.getUntrackedParameter<edm::InputTag>("caloJetTag")),
00081   chargedJetLabel_(parameters.getUntrackedParameter<edm::InputTag>("chargedJetTag")),
00082   trackLabel_(parameters.getUntrackedParameter<edm::InputTag>("trackTag")),
00083   vtxLabel_(parameters.getUntrackedParameter<edm::InputTag>("vtxTag")),
00084   bsLabel_(parameters.getParameter<edm::InputTag>("beamSpotTag")) 
00085 {
00086   // Constructor.
00087   std::vector<std::string> quality = parameters.getParameter<std::vector<std::string> >("quality");
00088   for (unsigned int j=0;j<quality.size();j++) quality_.push_back(reco::TrackBase::qualityByName(quality[j])); 
00089   std::vector<std::string> algorithm = parameters.getParameter<std::vector<std::string> >("algorithm");
00090   for (unsigned int j=0;j<algorithm.size();j++) algorithm_.push_back(reco::TrackBase::algoByName(algorithm[j])); 
00091 
00092   if (parameters.exists("hltTrgNames"))
00093     hltTrgNames_ = parameters.getUntrackedParameter<vector<string> >("hltTrgNames");
00094 
00095   if (parameters.exists("hltProcNames"))
00096      hltProcNames_ = parameters.getUntrackedParameter<vector<string> >("hltProcNames");
00097   else {
00098     //     hltProcNames_.push_back("FU");
00099      hltProcNames_.push_back("HLT");
00100   }
00101 
00102  
00103 }
00104 
00105 //--------------------------------------------------------------------------------------------------
00106 QcdUeDQM::~QcdUeDQM()
00107 {
00108   // Destructor.
00109 
00110 
00111 }
00112 
00113 //--------------------------------------------------------------------------------------------------
00114 void QcdUeDQM::analyze(const Event &iEvent, const EventSetup &iSetup) 
00115 {
00116   if( ! isHltConfigSuccessful_ ) return;
00117 
00118   // Analyze the given event.
00119    
00120    edm::Handle<reco::BeamSpot> beamSpot;
00121    bool ValidBS_ = iEvent.getByLabel(bsLabel_,beamSpot);
00122    if(!ValidBS_)return;
00123 
00124    edm::Handle<reco::TrackCollection>tracks ;
00125    bool ValidTrack_ = iEvent.getByLabel(trackLabel_,tracks);
00126    if(!ValidTrack_)return;
00127 
00128    //edm::Handle<reco::CandidateView> trkJets;
00129    //bool ValidTrackJet_ = iEvent.getByLabel (chargedJetLabel_,trkJets);
00130    //if(!ValidTrackJet_)return;
00131 
00132    edm::Handle< reco::TrackJetCollection > trkJets ;
00133    bool ValidTrackJet_ =  iEvent.getByLabel(chargedJetLabel_,trkJets );
00134    if(!ValidTrackJet_)return;
00135 
00136    
00137    edm::Handle<reco::CaloJetCollection> calJets;
00138    bool ValidCaloJet_ = iEvent.getByLabel (caloJetLabel_,calJets);
00139    if(!ValidCaloJet_)return;
00140  
00141    edm::Handle< reco::VertexCollection > vertexColl;
00142    bool ValidVtxColl_ = iEvent.getByLabel (vtxLabel_, vertexColl);
00143    if(!ValidVtxColl_)return;
00144 
00145    reco::TrackCollection tracks_sort = *tracks;
00146    std::sort(tracks_sort.begin(), tracks_sort.end(), PtSorter()); 
00147 
00148   // get tracker geometry
00149 /*  ESHandle<TrackerGeometry> trackerHandle;
00150   iSetup.get<TrackerDigiGeometryRecord>().get(trackerHandle);
00151   tgeo_ = trackerHandle.product();
00152   if (!tgeo_)return;
00153 */
00154   selected_.clear(); 
00155   fillHltBits(iEvent,iSetup);
00156   // select good tracks
00157   if(fillVtxPlots(beamSpot.product(),vertexColl))
00158   {
00159   fill1D(hNevts_,1);
00160   for(reco::TrackCollection::const_iterator Trk = tracks_sort.begin(); Trk != tracks_sort.end(); ++Trk)
00161    {
00162     
00163    if ( trackSelection(*Trk,beamSpot.product(),vtx1,vertexColl->size()) ) selected_.push_back( & * Trk );   
00164    }
00165  
00166       
00167     fillpTMaxRelated(selected_);
00168     fillChargedJetSpectra(trkJets);  
00169 //    fillCaloJetSpectra(calJets);
00170     fillUE_with_MaxpTtrack(selected_);
00171     if(trkJets->size() > 0)fillUE_with_ChargedJets(selected_,trkJets); 
00172    // if(calJets->size()>0)fillUE_with_CaloJets(selected_,calJets);
00173   
00174  }
00175 
00176 }
00177 
00178 //--------------------------------------------------------------------------------------------------
00179 void QcdUeDQM::beginJob() 
00180 {
00181   // Begin job and setup the DQM store.
00182 
00183   theDbe_ = Service<DQMStore>().operator->();
00184   if (!theDbe_)return;
00185   
00186   //  theDbe_->setCurrentFolder("Physics/QcdUe");
00187   
00188 }
00189 
00190 //--------------------------------------------------------------------------------------------------
00191 void QcdUeDQM::beginLuminosityBlock(const LuminosityBlock &l, 
00192                                        const EventSetup &iSetup)
00193 {
00194   if( ! isHltConfigSuccessful_ ) return;
00195 
00196   // At the moment, nothing needed to be done.
00197 }
00198 
00199 //--------------------------------------------------------------------------------------------------
00200 void QcdUeDQM::beginRun(const Run &run, const EventSetup &iSetup)
00201 {
00202 
00203  // indicating change of HLT cfg at run boundries
00204  // for HLTConfigProvider::init()
00205  bool isHltCfgChange = false;
00206  isHltConfigSuccessful_ = false; // init
00207 
00208  string teststr;
00209  for(size_t i=0; i<hltProcNames_.size(); ++i) {
00210    if (i>0) 
00211      teststr += ", ";
00212    teststr += hltProcNames_.at(i);
00213    if ( hltConfig.init( run, iSetup, hltProcNames_.at(i), isHltCfgChange ) ) {
00214      isHltConfigSuccessful_ = true;
00215      hltUsedResName_ = hltResName_;
00216      if (hltResName_.find(':')==string::npos)
00217        hltUsedResName_ += "::";
00218      else 
00219        hltUsedResName_ += ":";
00220      hltUsedResName_ += hltProcNames_.at(i);
00221      break;
00222    }
00223  }
00224  
00225  if ( ! isHltConfigSuccessful_ )return;
00226 
00227   // setup "Any" bit
00228   hltTrgBits_.clear();
00229   hltTrgBits_.push_back(-1);
00230   hltTrgDeci_.clear();
00231   hltTrgDeci_.push_back(true);
00232   hltTrgUsedNames_.clear();
00233   hltTrgUsedNames_.push_back("Any");
00234 
00235   // figure out relation of trigger name to trigger bit and store used trigger names/bits
00236   for(size_t i=0;i<hltTrgNames_.size();++i) {
00237     const string &n1(hltTrgNames_.at(i));
00238     //unsigned int hlt_prescale = hltConfig.prescaleValue(iSetup, n1);
00239     //cout<<"trigger=="<<n1<<"presc=="<<hlt_prescale<<endl;
00240     bool found = 0;
00241     for(size_t j=0;j<hltConfig.size();++j) {
00242       const string &n2(hltConfig.triggerName(j));
00243       if (n2==n1) {
00244         hltTrgBits_.push_back(j);
00245         hltTrgUsedNames_.push_back(n1);
00246         hltTrgDeci_.push_back(false);
00247         found = 1;
00248         break;
00249       }
00250     }      
00251     if (!found) {
00252       CP(2) cout<<"Could not find trigger bit"<<endl ;
00253     }
00254   }
00255  
00256   // book monitoring histograms
00257   createHistos();
00258   isHltConfigSuccessful_ = true;
00259 
00260 }
00261 
00262 //--------------------------------------------------------------------------------------------------
00263 void QcdUeDQM::book1D(std::vector<MonitorElement*> &mes, 
00264                          const std::string &name, const std::string &title, 
00265                          int nx, double x1, double x2, bool sumw2, bool sbox)
00266 {
00267   // Book 1D histos.
00268 
00269   for(size_t i=0;i<hltTrgUsedNames_.size();++i) {
00270     std::string folderName = "Physics/QcdUe/" + hltTrgUsedNames_.at(i);
00271     theDbe_->setCurrentFolder(folderName);
00272     MonitorElement *e = theDbe_->book1D(Form("%s_%s",name.c_str(),hltTrgUsedNames_.at(i).c_str()),
00273                                         Form("%s: %s",hltTrgUsedNames_.at(i).c_str(), title.c_str()), 
00274                                         nx, x1, x2);
00275     TH1 *h1 = e->getTH1();
00276     if (sumw2) {
00277       if( 0 == h1->GetSumw2N() ) { // protect against re-summing (would cause exception)
00278         h1->Sumw2();
00279       }
00280     }
00281     h1->SetStats(sbox);
00282     mes.push_back(e);
00283   }
00284 }
00285 
00286 //--------------------------------------------------------------------------------------------------
00287 void QcdUeDQM::book2D(std::vector<MonitorElement*> &mes, 
00288                          const std::string &name, const std::string &title, 
00289                          int nx, double x1, double x2, int ny, double y1, double y2, 
00290                          bool sumw2, bool sbox)
00291 {
00292   // Book 2D histos.
00293 
00294   for(size_t i=0;i<hltTrgUsedNames_.size();++i) {
00295     std::string folderName = "Physics/QcdUe/" + hltTrgUsedNames_.at(i);
00296     theDbe_->setCurrentFolder(folderName);
00297     MonitorElement *e = theDbe_->book2D(Form("%s_%s",name.c_str(),hltTrgUsedNames_.at(i).c_str()),
00298                                         Form("%s: %s",hltTrgUsedNames_.at(i).c_str(), title.c_str()), 
00299                                         nx, x1, x2, ny, y1, y2);
00300     TH1 *h1 = e->getTH1();
00301     if (sumw2) {
00302       if( 0 == h1->GetSumw2N() ) { // protect against re-summing (would cause exception)
00303         h1->Sumw2();
00304       }
00305     }
00306     h1->SetStats(sbox);
00307     mes.push_back(e);
00308   }
00309 }
00310 
00311 //--------------------------------------------------------------------------------------------------
00312 void QcdUeDQM::bookProfile(std::vector<MonitorElement*> &mes, 
00313                          const std::string &name, const std::string &title, 
00314                          int nx, double x1, double x2,  double y1, double y2, 
00315                          bool sumw2, bool sbox)
00316 {
00317   // Book Profile histos.
00318 
00319   for(size_t i=0;i<hltTrgUsedNames_.size();++i) {
00320     std::string folderName = "Physics/QcdUe/" + hltTrgUsedNames_.at(i);
00321     theDbe_->setCurrentFolder(folderName);
00322     MonitorElement *e = theDbe_->bookProfile(Form("%s_%s",name.c_str(),hltTrgUsedNames_.at(i).c_str()),
00323                                         Form("%s: %s",hltTrgUsedNames_.at(i).c_str(), title.c_str()), 
00324                                         nx, x1, x2, y1, y2," ");
00325     mes.push_back(e);
00326   }
00327 }
00328 //--------------------------------------------------------------------------------------------------
00329 void QcdUeDQM::createHistos()
00330 {
00331   // Book histograms if needed.
00332 
00333 
00334 /*  if (1) {
00335     theDbe_->setCurrentFolder("Physics/EventInfo/");
00336     repSumMap_  = theDbe_->book2D("reportSummaryMap","reportSummaryMap",1,0,1,1,0,1);
00337     repSummary_ = theDbe_->bookFloat("reportSummary");
00338   }
00339   */ 
00340    theDbe_->setCurrentFolder("Physics/QcdUe");
00341 
00342   if (1) {
00343     const int Nx = hltTrgUsedNames_.size();
00344     const double x1 = -0.5;
00345     const double x2 = Nx-0.5;
00346     h2TrigCorr_ = theDbe_->book2D("h2TriCorr","Trigger bit x vs y;y&&!x;x&&y",Nx,x1,x2,Nx,x1,x2);
00347     for(size_t i=1;i<=hltTrgUsedNames_.size();++i) {
00348       h2TrigCorr_->setBinLabel(i,hltTrgUsedNames_.at(i-1),1);
00349       h2TrigCorr_->setBinLabel(i,hltTrgUsedNames_.at(i-1),2);
00350     }
00351     TH1 *h = h2TrigCorr_->getTH1();
00352     if (h)
00353       h->SetStats(0);
00354   }
00355   book1D(hNevts_,"hNevts","number of events",2,0,2);
00356   book1D(hNtrackerLayer_,"hNtrackerLayer","number of tracker layers;multiplicity",20,-0.5,19.5 );
00357   book1D(hNtrackerPixelLayer_,"hNtrackerPixelLayer","number of pixel layers;multiplicity",10,-0.5,9.5 );
00358   book1D(hNtrackerStripPixelLayer_,"hNtrackerStripPixelLayer","number of strip + pixel layers;multiplicity",30,-0.5,39.5 );
00359   book1D(hRatioPtErrorPt_,"hRatioPtErrorPt","ratio of pT error and track pT",25,0.,5.);
00360   book1D(hTrkPt_,"hTrkPt","pT of all tracks",50,0.,50.);
00361   book1D(hTrkEta_,"hTrkEta","eta of all tracks",40,-4.,4.);
00362   book1D(hTrkPhi_,"hTrkPhi","phi of all tracks",40,-4.,4.);
00363   book1D(hRatioDxySigmaDxyBS_,"hRatioDxySigmaDxyBS","ratio of transverse impact parameter and its significance wrt beam spot",60,-10.,10);
00364   book1D(hRatioDxySigmaDxyPV_,"hRatioDxySigmaDxyPV","ratio of transverse impact parameter and its significance wrt PV",60,-10.,10);
00365   book1D(hRatioDzSigmaDzBS_,"hRatioDzSigmaDzBS","ratio of longitudinal impact parameter and its significance wrt beam spot",80,-20.,20);
00366   book1D(hRatioDzSigmaDzPV_,"hRatioDzSigmaDzPV","ratio of longitudinal impact parameter and its significance wrt PV",80,-20.,20);
00367   book1D(hTrkChi2_,"hTrkChi2","track chi2",30,0.,30);
00368   book1D(hTrkNdof_,"hTrkNdof","track NDOF",100,0,100);
00369 
00370   book1D(hNgoodTrk_,"hNgoodTrk","number of good tracks",50,-0.5,49.5);
00371 
00372   book1D(hGoodTrkPt500_,"hGoodTrkPt500","pT of all good tracks with pT > 500 MeV",50,0.,50.);
00373   book1D(hGoodTrkEta500_,"hGoodTrkEta500","eta of all good tracks pT > 500 MeV",40,-4.,4.);
00374   book1D(hGoodTrkPhi500_,"hGoodTrkPhi500","phi of all good tracks pT > 500 MeV",40,-4.,4.);
00375 
00376   book1D(hGoodTrkPt900_,"hGoodTrkPt900","pT of all good tracks with pT > 900 MeV",50,0.,50.);
00377   book1D(hGoodTrkEta900_,"hGoodTrkEta900","eta of all good tracks pT > 900 MeV",40,-4.,4.);
00378   book1D(hGoodTrkPhi900_,"hGoodTrkPhi900","phi of all good tracks pT > 900 MeV",40,-4.,4.);
00379 
00380   book1D(hNvertices_,"hNvertices","number of vertices",5,-0.5,4.5);
00381   book1D(hVertex_z_,"hVertex_z","z position of vertex; z[cm]",200,-50,50);
00382   book1D(hVertex_y_,"hVertex_y","y position of vertex; y[cm]",100,-5,5);
00383   book1D(hVertex_x_,"hVertex_x","x position of vertex; x[cm]",100,-5,5);
00384   book1D(hVertex_ndof_,"hVertex_ndof","ndof of vertex",100,0,100);
00385   book1D(hVertex_rho_,"hVertex_rho","rho of vertex",100,0,5);
00386   book1D(hVertex_z_bs_,"hVertex_z_bs","z position of vertex from beamspot; z[cm]",200,-50,50);
00387 
00388   book1D(hBeamSpot_z_,"hBeamSpot_z","z position of beamspot; z[cm]",100,-20,20);
00389   book1D(hBeamSpot_y_,"hBeamSpot_y","y position of beamspot; y[cm]",50,-10,10);
00390   book1D(hBeamSpot_x_,"hBeamSpot_x","x position of beamspot; x[cm]",50,-10,10);
00391 
00392 
00393   if (1) {
00394     const int Nx = 25;
00395     const double x1 = 0.0;
00396     const double x2 = 50.0;
00397     book1D(hLeadingTrack_pTSpectrum_,"hLeadingTrack_pTSpectrum","pT spectrum of leading track;pT(GeV/c)",Nx,x1,x2);
00398 //    book1D(hLeadingCaloJet_pTSpectrum_,"hLeadingCalo_pTSpectrum","pT spectrum of leading calo jet;pT(GeV/c)",Nx,x1,x2);
00399     book1D(hLeadingChargedJet_pTSpectrum_,"hLeadingChargedJet_pTSpectrum","pT spectrum of leading track jet;pT(GeV/c)",Nx,x1,x2);
00400     
00401   }
00402   
00403   if (1) {
00404     const int Nx = 24;
00405     const double x1 = -4.;
00406     const double x2 =  4.;
00407     book1D(hLeadingTrack_phiSpectrum_,"hLeadingTrack_phiSpectrum","#phi spectrum of leading track;#phi",Nx,x1,x2);
00408   //  book1D(hLeadingCaloJet_phiSpectrum_,"hLeadingCaloJet_phiSpectrum","#phi spectrum of leading calo jet;#phi",Nx,x1,x2);
00409     book1D(hLeadingChargedJet_phiSpectrum_,"hLeadingChargedJet_phiSpectrum","#phi spectrum of leading track jet;#phi",Nx,x1,x2);
00410 
00411   }
00412   
00413   if (1) {
00414     const int Nx = 24;
00415     const double x1 = -4.;
00416     const double x2 =  4.;
00417     book1D(hLeadingTrack_etaSpectrum_,"hLeadingTrack_etaSpectrum","#eta spectrum of leading track;#eta",Nx,x1,x2);
00418     //book1D(hLeadingCaloJet_etaSpectrum_,"hLeadingCaloJet_etaSpectrum","#eta spectrum of leading calo jet;#eta",Nx,x1,x2);
00419     book1D(hLeadingChargedJet_etaSpectrum_,"hLeadingChargedJet_etaSpectrum","#eta spectrum of leading track jet;#eta",Nx,x1,x2);
00420 
00421   }
00422 
00423 
00424 if (1) {
00425     const int Nx = 75;
00426     const double x1 = 0.0;
00427     const double x2 = 75.0;
00428     const double y1 = 0.;
00429     const double y2 = 10.;
00430     bookProfile(hdNdEtadPhi_pTMax_Toward500_,"hdNdEtadPhi_pTMax_Toward500", 
00431                  "Average number of tracks (pT > 500 MeV) in toward region vs leading track pT;pT(GeV/c);dN/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00432     bookProfile(hdNdEtadPhi_pTMax_Transverse500_,"hdNdEtadPhi_pTMax_Transverse500", 
00433                  "Average number of tracks (pT > 500 MeV) in transverse region vs leading track pT;pT(GeV/c);dN/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00434     bookProfile(hdNdEtadPhi_pTMax_Away500_,"hdNdEtadPhi_pTMax_Away500", 
00435                  "Average number of tracks (pT > 500 MeV) in away region vs leading track pT;pT(GeV/c);dN/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00436  /*
00437     bookProfile(hdNdEtadPhi_caloJet_Toward500_,"hdNdEtadPhi_caloJet_Toward500", 
00438                  "Average number of tracks (pT > 500 MeV) in toward region vs leading calo jet pT;pT(GeV/c);dN/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00439     bookProfile(hdNdEtadPhi_caloJet_Transverse500_,"hdNdEtadPhi_caloJet_Transverse500", 
00440                  "Average number of tracks (pT > 500 MeV) in transverse region vs leading calo jet pT;pT(GeV/c);dN/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00441     bookProfile(hdNdEtadPhi_caloJet_Away500_,"hdNdEtadPhi_caloJet_Away500", 
00442                  "Average number of tracks (pT > 500 MeV) in away region vs leading calo jet pT;pT(GeV/c);dN/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);   
00443   */
00444 
00445     bookProfile(hdNdEtadPhi_trackJet_Toward500_,"hdNdEtadPhi_trackJet_Toward500", 
00446                  "Average number of tracks (pT > 500 MeV) in toward region vs leading track jet pT;pT(GeV/c);dN/d#eta d#phi",Nx,x1,x2,y1,y2);
00447     bookProfile(hdNdEtadPhi_trackJet_Transverse500_,"hdNdEtadPhi_trackJet_Transverse500", 
00448                  "Average number of tracks (pT > 500 MeV) in transverse region vs leading track jet pT;pT(GeV/c);dN/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00449     bookProfile(hdNdEtadPhi_trackJet_Away500_,"hdNdEtadPhi_trackJet_Away500", 
00450                  "Average number of tracks (pT > 500 MeV) in away region vs leading track jet pT;pT(GeV/c);dN/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00451 
00452 
00453  
00454     bookProfile(hpTSumdEtadPhi_pTMax_Toward500_,"hpTSumdEtadPhi_pTMax_Toward500", 
00455                  "Average number of tracks (pT > 500 MeV) in toward region vs leading track pT;pT(GeV/c);dpTSum/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00456     bookProfile(hpTSumdEtadPhi_pTMax_Transverse500_,"hpTSumdEtadPhi_pTMax_Transverse500", 
00457                  "Average number of tracks (pT > 500 MeV) in transverse region vs leading track pT;pT(GeV/c);dpTSum/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00458     bookProfile(hpTSumdEtadPhi_pTMax_Away500_,"hpTSumdEtadPhi_pTMax_Away500", 
00459                  "Average number of tracks (pT > 500 MeV) in away region vs leading track pT;pT(GeV/c);dpTSum/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00460 /* 
00461     bookProfile(hpTSumdEtadPhi_caloJet_Toward500_,"hpTSumdEtadPhi_caloJet_Toward500", 
00462                  "Average number of tracks (pT > 500 MeV) in toward region vs leading calo jet pT;pT(GeV/c);dpTSum/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00463     bookProfile(hpTSumdEtadPhi_caloJet_Transverse500_,"hpTSumdEtadPhi_caloJet_Transverse500", 
00464                  "Average number of tracks (pT > 500 MeV) in transverse region vs leading calo jet pT;pT(GeV/c);dpTSum/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00465     bookProfile(hpTSumdEtadPhi_caloJet_Away500_,"hpTSumdEtadPhi_caloJet_Away500", 
00466                  "Average number of tracks (pT > 500 MeV) in away region vs leading calo jet pT;pT(GeV/c);dpTSum/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);   
00467 */
00468   
00469     bookProfile(hpTSumdEtadPhi_trackJet_Toward500_,"hpTSumdEtadPhi_trackJet_Toward500", 
00470                  "Average number of tracks (pT > 500 MeV) in toward region vs leading track jet pT;pT(GeV/c);dpTSum/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00471     bookProfile(hpTSumdEtadPhi_trackJet_Transverse500_,"hpTSumdEtadPhi_trackJet_Transverse500", 
00472                  "Average number of tracks (pT > 500 MeV) in transverse region vs leading track jet pT;pT(GeV/c);dpTSum/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00473     bookProfile(hpTSumdEtadPhi_trackJet_Away500_,"hpTSumdEtadPhi_trackJet_Away500", 
00474                  "Average number of tracks (pT > 500 MeV) in away region vs leading track jet pT;pT(GeV/c);dpTSum/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00475 
00476    
00477    bookProfile(hdNdEtadPhi_pTMax_Toward900_,"hdNdEtadPhi_pTMax_Toward900",
00478                  "Average number of tracks (pT > 900 MeV) in toward region vs leading track pT;pT(GeV/c);dN/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00479     bookProfile(hdNdEtadPhi_pTMax_Transverse900_,"hdNdEtadPhi_pTMax_Transverse900",
00480                  "Average number of tracks (pT > 900 MeV) in transverse region vs leading track pT;pT(GeV/c);dN/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00481     bookProfile(hdNdEtadPhi_pTMax_Away900_,"hdNdEtadPhi_pTMax_Away900",
00482                  "Average number of tracks (pT > 900 MeV) in away region vs leading track pT;pT(GeV/c);dN/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00483 /*
00484     bookProfile(hdNdEtadPhi_caloJet_Toward900_,"hdNdEtadPhi_caloJet_Toward900",
00485                  "Average number of tracks (pT > 900 MeV) in toward region vs leading calo jet pT;pT(GeV/c);dN/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00486     bookProfile(hdNdEtadPhi_caloJet_Transverse900_,"hdNdEtadPhi_caloJet_Transverse900",
00487                  "Average number of tracks (pT > 900 MeV) in transverse region vs leading calo jet pT;pT(GeV/c);dN/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00488     bookProfile(hdNdEtadPhi_caloJet_Away900_,"hdNdEtadPhi_caloJet_Away900",
00489                  "Average number of tracks (pT > 900 MeV) in away region vs leading calo jet pT;pT(GeV/c);dN/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00490 */
00491 
00492     bookProfile(hdNdEtadPhi_trackJet_Toward900_,"hdNdEtadPhi_trackJet_Toward900",
00493                  "Average number of tracks (pT > 900 MeV) in toward region vs leading track jet pT;pT(GeV/c);dN/d#eta d#phi",Nx,x1,x2,y1,y2);
00494     bookProfile(hdNdEtadPhi_trackJet_Transverse900_,"hdNdEtadPhi_trackJet_Transverse900",
00495                  "Average number of tracks (pT > 900 MeV) in transverse region vs leading track jet pT;pT(GeV/c);dN/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00496     bookProfile(hdNdEtadPhi_trackJet_Away900_,"hdNdEtadPhi_trackJet_Away900",
00497                  "Average number of tracks (pT > 900 MeV) in away region vs leading track jet pT;pT(GeV/c);dN/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00498 
00499 
00500 
00501     bookProfile(hpTSumdEtadPhi_pTMax_Toward900_,"hpTSumdEtadPhi_pTMax_Toward900",
00502                  "Average number of tracks (pT > 900 MeV) in toward region vs leading track pT;pT(GeV/c);dpTSum/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00503     bookProfile(hpTSumdEtadPhi_pTMax_Transverse900_,"hpTSumdEtadPhi_pTMax_Transverse900",
00504                  "Average number of tracks (pT > 900 MeV) in transverse region vs leading track pT;pT(GeV/c);dpTSum/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00505     bookProfile(hpTSumdEtadPhi_pTMax_Away900_,"hpTSumdEtadPhi_pTMax_Away900",
00506                  "Average number of tracks (pT > 900 MeV) in away region vs leading track pT;pT(GeV/c);dpTSum/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00507 /*
00508     bookProfile(hpTSumdEtadPhi_caloJet_Toward900_,"hpTSumdEtadPhi_caloJet_Toward900",
00509                  "Average number of tracks (pT > 900 MeV) in toward region vs leading calo jet pT;pT(GeV/c);dpTSum/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00510     bookProfile(hpTSumdEtadPhi_caloJet_Transverse900_,"hpTSumdEtadPhi_caloJet_Transverse900",
00511                  "Average number of tracks (pT > 900 MeV) in transverse region vs leading calo jet pT;pT(GeV/c);dpTSum/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00512     bookProfile(hpTSumdEtadPhi_caloJet_Away900_,"hpTSumdEtadPhi_caloJet_Away900",
00513                  "Average number of tracks (pT > 900 MeV) in away region vs leading calo jet pT;pT(GeV/c);dpTSum/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00514 */
00515     bookProfile(hpTSumdEtadPhi_trackJet_Toward900_,"hpTSumdEtadPhi_trackJet_Toward900",
00516                  "Average number of tracks (pT > 900 MeV) in toward region vs leading track jet pT;pT(GeV/c);dpTSum/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00517     bookProfile(hpTSumdEtadPhi_trackJet_Transverse900_,"hpTSumdEtadPhi_trackJet_Transverse900",
00518                  "Average number of tracks (pT > 900 MeV) in transverse region vs leading track jet pT;pT(GeV/c);dpTSum/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00519     bookProfile(hpTSumdEtadPhi_trackJet_Away900_,"hpTSumdEtadPhi_trackJet_Away900",
00520                  "Average number of tracks (pT > 900 MeV) in away region vs leading track jet pT;pT(GeV/c);dpTSum/d#eta d#phi",Nx,x1,x2,y1,y2,0,0);
00521  
00522    
00523   }
00524 
00525 if (1) {
00526     const int Nx = 20;
00527     const double x1 = 0.0;
00528     const double x2 = 20.0;
00529 
00530         book1D(hChargedJetMulti_,"hChargedJetMulti","Charged jet multiplicity;multiplicities",Nx,x1,x2);
00531         //book1D(hCaloJetMulti_,"hCaloJetMulti","Calo jet multiplicity;multiplicities",Nx,x1,x2);
00532 
00533   }
00534 
00535 
00536 if (1) {
00537     const int Nx = 60;
00538     const double x1 = -180.0;
00539     const double x2 = 180.0;
00540 
00541         book1D(hdPhi_maxpTTrack_tracks_,"hdPhi_maxpTTrack_tracks","delta phi between leading tracks and other tracks;#Delta#phi(leading track-track)",Nx,x1,x2);
00542 //        book1D(hdPhi_caloJet_tracks_,"hdPhi_caloJet_tracks","delta phi between leading calo jet  and tracks;#Delta#phi(leading calo jet-track)",Nx,x1,x2);
00543         book1D(hdPhi_chargedJet_tracks_,"hdPhi_chargedJet_tracks","delta phi between leading charged jet  and tracks;#Delta#phi(leading charged jet-track)",Nx,x1,x2);
00544 
00545 }
00546             
00547 
00548 }
00549 
00550 //--------------------------------------------------------------------------------------------------
00551 void QcdUeDQM::endJob(void) 
00552 {
00553 }
00554 
00555 //--------------------------------------------------------------------------------------------------
00556 
00557 void QcdUeDQM::endLuminosityBlock(const LuminosityBlock &l, 
00558                                      const EventSetup &iSetup)
00559 {
00560   if( ! isHltConfigSuccessful_ ) return;
00561 
00562   // Update various histograms.
00563 
00564   //repSummary_->Fill(1.);
00565  // repSumMap_->Fill(0.5,0.5,1.);
00566 
00567 }
00568 
00569 //--------------------------------------------------------------------------------------------------
00570 
00571 void QcdUeDQM::endRun(const Run &, const EventSetup &)
00572 {
00573   if( ! isHltConfigSuccessful_ ) return;
00574 
00575   // End run, cleanup. TODO: can this be called several times in DQM???
00576 
00577 }
00578 
00579 //--------------------------------------------------------------------------------------------------
00580 void QcdUeDQM::fill1D(std::vector<TH1F*> &hs, double val, double w)
00581 {
00582   // Loop over histograms and fill if trigger has fired.
00583 
00584   for(size_t i=0;i<hs.size();++i) {
00585     if (!hltTrgDeci_.at(i))
00586       continue;
00587     hs.at(i)->Fill(val,w);
00588   }
00589 }
00590 
00591 //--------------------------------------------------------------------------------------------------
00592 void QcdUeDQM::fill1D(std::vector<MonitorElement*> &mes, double val, double w)
00593 {
00594   // Loop over histograms and fill if trigger has fired.
00595 
00596   for(size_t i=0;i<mes.size();++i) {
00597     if (!hltTrgDeci_.at(i))
00598       continue;
00599     mes.at(i)->Fill(val,w);
00600   }
00601 }
00602 
00603 //--------------------------------------------------------------------------------------------------
00604 void QcdUeDQM::setLabel1D(std::vector<MonitorElement*> &mes)
00605 {
00606   // Loop over histograms and fill if trigger has fired.
00607   string cut[5] = {"Nevt","vtx!=bmspt","Zvtx<10cm","pT>1GeV","trackFromVtx"};
00608   for(size_t i=0;i<mes.size();++i) {
00609     if (!hltTrgDeci_.at(i))
00610       continue;
00611     for(size_t j = 1;j < 6;j++)mes.at(i)->setBinLabel(j,cut[j-1],1);
00612   }
00613 }
00614 
00615 //--------------------------------------------------------------------------------------------------
00616 void QcdUeDQM::fill2D(std::vector<TH2F*> &hs, double valx, double valy, double w)
00617 {
00618   // Loop over histograms and fill if trigger has fired.
00619 
00620   for(size_t i=0;i<hs.size();++i) {
00621     if (!hltTrgDeci_.at(i))
00622       continue;
00623     hs.at(i)->Fill(valx, valy ,w);
00624   }
00625 }
00626 
00627 //--------------------------------------------------------------------------------------------------
00628 void QcdUeDQM::fill2D(std::vector<MonitorElement*> &mes, double valx, double valy, double w)
00629 {
00630   // Loop over histograms and fill if trigger has fired.
00631 
00632   for(size_t i=0;i<mes.size();++i) {
00633     if (!hltTrgDeci_.at(i))
00634       continue;
00635     mes.at(i)->Fill(valx, valy ,w);
00636   }
00637 }
00638 //--------------------------------------------------------------------------------------------------
00639 void QcdUeDQM::fillProfile(std::vector<TProfile*> &hs, double valx, double valy, double w)
00640 {
00641   // Loop over histograms and fill if trigger has fired.
00642 
00643   for(size_t i=0;i<hs.size();++i) {
00644     if (!hltTrgDeci_.at(i))
00645       continue;
00646     hs.at(i)->Fill(valx, valy ,w);
00647   }
00648 }
00649 
00650 //--------------------------------------------------------------------------------------------------
00651 void QcdUeDQM::fillProfile(std::vector<MonitorElement*> &mes, double valx, double valy, double w)
00652 {
00653   // Loop over histograms and fill if trigger has fired.
00654 
00655   for(size_t i=0;i<mes.size();++i) {
00656     if (!hltTrgDeci_.at(i))
00657       continue;
00658    const double y = valy*w; 
00659     mes.at(i)->Fill(valx, y);
00660   }
00661 }
00662 
00663 //--------------------------------------------------------------------------------------------------
00664 bool QcdUeDQM::trackSelection(const reco::Track &trk, const reco::BeamSpot* bs, const reco::Vertex vtx, int sizevtx )
00665 {
00666 //-------------Fill basic histograms---------
00667 
00668 
00669 //-------------------------------------------
00670 
00671  bool goodTrk = false;
00672 
00673  if(sizevtx!=1) return 0;    //selection events with only a vertex
00674 
00675  //Fill basic information of all the tracks
00676   fill1D(hNtrackerLayer_,trk.hitPattern().trackerLayersWithMeasurement());
00677   fill1D(hNtrackerPixelLayer_,trk.hitPattern().pixelLayersWithMeasurement());
00678   fill1D(hNtrackerStripPixelLayer_,(trk.hitPattern().pixelLayersWithMeasurement() +  trk.hitPattern().numberOfValidStripLayersWithMonoAndStereo()));
00679   fill1D(hRatioPtErrorPt_,(trk.ptError()/trk.pt()));
00680   fill1D(hTrkPt_,trk.pt());
00681   fill1D(hTrkEta_,trk.eta());  
00682   fill1D(hTrkPhi_,trk.phi());
00683   fill1D(hRatioDxySigmaDxyBS_,(trk.dxy(bs->position())/trk.dxyError()));
00684   fill1D(hRatioDxySigmaDxyPV_,(trk.dxy(vtx.position())/trk.dxyError()));
00685   fill1D(hRatioDzSigmaDzBS_,(trk.dz(bs->position())/trk.dzError()));
00686   fill1D(hRatioDzSigmaDzPV_,(trk.dz(vtx.position())/trk.dzError()));
00687   fill1D(hTrkChi2_,trk.normalizedChi2());
00688   fill1D(hTrkNdof_,trk.ndof());   
00689 
00690   fill1D(hBeamSpot_x_,bs->x0());
00691   fill1D(hBeamSpot_y_,bs->y0()); 
00692   fill1D(hBeamSpot_z_,bs->z0());
00693  
00694  //number of layers
00695  bool layerMinCutbool=false;
00696     if (trk.hitPattern().trackerLayersWithMeasurement() >= minHit_ ||
00697                 (trk.hitPattern().trackerLayersWithMeasurement()==3 && trk.hitPattern().pixelLayersWithMeasurement()==3 && allowTriplets_))
00698       layerMinCutbool=true;
00699 
00700 
00701   //number of pixel layers
00702  bool pxlLayerMinCutbool=false;
00703     if (trk.hitPattern().pixelLayersWithMeasurement() >=pxlLayerMinCut_) pxlLayerMinCutbool=true;
00704 
00705 
00706  // cut on the hits in pixel layers
00707  bool hasPIX1 = false;
00708     if (requirePIX1_) {
00709       const reco::HitPattern& p = trk.hitPattern();
00710       for (int i=0; i<p.numberOfHits(); i++) {
00711         uint32_t hit = p.getHitPattern(i);
00712         if (p.validHitFilter(hit) && p.pixelHitFilter(hit) && p.getLayer(hit)==1) hasPIX1 = true;
00713       }
00714     }else hasPIX1 = true;
00715  
00716  // cut on the pT error
00717  bool ptErrorbool=false; 
00718     if (trk.ptError()/trk.pt() < ptErr_pt_ || 
00719         (trk.hitPattern().trackerLayersWithMeasurement()==3 && trk.hitPattern().pixelLayersWithMeasurement()==3 && allowTriplets_)) ptErrorbool=true; 
00720  // quality cut
00721  bool quality_ok = true;  
00722  if (quality_.size()!=0) {
00723       quality_ok = false;
00724       for (unsigned int i = 0; i<quality_.size();++i) {
00725         if (trk.quality(quality_[i])){
00726           quality_ok = true;
00727           break;
00728         }
00729       }
00730     }
00731  //-----
00732  bool algo_ok = true;
00733     if (algorithm_.size()!=0) {
00734       if (std::find(algorithm_.begin(),algorithm_.end(),trk.algo())==algorithm_.end()) algo_ok = false;
00735     }
00736  
00737 
00738  if(bsuse_==1)
00739       {
00740     if(hasPIX1 &&  pxlLayerMinCutbool && layerMinCutbool &&  (trk.hitPattern().pixelLayersWithMeasurement() +  trk.hitPattern().numberOfValidStripLayersWithMonoAndStereo()) >= min3DHit_ && ptErrorbool && fabs(trk.pt()) >= ptMin_ &&  trk.eta() >= minRapidity_ && trk.eta() <= maxRapidity_ &&  fabs(trk.dxy(bs->position())/trk.dxyError()) < tip_ &&  fabs(trk.dz(bs->position())/trk.dzError()) < lip_  &&  trk.normalizedChi2()<=maxChi2_ &&  quality_ok &&  algo_ok)goodTrk=true ;
00741       }
00742 
00743     if(bsuse_==0)
00744      {
00745        if(hasPIX1 &&  pxlLayerMinCutbool &&  layerMinCutbool &&  (trk.hitPattern().pixelLayersWithMeasurement() +  trk.hitPattern().numberOfValidStripLayersWithMonoAndStereo())  >= min3DHit_ && ptErrorbool && fabs(trk.pt()) >= ptMin_ && trk.eta() >= minRapidity_ && trk.eta() <= maxRapidity_ && fabs(trk.dxy(vtx.position())/trk.dxyError()) < tip_ && fabs(trk.dz(vtx.position())/trk.dzError()) < lip_  && trk.normalizedChi2()<=maxChi2_ && quality_ok &&  algo_ok)goodTrk=true;
00746      
00747      }
00748 
00749   return goodTrk;
00750 
00751  }
00752 //--------------------------------------------------------------------------------------------------
00753 bool  QcdUeDQM::fillVtxPlots( const reco::BeamSpot* bs, const edm::Handle< reco::VertexCollection > vtxColl)
00754 {
00755   const reco::VertexCollection theVertices = *(vtxColl.product());
00756   bool goodVtx = false;
00757   fill1D(hNvertices_,theVertices.size()); 
00758     for (reco::VertexCollection::const_iterator vertexIt = theVertices.begin(); vertexIt != theVertices.end(); ++vertexIt) 
00759       {
00760         fill1D(hVertex_z_,vertexIt->z());
00761         fill1D(hVertex_y_,vertexIt->y());
00762         fill1D(hVertex_x_,vertexIt->x());
00763         fill1D(hVertex_ndof_,vertexIt->ndof());
00764         fill1D(hVertex_rho_,vertexIt->position().rho()); 
00765         fill1D(hVertex_z_bs_,(vertexIt->z() -bs->z0()));
00766 
00767         if(fabs(vertexIt->z() -bs->z0()) < diffvtxbs_  && vertexIt->ndof() >= 4 && vertexIt->position().rho()<= 2.0)
00768          {
00769          goodVtx = true;
00770          vtx1=(*vertexIt);
00771          
00772          break;
00773          }
00774       } // Loop over vertcies
00775    return goodVtx;
00776 }
00777 //--------------------------------------------------------------------------------------------------
00778 void QcdUeDQM::fillpTMaxRelated(const std::vector<const reco::Track *> &track)
00779  {
00780    fill1D(hNgoodTrk_,track.size());
00781    if(track.size()>0)
00782    {
00783    fill1D(hLeadingTrack_pTSpectrum_,track[0]->pt());
00784    fill1D(hLeadingTrack_phiSpectrum_,track[0]->phi());
00785    fill1D(hLeadingTrack_etaSpectrum_,track[0]->eta());
00786    }
00787      for(size_t i = 0; i < track.size(); i++)
00788        {
00789         fill1D(hGoodTrkPt500_,track[i]->pt());
00790         fill1D(hGoodTrkEta500_,track[i]->eta());
00791         fill1D(hGoodTrkPhi500_,track[i]->phi());
00792         if(track[i]->pt() > 0.9)
00793         {
00794         fill1D(hGoodTrkPt900_,track[i]->pt());
00795         fill1D(hGoodTrkEta900_,track[i]->eta());
00796         fill1D(hGoodTrkPhi900_,track[i]->phi());
00797         }
00798         }
00799                      
00800  }
00801 
00802 
00803 void QcdUeDQM::fillChargedJetSpectra(const edm::Handle<reco::TrackJetCollection> trackJets)
00804 {
00805   fill1D(hChargedJetMulti_,trackJets->size());
00806   for( reco::TrackJetCollection::const_iterator f  = trackJets->begin();  f != trackJets->end(); f++) 
00807     {
00808       if(f != trackJets->begin())continue;
00809       fill1D(hLeadingChargedJet_pTSpectrum_,f->pt());
00810       fill1D(hLeadingChargedJet_etaSpectrum_,f->eta());
00811       fill1D(hLeadingChargedJet_phiSpectrum_,f->phi());
00812     } 
00813         
00814 }
00815 
00816 /*
00817 void QcdUeDQM::fillCaloJetSpectra(const edm::Handle<reco::CaloJetCollection> caloJets)
00818 {
00819   fill1D(hCaloJetMulti_,caloJets->size());
00820    for( reco::CaloJetCollection::const_iterator f  = caloJets->begin();  f != caloJets->end(); f++)
00821      {
00822        if(f != caloJets->begin())continue;
00823        fill1D(hLeadingCaloJet_pTSpectrum_,f->pt()); 
00824        fill1D(hLeadingCaloJet_etaSpectrum_,f->eta());
00825        fill1D(hLeadingCaloJet_phiSpectrum_,f->phi());
00826      }
00827    
00828 }
00829 
00830 
00831 */
00832 
00833 
00834 /*
00835  weight for transverse/toward/away region = 0.12
00836  
00837 
00838 */
00839 
00840 void QcdUeDQM::fillUE_with_MaxpTtrack(const std::vector<const reco::Track*>  &track)
00841 {
00842 double w = 0.119;          
00843 //double w = 1.;
00844 double nTrk500_TransReg = 0;
00845 double nTrk500_AwayReg = 0;
00846 double nTrk500_TowardReg = 0;
00847  
00848 double pTSum500_TransReg = 0;
00849 double pTSum500_AwayReg = 0;
00850 double pTSum500_TowardReg = 0;
00851 
00852 
00853 double nTrk900_TransReg = 0;
00854 double nTrk900_AwayReg = 0;
00855 double nTrk900_TowardReg = 0;
00856 
00857 double pTSum900_TransReg = 0;
00858 double pTSum900_AwayReg = 0;
00859 double pTSum900_TowardReg = 0;
00860    if(track.size() > 0) 
00861     {
00862      if(track[0]->pt() > 1.)
00863          {
00864             for(size_t i = 1; i < track.size();i++)
00865                  {
00866                         
00867                        double dphi = (180./PI)*(deltaPhi(track[0]->phi(),track[i]->phi()));
00868                        fill1D(hdPhi_maxpTTrack_tracks_,dphi);
00869                        if(fabs(dphi)>60. && fabs(dphi)<120.)
00870                          {
00871                                pTSum500_TransReg =  pTSum500_TransReg + track[i]->pt();     
00872                                nTrk500_TransReg++;
00873                                if(track[i]->pt() > 0.9)
00874                                {
00875                                pTSum900_TransReg =  pTSum900_TransReg + track[i]->pt();
00876                                nTrk900_TransReg++;
00877                                }
00878                          }            
00879                         
00880                        if(fabs(dphi)>120. && fabs(dphi)<180.)
00881                          {
00882                                pTSum500_AwayReg =  pTSum500_AwayReg + track[i]->pt();   
00883                                nTrk500_AwayReg++;
00884                                 if(track[i]->pt() > 0.9)
00885                                 {
00886                                 pTSum900_AwayReg =  pTSum900_AwayReg + track[i]->pt();
00887                                 nTrk900_AwayReg++;
00888 
00889                                 } 
00890                          } 
00891                        
00892                        if(fabs(dphi)<60.)
00893                          {
00894                                pTSum500_TowardReg =  pTSum500_TowardReg + track[i]->pt();
00895                                nTrk500_TowardReg++;
00896                                if(track[i]->pt() > 0.9)
00897                                {
00898                                pTSum900_TowardReg =  pTSum900_TowardReg + track[i]->pt();
00899                                nTrk900_TowardReg++;
00900                                } 
00901                          }           
00902                      } // track loop 
00903                  }// leading track
00904              // non empty collection
00905                fillProfile(hdNdEtadPhi_pTMax_Toward500_, track[0]->pt(),nTrk500_TowardReg,w);
00906                fillProfile(hdNdEtadPhi_pTMax_Transverse500_, track[0]->pt(),nTrk500_TransReg,w);
00907                fillProfile(hdNdEtadPhi_pTMax_Away500_, track[0]->pt(),nTrk500_AwayReg,w);
00908 
00909                fillProfile(hpTSumdEtadPhi_pTMax_Toward500_,track[0]->pt() ,pTSum500_TowardReg,w);
00910                fillProfile(hpTSumdEtadPhi_pTMax_Transverse500_,track[0]->pt(),pTSum500_TransReg,w);
00911                fillProfile(hpTSumdEtadPhi_pTMax_Away500_, track[0]->pt(),pTSum500_AwayReg,w);
00912 
00913                fillProfile(hdNdEtadPhi_pTMax_Toward900_, track[0]->pt(),nTrk900_TowardReg,w);
00914                fillProfile(hdNdEtadPhi_pTMax_Transverse900_, track[0]->pt(),nTrk900_TransReg,w);
00915                fillProfile(hdNdEtadPhi_pTMax_Away900_, track[0]->pt(),nTrk900_AwayReg,w);
00916 
00917                fillProfile(hpTSumdEtadPhi_pTMax_Toward900_,track[0]->pt() ,pTSum900_TowardReg,w);
00918                fillProfile(hpTSumdEtadPhi_pTMax_Transverse900_,track[0]->pt(),pTSum900_TransReg,w);
00919                fillProfile(hpTSumdEtadPhi_pTMax_Away900_, track[0]->pt(),pTSum900_AwayReg,w);
00920      }
00921 }
00922 
00923 void QcdUeDQM::fillUE_with_ChargedJets(const std::vector<const reco::Track *>  &track, const edm::Handle<reco::TrackJetCollection> &trackJets)
00924 {
00925 double w = 0.119;
00926 double nTrk500_TransReg = 0;
00927 double nTrk500_AwayReg = 0;
00928 double nTrk500_TowardReg = 0;
00929   
00930 double pTSum500_TransReg = 0;
00931 double pTSum500_AwayReg = 0;
00932 double pTSum500_TowardReg = 0;
00933 
00934 
00935 double nTrk900_TransReg = 0;
00936 double nTrk900_AwayReg = 0;
00937 double nTrk900_TowardReg = 0;
00938 
00939 double pTSum900_TransReg = 0;
00940 double pTSum900_AwayReg = 0;
00941 double pTSum900_TowardReg = 0;
00942 
00943    if(!(trackJets->empty()) && (trackJets->begin())->pt() > 1.)
00944          {
00945          double jetPhi = (trackJets->begin())->phi();
00946            for(size_t i = 0; i < track.size();i++)
00947                    {
00948                          double dphi = (180./PI)*(deltaPhi(jetPhi,track[i]->phi()));
00949                          fill1D(hdPhi_chargedJet_tracks_,dphi);
00950                          if(fabs(dphi)>60. && fabs(dphi)<120.)
00951                            {
00952                                  pTSum500_TransReg =  pTSum500_TransReg + track[i]->pt();
00953                                  nTrk500_TransReg++;
00954                                  if(track[i]->pt() > 0.9)
00955                                  {
00956                                  pTSum900_TransReg =  pTSum900_TransReg + track[i]->pt();
00957                                  nTrk900_TransReg++;
00958                                  }
00959                            }
00960                         
00961                          if(fabs(dphi)>120. && fabs(dphi)<180.)
00962                            {
00963                                  pTSum500_AwayReg =  pTSum500_AwayReg + track[i]->pt();
00964                                  nTrk500_AwayReg++;
00965                                   if(track[i]->pt() > 0.9)
00966                                  {
00967                                  pTSum900_AwayReg =  pTSum900_AwayReg + track[i]->pt();
00968                                  nTrk900_AwayReg++;
00969                                  }
00970                            }
00971                          if(fabs(dphi)<60.)
00972                            {
00973                                  pTSum500_TowardReg =  pTSum500_TowardReg + track[i]->pt();
00974                                  nTrk500_TowardReg++;
00975                                   if(track[i]->pt() > 0.9)
00976                                  {
00977                                  pTSum900_TowardReg =  pTSum900_TowardReg + track[i]->pt();
00978                                  nTrk900_TowardReg++;
00979                                  } 
00980                            }
00981                        }// tracks loop 
00982 
00983                    }// leading track jet
00984                  
00985                  fillProfile(hdNdEtadPhi_trackJet_Toward500_, (trackJets->begin())->pt(),nTrk500_TowardReg,w);
00986                  fillProfile(hdNdEtadPhi_trackJet_Transverse500_, (trackJets->begin())->pt(),nTrk500_TransReg,w);
00987                  fillProfile(hdNdEtadPhi_trackJet_Away500_, (trackJets->begin())->pt(),nTrk500_AwayReg,w);
00988                  
00989                  fillProfile(hpTSumdEtadPhi_trackJet_Toward500_, (trackJets->begin())->pt(),pTSum500_TowardReg,w);
00990                  fillProfile(hpTSumdEtadPhi_trackJet_Transverse500_, (trackJets->begin())->pt(),pTSum500_TransReg,w);
00991                  fillProfile(hpTSumdEtadPhi_trackJet_Away500_, (trackJets->begin())->pt(),pTSum500_AwayReg,w);
00992 
00993                  fillProfile(hdNdEtadPhi_trackJet_Toward900_, (trackJets->begin())->pt(),nTrk900_TowardReg,w);
00994                  fillProfile(hdNdEtadPhi_trackJet_Transverse900_, (trackJets->begin())->pt(),nTrk900_TransReg,w);
00995                  fillProfile(hdNdEtadPhi_trackJet_Away900_, (trackJets->begin())->pt(),nTrk900_AwayReg,w);
00996 
00997                  fillProfile(hpTSumdEtadPhi_trackJet_Toward900_, (trackJets->begin())->pt(),pTSum900_TowardReg,w);
00998                  fillProfile(hpTSumdEtadPhi_trackJet_Transverse900_, (trackJets->begin())->pt(),pTSum900_TransReg,w);
00999                  fillProfile(hpTSumdEtadPhi_trackJet_Away900_, (trackJets->begin())->pt(),pTSum900_AwayReg,w);  
01000 }
01001 
01002 /*
01003 void QcdUeDQM:: fillUE_with_CaloJets(const std::vector<const reco::Track *>  &track, const edm::Handle<reco::CaloJetCollection> &caloJets)
01004 {
01005 double w = 0.119;
01006 double nTrk500_TransReg = 0;
01007 double nTrk500_AwayReg = 0;
01008 double nTrk500_TowardReg = 0;
01009                  
01010 double pTSum500_TransReg = 0;
01011 double pTSum500_AwayReg = 0;
01012 double pTSum500_TowardReg = 0;
01013 
01014 double nTrk900_TransReg = 0;
01015 double nTrk900_AwayReg = 0;
01016 double nTrk900_TowardReg = 0;
01017 
01018 double pTSum900_TransReg = 0;
01019 double pTSum900_AwayReg = 0;
01020 double pTSum900_TowardReg = 0;
01021     if(!(caloJets->empty()) && (caloJets->begin())->pt() > 1.)
01022           {
01023             double jetPhi = (caloJets->begin())->phi();
01024              for(size_t i = 0; i < track.size();i++)
01025                    {
01026                          double dphi = (180./PI)*(deltaPhi(jetPhi,track[i]->phi()));
01027                          fill1D(hdPhi_caloJet_tracks_,dphi);
01028                          if(fabs(dphi)>60. && fabs(dphi)<120.)
01029                            {
01030                                  pTSum500_TransReg =  pTSum500_TransReg + track[i]->pt();
01031                                  nTrk500_TransReg++;
01032                                  if(track[i]->pt() > 0.9)
01033                                  {
01034                                  pTSum900_TransReg =  pTSum900_TransReg + track[i]->pt();
01035                                  nTrk900_TransReg++;
01036                                  }
01037                            }
01038                          if(fabs(dphi)>120. && fabs(dphi)<180.)
01039                            {
01040                                  pTSum500_AwayReg =  pTSum500_AwayReg + track[i]->pt();
01041                                  nTrk500_AwayReg++;
01042                                  if(track[i]->pt() > 0.9)
01043                                  {
01044                                  pTSum900_AwayReg =  pTSum900_AwayReg + track[i]->pt();
01045                                  nTrk900_AwayReg++;
01046                                  }
01047                            }
01048                          if(fabs(dphi)<60.)
01049                            {
01050                                  pTSum500_TowardReg =  pTSum500_TowardReg + track[i]->pt();
01051                                  nTrk500_TowardReg++;
01052                                  if(track[i]->pt() > 0.9)
01053                                  {
01054                                  pTSum900_TowardReg =  pTSum900_TowardReg + track[i]->pt();
01055                                  nTrk900_TowardReg++; 
01056                                  }
01057                            }
01058                        }// tracks loop 
01059 
01060                    }// leading calo jet
01061                  fillProfile(hdNdEtadPhi_caloJet_Toward500_, (caloJets->begin())->pt(),nTrk500_TowardReg,w);
01062                  fillProfile(hdNdEtadPhi_caloJet_Transverse500_, (caloJets->begin())->pt(),nTrk500_TransReg,w);
01063                  fillProfile(hdNdEtadPhi_caloJet_Away500_, (caloJets->begin())->pt(),nTrk500_AwayReg,w);
01064                    
01065                  fillProfile(hpTSumdEtadPhi_caloJet_Toward500_, (caloJets->begin())->pt(),pTSum500_TowardReg,w);
01066                  fillProfile(hpTSumdEtadPhi_caloJet_Transverse500_, (caloJets->begin())->pt(),pTSum500_TransReg,w);
01067                  fillProfile(hpTSumdEtadPhi_caloJet_Away500_, (caloJets->begin())->pt(),pTSum500_AwayReg,w);
01068 
01069                  fillProfile(hdNdEtadPhi_caloJet_Toward900_, (caloJets->begin())->pt(),nTrk900_TowardReg,w);
01070                  fillProfile(hdNdEtadPhi_caloJet_Transverse900_, (caloJets->begin())->pt(),nTrk900_TransReg,w);
01071                  fillProfile(hdNdEtadPhi_caloJet_Away900_, (caloJets->begin())->pt(),nTrk900_AwayReg,w);
01072 
01073                  fillProfile(hpTSumdEtadPhi_caloJet_Toward900_, (caloJets->begin())->pt(),pTSum900_TowardReg,w);
01074                  fillProfile(hpTSumdEtadPhi_caloJet_Transverse900_, (caloJets->begin())->pt(),pTSum900_TransReg,w);
01075                  fillProfile(hpTSumdEtadPhi_caloJet_Away900_, (caloJets->begin())->pt(),pTSum900_AwayReg,w);
01076 
01077                  
01078 }
01079 
01080 */
01081 void QcdUeDQM::fillHltBits(const Event &iEvent,const EventSetup &iSetup)
01082 {
01083   // Fill HLT trigger bits.
01084 
01085   Handle<TriggerResults> triggerResultsHLT;
01086   getProduct(hltUsedResName_, triggerResultsHLT, iEvent);
01087   
01088 /*  const unsigned int ntrigs(triggerResultsHLT.product()->size());
01089   if( ntrigs != 0 ) {
01090   const edm::TriggerNames & triggerNames = iEvent.triggerNames(*triggerResultsHLT);
01091   for (unsigned int j=0; j!=ntrigs; j++) {
01092 //    if(triggerResultsHLT->accept(j)){
01093 //    unsigned int hlt_prescale = hltConfig_.prescaleValue(iEvent, iSetup, triggerName);
01094     cout<<"trigger fired"<<triggerNames.triggerName(j)<<endl;
01095     for(unsigned int itrigName = 0; itrigName < hltTrgNames_.size(); itrigName++ ) {
01096      unsigned int hlt_prescale = hltConfig.prescaleValue(iEvent, iSetup, hltTrgNames_[itrigName]);
01097 
01098      if(triggerNames.triggerIndex(hltTrgNames_[itrigName]) >= (unsigned int)ntrigs ) continue; 
01099 //    if( triggerResultsHLT->accept(triggerNames.triggerIndex(hltTrgNames_[itrigName])) )cout<<hltTrgNames_[itrigName]<<endl;
01100 
01101  }
01102     }
01103    } 
01104  // }
01105 */
01106   for(size_t i=0;i<hltTrgBits_.size();++i) {
01107     if (hltTrgBits_.at(i)<0) 
01108       continue; //ignore unknown trigger 
01109     size_t tbit = hltTrgBits_.at(i);
01110     if (tbit<triggerResultsHLT->size()) {
01111       hltTrgDeci_[i] = triggerResultsHLT->accept(tbit);
01112   }
01113  }
01114   // fill correlation histogram
01115   for(size_t i=0;i<hltTrgBits_.size();++i) {
01116     if (hltTrgDeci_.at(i))
01117       h2TrigCorr_->Fill(i,i);
01118     for(size_t j=i+1;j<hltTrgBits_.size();++j) {
01119       if (hltTrgDeci_.at(i) && hltTrgDeci_.at(j))
01120         h2TrigCorr_->Fill(i,j);
01121       if (hltTrgDeci_.at(i) && !hltTrgDeci_.at(j))
01122         h2TrigCorr_->Fill(j,i);
01123     }
01124   }
01125 }
01126 
01127 
01128 
01129 //--------------------------------------------------------------------------------------------------
01130