CMS 3D CMS Logo

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