Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:08:07

0001 /*
0002     This is the DQM code for UE physics plots
0003     11/12/2009 Sunil Bansal
0004 */
0005 #include "DQM/Physics/src/QcdUeDQM.h"
0006 #include "DataFormats/Common/interface/TriggerResults.h"
0007 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
0008 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
0009 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHitCollection.h"
0010 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0011 #include "DataFormats/VertexReco/interface/Vertex.h"
0012 #include "DQMServices/Core/interface/DQMStore.h"
0013 #include "FWCore/Framework/interface/Event.h"
0014 #include "FWCore/Framework/interface/EventSetup.h"
0015 #include "FWCore/Framework/interface/ESHandle.h"
0016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0017 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0018 #include "FWCore/ServiceRegistry/interface/Service.h"
0019 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
0020 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0021 #include "Geometry/CommonDetUnit/interface/PixelGeomDetType.h"
0022 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0023 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0024 #include "Geometry/TrackerGeometryBuilder/interface/RectangularPixelTopology.h"
0025 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
0026 #include "AnalysisDataFormats/TrackInfo/interface/TrackInfo.h"
0027 #include "AnalysisDataFormats/TrackInfo/interface/TrackInfoTrackAssociation.h"
0028 #include "DataFormats/TrackReco/interface/Track.h"
0029 #include "DataFormats/Candidate/interface/CandidateFwd.h"
0030 #include "DataFormats/Candidate/interface/CandMatchMap.h"
0031 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
0032 #include "DataFormats/Math/interface/deltaPhi.h"
0033 #include <TString.h>
0034 #include <TMath.h>
0035 #include <TH1F.h>
0036 #include <TH2F.h>
0037 #include <TH3F.h>
0038 #include <TProfile.h>
0039 using namespace std;
0040 using namespace edm;
0041 
0042 #define CP(level) if (level >= verbose_)
0043 
0044 struct deleter {
0045   void operator()(TH3F *&h) {
0046     delete h;
0047     h = nullptr;
0048   }
0049 };
0050 
0051 QcdUeDQM::QcdUeDQM(const ParameterSet &parameters)
0052     : hltResName_(parameters.getUntrackedParameter<string>("hltTrgResults")),
0053       verbose_(parameters.getUntrackedParameter<int>("verbose", 3)),
0054       tgeo_(nullptr),
0055       repSumMap_(nullptr),
0056       repSummary_(nullptr),
0057       h2TrigCorr_(nullptr),
0058       ptMin_(parameters.getParameter<double>("ptMin")),
0059       minRapidity_(parameters.getParameter<double>("minRapidity")),
0060       maxRapidity_(parameters.getParameter<double>("maxRapidity")),
0061       tip_(parameters.getParameter<double>("tip")),
0062       lip_(parameters.getParameter<double>("lip")),
0063       diffvtxbs_(parameters.getParameter<double>("diffvtxbs")),
0064       ptErr_pt_(parameters.getParameter<double>("ptErr_pt")),
0065       vtxntk_(parameters.getParameter<double>("vtxntk")),
0066       minHit_(parameters.getParameter<int>("minHit")),
0067       pxlLayerMinCut_(parameters.getParameter<double>("pxlLayerMinCut")),
0068       requirePIX1_(parameters.getParameter<bool>("requirePIX1")),
0069       min3DHit_(parameters.getParameter<int>("min3DHit")),
0070       maxChi2_(parameters.getParameter<double>("maxChi2")),
0071       bsuse_(parameters.getParameter<bool>("bsuse")),
0072       allowTriplets_(parameters.getParameter<bool>("allowTriplets")),
0073       bsPos_(parameters.getParameter<double>("bsPos")),
0074       caloJetLabel_(consumes<reco::CaloJetCollection>(parameters.getUntrackedParameter<edm::InputTag>("caloJetTag"))),
0075       chargedJetLabel_(
0076           consumes<reco::TrackJetCollection>(parameters.getUntrackedParameter<edm::InputTag>("chargedJetTag"))),
0077       trackLabel_(consumes<reco::TrackCollection>(parameters.getUntrackedParameter<edm::InputTag>("trackTag"))),
0078       vtxLabel_(consumes<reco::VertexCollection>(parameters.getUntrackedParameter<edm::InputTag>("vtxTag"))),
0079       bsLabel_(consumes<reco::BeamSpot>(parameters.getParameter<edm::InputTag>("beamSpotTag"))) {
0080   // Constructor.
0081   std::vector<std::string> quality = parameters.getParameter<std::vector<std::string> >("quality");
0082   for (unsigned int j = 0; j < quality.size(); j++)
0083     quality_.push_back(reco::TrackBase::qualityByName(quality[j]));
0084   std::vector<std::string> algorithm = parameters.getParameter<std::vector<std::string> >("algorithm");
0085   for (unsigned int j = 0; j < algorithm.size(); j++)
0086     algorithm_.push_back(reco::TrackBase::algoByName(algorithm[j]));
0087 
0088   if (parameters.exists("hltTrgNames"))
0089     hltTrgNames_ = parameters.getUntrackedParameter<vector<string> >("hltTrgNames");
0090 
0091   if (parameters.exists("hltProcNames"))
0092     hltProcNames_ = parameters.getUntrackedParameter<vector<string> >("hltProcNames");
0093   else {
0094     //     hltProcNames_.push_back("FU");
0095     hltProcNames_.push_back("HLT");
0096   }
0097 
0098   isHltConfigSuccessful_ = false;  // init
0099 }
0100 
0101 QcdUeDQM::~QcdUeDQM() {
0102   // Destructor.
0103 }
0104 
0105 void QcdUeDQM::dqmBeginRun(const Run &run, const EventSetup &iSetup) {
0106   // indicating change of HLT cfg at run boundries
0107   // for HLTConfigProvider::init()
0108   bool isHltCfgChange = false;
0109   isHltConfigSuccessful_ = false;  // init
0110 
0111   string teststr;
0112   for (size_t i = 0; i < hltProcNames_.size(); ++i) {
0113     if (i > 0)
0114       teststr += ", ";
0115     teststr += hltProcNames_.at(i);
0116     if (hltConfig.init(run, iSetup, hltProcNames_.at(i), isHltCfgChange)) {
0117       isHltConfigSuccessful_ = true;
0118       hltUsedResName_ = hltResName_;
0119       if (hltResName_.find(':') == string::npos)
0120         hltUsedResName_ += "::";
0121       else
0122         hltUsedResName_ += ":";
0123       hltUsedResName_ += hltProcNames_.at(i);
0124       break;
0125     }
0126   }
0127 
0128   if (!isHltConfigSuccessful_)
0129     return;
0130 
0131   // setup "Any" bit
0132   hltTrgBits_.clear();
0133   hltTrgBits_.push_back(-1);
0134   hltTrgDeci_.clear();
0135   hltTrgDeci_.push_back(true);
0136   hltTrgUsedNames_.clear();
0137   hltTrgUsedNames_.push_back("Any");
0138 
0139   // figure out relation of trigger name to trigger bit and store used trigger
0140   // names/bits
0141   for (size_t i = 0; i < hltTrgNames_.size(); ++i) {
0142     const string &n1(hltTrgNames_.at(i));
0143     // unsigned int hlt_prescale = hltConfig.prescaleValue(iSetup, n1);
0144     // cout<<"trigger=="<<n1<<"presc=="<<hlt_prescale<<endl;
0145     bool found = false;
0146     for (size_t j = 0; j < hltConfig.size(); ++j) {
0147       const string &n2(hltConfig.triggerName(j));
0148       if (n2 == n1) {
0149         hltTrgBits_.push_back(j);
0150         hltTrgUsedNames_.push_back(n1);
0151         hltTrgDeci_.push_back(false);
0152         found = true;
0153         break;
0154       }
0155     }
0156     if (!found) {
0157       CP(2) cout << "Could not find trigger bit" << endl;
0158     }
0159   }
0160   isHltConfigSuccessful_ = true;
0161 }
0162 
0163 void QcdUeDQM::bookHistograms(DQMStore::IBooker &iBooker, edm::Run const &, edm::EventSetup const &) {
0164   // Book histograms if needed.
0165   iBooker.setCurrentFolder("Physics/QcdUe");
0166 
0167   if (true) {
0168     const int Nx = hltTrgUsedNames_.size();
0169     const double x1 = -0.5;
0170     const double x2 = Nx - 0.5;
0171     h2TrigCorr_ = iBooker.book2D("h2TriCorr", "Trigger bit x vs y;y&&!x;x&&y", Nx, x1, x2, Nx, x1, x2);
0172     for (size_t i = 1; i <= hltTrgUsedNames_.size(); ++i) {
0173       h2TrigCorr_->setBinLabel(i, hltTrgUsedNames_.at(i - 1), 1);
0174       h2TrigCorr_->setBinLabel(i, hltTrgUsedNames_.at(i - 1), 2);
0175     }
0176     TH1 *h = h2TrigCorr_->getTH1();
0177     if (h)
0178       h->SetStats(false);
0179   }
0180   book1D(iBooker, hNevts_, "hNevts", "number of events", 2, 0, 2);
0181   book1D(iBooker, hNtrackerLayer_, "hNtrackerLayer", "number of tracker layers;multiplicity", 20, -0.5, 19.5);
0182   book1D(iBooker, hNtrackerPixelLayer_, "hNtrackerPixelLayer", "number of pixel layers;multiplicity", 10, -0.5, 9.5);
0183   book1D(iBooker,
0184          hNtrackerStripPixelLayer_,
0185          "hNtrackerStripPixelLayer",
0186          "number of strip + pixel layers;multiplicity",
0187          30,
0188          -0.5,
0189          39.5);
0190   book1D(iBooker, hRatioPtErrorPt_, "hRatioPtErrorPt", "ratio of pT error and track pT", 25, 0., 5.);
0191   book1D(iBooker, hTrkPt_, "hTrkPt", "pT of all tracks", 50, 0., 50.);
0192   book1D(iBooker, hTrkEta_, "hTrkEta", "eta of all tracks", 40, -4., 4.);
0193   book1D(iBooker, hTrkPhi_, "hTrkPhi", "phi of all tracks", 40, -4., 4.);
0194   book1D(iBooker,
0195          hRatioDxySigmaDxyBS_,
0196          "hRatioDxySigmaDxyBS",
0197          "ratio of transverse impact parameter and its significance wrt beam spot",
0198          60,
0199          -10.,
0200          10);
0201   book1D(iBooker,
0202          hRatioDxySigmaDxyPV_,
0203          "hRatioDxySigmaDxyPV",
0204          "ratio of transverse impact parameter and its significance wrt PV",
0205          60,
0206          -10.,
0207          10);
0208   book1D(iBooker,
0209          hRatioDzSigmaDzBS_,
0210          "hRatioDzSigmaDzBS",
0211          "ratio of longitudinal impact parameter and its significance wrt beam "
0212          "spot",
0213          80,
0214          -20.,
0215          20);
0216   book1D(iBooker,
0217          hRatioDzSigmaDzPV_,
0218          "hRatioDzSigmaDzPV",
0219          "ratio of longitudinal impact parameter and its significance wrt PV",
0220          80,
0221          -20.,
0222          20);
0223   book1D(iBooker, hTrkChi2_, "hTrkChi2", "track chi2", 30, 0., 30);
0224   book1D(iBooker, hTrkNdof_, "hTrkNdof", "track NDOF", 100, 0, 100);
0225 
0226   book1D(iBooker, hNgoodTrk_, "hNgoodTrk", "number of good tracks", 50, -0.5, 49.5);
0227 
0228   book1D(iBooker, hGoodTrkPt500_, "hGoodTrkPt500", "pT of all good tracks with pT > 500 MeV", 50, 0., 50.);
0229   book1D(iBooker, hGoodTrkEta500_, "hGoodTrkEta500", "eta of all good tracks pT > 500 MeV", 40, -4., 4.);
0230   book1D(iBooker, hGoodTrkPhi500_, "hGoodTrkPhi500", "phi of all good tracks pT > 500 MeV", 40, -4., 4.);
0231 
0232   book1D(iBooker, hGoodTrkPt900_, "hGoodTrkPt900", "pT of all good tracks with pT > 900 MeV", 50, 0., 50.);
0233   book1D(iBooker, hGoodTrkEta900_, "hGoodTrkEta900", "eta of all good tracks pT > 900 MeV", 40, -4., 4.);
0234   book1D(iBooker, hGoodTrkPhi900_, "hGoodTrkPhi900", "phi of all good tracks pT > 900 MeV", 40, -4., 4.);
0235 
0236   book1D(iBooker, hNvertices_, "hNvertices", "number of vertices", 5, -0.5, 4.5);
0237   book1D(iBooker, hVertex_z_, "hVertex_z", "z position of vertex; z[cm]", 200, -50, 50);
0238   book1D(iBooker, hVertex_y_, "hVertex_y", "y position of vertex; y[cm]", 100, -5, 5);
0239   book1D(iBooker, hVertex_x_, "hVertex_x", "x position of vertex; x[cm]", 100, -5, 5);
0240   book1D(iBooker, hVertex_ndof_, "hVertex_ndof", "ndof of vertex", 100, 0, 100);
0241   book1D(iBooker, hVertex_rho_, "hVertex_rho", "rho of vertex", 100, 0, 5);
0242   book1D(iBooker, hVertex_z_bs_, "hVertex_z_bs", "z position of vertex from beamspot; z[cm]", 200, -50, 50);
0243 
0244   book1D(iBooker, hBeamSpot_z_, "hBeamSpot_z", "z position of beamspot; z[cm]", 100, -20, 20);
0245   book1D(iBooker, hBeamSpot_y_, "hBeamSpot_y", "y position of beamspot; y[cm]", 50, -10, 10);
0246   book1D(iBooker, hBeamSpot_x_, "hBeamSpot_x", "x position of beamspot; x[cm]", 50, -10, 10);
0247 
0248   if (true) {
0249     const int Nx = 25;
0250     const double x1 = 0.0;
0251     const double x2 = 50.0;
0252     book1D(iBooker,
0253            hLeadingTrack_pTSpectrum_,
0254            "hLeadingTrack_pTSpectrum",
0255            "pT spectrum of leading track;pT(GeV/c)",
0256            Nx,
0257            x1,
0258            x2);
0259     book1D(iBooker,
0260            hLeadingChargedJet_pTSpectrum_,
0261            "hLeadingChargedJet_pTSpectrum",
0262            "pT spectrum of leading track jet;pT(GeV/c)",
0263            Nx,
0264            x1,
0265            x2);
0266   }
0267 
0268   if (true) {
0269     const int Nx = 24;
0270     const double x1 = -4.;
0271     const double x2 = 4.;
0272     book1D(iBooker,
0273            hLeadingTrack_phiSpectrum_,
0274            "hLeadingTrack_phiSpectrum",
0275            "#phi spectrum of leading track;#phi",
0276            Nx,
0277            x1,
0278            x2);
0279     book1D(iBooker,
0280            hLeadingChargedJet_phiSpectrum_,
0281            "hLeadingChargedJet_phiSpectrum",
0282            "#phi spectrum of leading track jet;#phi",
0283            Nx,
0284            x1,
0285            x2);
0286   }
0287 
0288   if (true) {
0289     const int Nx = 24;
0290     const double x1 = -4.;
0291     const double x2 = 4.;
0292     book1D(iBooker,
0293            hLeadingTrack_etaSpectrum_,
0294            "hLeadingTrack_etaSpectrum",
0295            "#eta spectrum of leading track;#eta",
0296            Nx,
0297            x1,
0298            x2);
0299     book1D(iBooker,
0300            hLeadingChargedJet_etaSpectrum_,
0301            "hLeadingChargedJet_etaSpectrum",
0302            "#eta spectrum of leading track jet;#eta",
0303            Nx,
0304            x1,
0305            x2);
0306   }
0307 
0308   if (true) {
0309     const int Nx = 75;
0310     const double x1 = 0.0;
0311     const double x2 = 75.0;
0312     const double y1 = 0.;
0313     const double y2 = 10.;
0314     bookProfile(iBooker,
0315                 hdNdEtadPhi_pTMax_Toward500_,
0316                 "hdNdEtadPhi_pTMax_Toward500",
0317                 "Average number of tracks (pT > 500 MeV) in toward region vs "
0318                 "leading track pT;pT(GeV/c);dN/d#eta d#phi",
0319                 Nx,
0320                 x1,
0321                 x2,
0322                 y1,
0323                 y2,
0324                 false,
0325                 false);
0326     bookProfile(iBooker,
0327                 hdNdEtadPhi_pTMax_Transverse500_,
0328                 "hdNdEtadPhi_pTMax_Transverse500",
0329                 "Average number of tracks (pT > 500 MeV) in transverse region "
0330                 "vs leading track pT;pT(GeV/c);dN/d#eta d#phi",
0331                 Nx,
0332                 x1,
0333                 x2,
0334                 y1,
0335                 y2,
0336                 false,
0337                 false);
0338     bookProfile(iBooker,
0339                 hdNdEtadPhi_pTMax_Away500_,
0340                 "hdNdEtadPhi_pTMax_Away500",
0341                 "Average number of tracks (pT > 500 MeV) in away region vs "
0342                 "leading track pT;pT(GeV/c);dN/d#eta d#phi",
0343                 Nx,
0344                 x1,
0345                 x2,
0346                 y1,
0347                 y2,
0348                 false,
0349                 false);
0350     bookProfile(iBooker,
0351                 hdNdEtadPhi_trackJet_Toward500_,
0352                 "hdNdEtadPhi_trackJet_Toward500",
0353                 "Average number of tracks (pT > 500 MeV) in toward region vs "
0354                 "leading track jet pT;pT(GeV/c);dN/d#eta d#phi",
0355                 Nx,
0356                 x1,
0357                 x2,
0358                 y1,
0359                 y2);
0360     bookProfile(iBooker,
0361                 hdNdEtadPhi_trackJet_Transverse500_,
0362                 "hdNdEtadPhi_trackJet_Transverse500",
0363                 "Average number of tracks (pT > 500 MeV) in transverse region "
0364                 "vs leading track jet pT;pT(GeV/c);dN/d#eta d#phi",
0365                 Nx,
0366                 x1,
0367                 x2,
0368                 y1,
0369                 y2,
0370                 false,
0371                 false);
0372     bookProfile(iBooker,
0373                 hdNdEtadPhi_trackJet_Away500_,
0374                 "hdNdEtadPhi_trackJet_Away500",
0375                 "Average number of tracks (pT > 500 MeV) in away region vs "
0376                 "leading track jet pT;pT(GeV/c);dN/d#eta d#phi",
0377                 Nx,
0378                 x1,
0379                 x2,
0380                 y1,
0381                 y2,
0382                 false,
0383                 false);
0384 
0385     bookProfile(iBooker,
0386                 hpTSumdEtadPhi_pTMax_Toward500_,
0387                 "hpTSumdEtadPhi_pTMax_Toward500",
0388                 "Average number of tracks (pT > 500 MeV) in toward region vs "
0389                 "leading track pT;pT(GeV/c);dpTSum/d#eta d#phi",
0390                 Nx,
0391                 x1,
0392                 x2,
0393                 y1,
0394                 y2,
0395                 false,
0396                 false);
0397     bookProfile(iBooker,
0398                 hpTSumdEtadPhi_pTMax_Transverse500_,
0399                 "hpTSumdEtadPhi_pTMax_Transverse500",
0400                 "Average number of tracks (pT > 500 MeV) in transverse region "
0401                 "vs leading track pT;pT(GeV/c);dpTSum/d#eta d#phi",
0402                 Nx,
0403                 x1,
0404                 x2,
0405                 y1,
0406                 y2,
0407                 false,
0408                 false);
0409     bookProfile(iBooker,
0410                 hpTSumdEtadPhi_pTMax_Away500_,
0411                 "hpTSumdEtadPhi_pTMax_Away500",
0412                 "Average number of tracks (pT > 500 MeV) in away region vs "
0413                 "leading track pT;pT(GeV/c);dpTSum/d#eta d#phi",
0414                 Nx,
0415                 x1,
0416                 x2,
0417                 y1,
0418                 y2,
0419                 false,
0420                 false);
0421     bookProfile(iBooker,
0422                 hpTSumdEtadPhi_trackJet_Toward500_,
0423                 "hpTSumdEtadPhi_trackJet_Toward500",
0424                 "Average number of tracks (pT > 500 MeV) in toward region vs "
0425                 "leading track jet pT;pT(GeV/c);dpTSum/d#eta d#phi",
0426                 Nx,
0427                 x1,
0428                 x2,
0429                 y1,
0430                 y2,
0431                 false,
0432                 false);
0433     bookProfile(iBooker,
0434                 hpTSumdEtadPhi_trackJet_Transverse500_,
0435                 "hpTSumdEtadPhi_trackJet_Transverse500",
0436                 "Average number of tracks (pT > 500 MeV) in transverse region "
0437                 "vs leading track jet pT;pT(GeV/c);dpTSum/d#eta d#phi",
0438                 Nx,
0439                 x1,
0440                 x2,
0441                 y1,
0442                 y2,
0443                 false,
0444                 false);
0445     bookProfile(iBooker,
0446                 hpTSumdEtadPhi_trackJet_Away500_,
0447                 "hpTSumdEtadPhi_trackJet_Away500",
0448                 "Average number of tracks (pT > 500 MeV) in away region vs "
0449                 "leading track jet pT;pT(GeV/c);dpTSum/d#eta d#phi",
0450                 Nx,
0451                 x1,
0452                 x2,
0453                 y1,
0454                 y2,
0455                 false,
0456                 false);
0457 
0458     bookProfile(iBooker,
0459                 hdNdEtadPhi_pTMax_Toward900_,
0460                 "hdNdEtadPhi_pTMax_Toward900",
0461                 "Average number of tracks (pT > 900 MeV) in toward region vs "
0462                 "leading track pT;pT(GeV/c);dN/d#eta d#phi",
0463                 Nx,
0464                 x1,
0465                 x2,
0466                 y1,
0467                 y2,
0468                 false,
0469                 false);
0470     bookProfile(iBooker,
0471                 hdNdEtadPhi_pTMax_Transverse900_,
0472                 "hdNdEtadPhi_pTMax_Transverse900",
0473                 "Average number of tracks (pT > 900 MeV) in transverse region "
0474                 "vs leading track pT;pT(GeV/c);dN/d#eta d#phi",
0475                 Nx,
0476                 x1,
0477                 x2,
0478                 y1,
0479                 y2,
0480                 false,
0481                 false);
0482     bookProfile(iBooker,
0483                 hdNdEtadPhi_pTMax_Away900_,
0484                 "hdNdEtadPhi_pTMax_Away900",
0485                 "Average number of tracks (pT > 900 MeV) in away region vs "
0486                 "leading track pT;pT(GeV/c);dN/d#eta d#phi",
0487                 Nx,
0488                 x1,
0489                 x2,
0490                 y1,
0491                 y2,
0492                 false,
0493                 false);
0494     bookProfile(iBooker,
0495                 hdNdEtadPhi_trackJet_Toward900_,
0496                 "hdNdEtadPhi_trackJet_Toward900",
0497                 "Average number of tracks (pT > 900 MeV) in toward region vs "
0498                 "leading track jet pT;pT(GeV/c);dN/d#eta d#phi",
0499                 Nx,
0500                 x1,
0501                 x2,
0502                 y1,
0503                 y2);
0504     bookProfile(iBooker,
0505                 hdNdEtadPhi_trackJet_Transverse900_,
0506                 "hdNdEtadPhi_trackJet_Transverse900",
0507                 "Average number of tracks (pT > 900 MeV) in transverse region "
0508                 "vs leading track jet pT;pT(GeV/c);dN/d#eta d#phi",
0509                 Nx,
0510                 x1,
0511                 x2,
0512                 y1,
0513                 y2,
0514                 false,
0515                 false);
0516     bookProfile(iBooker,
0517                 hdNdEtadPhi_trackJet_Away900_,
0518                 "hdNdEtadPhi_trackJet_Away900",
0519                 "Average number of tracks (pT > 900 MeV) in away region vs "
0520                 "leading track jet pT;pT(GeV/c);dN/d#eta d#phi",
0521                 Nx,
0522                 x1,
0523                 x2,
0524                 y1,
0525                 y2,
0526                 false,
0527                 false);
0528 
0529     bookProfile(iBooker,
0530                 hpTSumdEtadPhi_pTMax_Toward900_,
0531                 "hpTSumdEtadPhi_pTMax_Toward900",
0532                 "Average number of tracks (pT > 900 MeV) in toward region vs "
0533                 "leading track pT;pT(GeV/c);dpTSum/d#eta d#phi",
0534                 Nx,
0535                 x1,
0536                 x2,
0537                 y1,
0538                 y2,
0539                 false,
0540                 false);
0541     bookProfile(iBooker,
0542                 hpTSumdEtadPhi_pTMax_Transverse900_,
0543                 "hpTSumdEtadPhi_pTMax_Transverse900",
0544                 "Average number of tracks (pT > 900 MeV) in transverse region "
0545                 "vs leading track pT;pT(GeV/c);dpTSum/d#eta d#phi",
0546                 Nx,
0547                 x1,
0548                 x2,
0549                 y1,
0550                 y2,
0551                 false,
0552                 false);
0553     bookProfile(iBooker,
0554                 hpTSumdEtadPhi_pTMax_Away900_,
0555                 "hpTSumdEtadPhi_pTMax_Away900",
0556                 "Average number of tracks (pT > 900 MeV) in away region vs "
0557                 "leading track pT;pT(GeV/c);dpTSum/d#eta d#phi",
0558                 Nx,
0559                 x1,
0560                 x2,
0561                 y1,
0562                 y2,
0563                 false,
0564                 false);
0565     bookProfile(iBooker,
0566                 hpTSumdEtadPhi_trackJet_Toward900_,
0567                 "hpTSumdEtadPhi_trackJet_Toward900",
0568                 "Average number of tracks (pT > 900 MeV) in toward region vs "
0569                 "leading track jet pT;pT(GeV/c);dpTSum/d#eta d#phi",
0570                 Nx,
0571                 x1,
0572                 x2,
0573                 y1,
0574                 y2,
0575                 false,
0576                 false);
0577     bookProfile(iBooker,
0578                 hpTSumdEtadPhi_trackJet_Transverse900_,
0579                 "hpTSumdEtadPhi_trackJet_Transverse900",
0580                 "Average number of tracks (pT > 900 MeV) in transverse region "
0581                 "vs leading track jet pT;pT(GeV/c);dpTSum/d#eta d#phi",
0582                 Nx,
0583                 x1,
0584                 x2,
0585                 y1,
0586                 y2,
0587                 false,
0588                 false);
0589     bookProfile(iBooker,
0590                 hpTSumdEtadPhi_trackJet_Away900_,
0591                 "hpTSumdEtadPhi_trackJet_Away900",
0592                 "Average number of tracks (pT > 900 MeV) in away region vs "
0593                 "leading track jet pT;pT(GeV/c);dpTSum/d#eta d#phi",
0594                 Nx,
0595                 x1,
0596                 x2,
0597                 y1,
0598                 y2,
0599                 false,
0600                 false);
0601   }
0602 
0603   if (true) {
0604     const int Nx = 20;
0605     const double x1 = 0.0;
0606     const double x2 = 20.0;
0607 
0608     book1D(iBooker, hChargedJetMulti_, "hChargedJetMulti", "Charged jet multiplicity;multiplicities", Nx, x1, x2);
0609   }
0610 
0611   if (true) {
0612     const int Nx = 60;
0613     const double x1 = -180.0;
0614     const double x2 = 180.0;
0615 
0616     book1D(iBooker,
0617            hdPhi_maxpTTrack_tracks_,
0618            "hdPhi_maxpTTrack_tracks",
0619            "delta phi between leading tracks and other "
0620            "tracks;#Delta#phi(leading track-track)",
0621            Nx,
0622            x1,
0623            x2);
0624     book1D(iBooker,
0625            hdPhi_chargedJet_tracks_,
0626            "hdPhi_chargedJet_tracks",
0627            "delta phi between leading charged jet  and "
0628            "tracks;#Delta#phi(leading charged jet-track)",
0629            Nx,
0630            x1,
0631            x2);
0632   }
0633 }
0634 
0635 void QcdUeDQM::analyze(const Event &iEvent, const EventSetup &iSetup) {
0636   if (!isHltConfigSuccessful_)
0637     return;
0638 
0639   // Analyze the given event.
0640 
0641   edm::Handle<reco::BeamSpot> beamSpot;
0642   bool ValidBS_ = iEvent.getByToken(bsLabel_, beamSpot);
0643   if (!ValidBS_)
0644     return;
0645 
0646   edm::Handle<reco::TrackCollection> tracks;
0647   bool ValidTrack_ = iEvent.getByToken(trackLabel_, tracks);
0648   if (!ValidTrack_)
0649     return;
0650 
0651   edm::Handle<reco::TrackJetCollection> trkJets;
0652   bool ValidTrackJet_ = iEvent.getByToken(chargedJetLabel_, trkJets);
0653   if (!ValidTrackJet_)
0654     return;
0655 
0656   edm::Handle<reco::CaloJetCollection> calJets;
0657   bool ValidCaloJet_ = iEvent.getByToken(caloJetLabel_, calJets);
0658   if (!ValidCaloJet_)
0659     return;
0660 
0661   edm::Handle<reco::VertexCollection> vertexColl;
0662   bool ValidVtxColl_ = iEvent.getByToken(vtxLabel_, vertexColl);
0663   if (!ValidVtxColl_)
0664     return;
0665 
0666   reco::TrackCollection tracks_sort = *tracks;
0667   std::sort(tracks_sort.begin(), tracks_sort.end(), PtSorter());
0668 
0669   // get tracker geometry
0670   /*  ESHandle<TrackerGeometry> trackerHandle;
0671    iSetup.get<TrackerDigiGeometryRecord>().get(trackerHandle);
0672    tgeo_ = trackerHandle.product();
0673    if (!tgeo_)return;
0674  */
0675   selected_.clear();
0676   fillHltBits(iEvent, iSetup);
0677   // select good tracks
0678   if (fillVtxPlots(beamSpot.product(), vertexColl)) {
0679     fill1D(hNevts_, 1);
0680     for (reco::TrackCollection::const_iterator Trk = tracks_sort.begin(); Trk != tracks_sort.end(); ++Trk) {
0681       if (trackSelection(*Trk, beamSpot.product(), vtx1, vertexColl->size()))
0682         selected_.push_back(&*Trk);
0683     }
0684 
0685     fillpTMaxRelated(selected_);
0686     fillChargedJetSpectra(trkJets);
0687     //    fillCaloJetSpectra(calJets);
0688     fillUE_with_MaxpTtrack(selected_);
0689     if (!trkJets->empty())
0690       fillUE_with_ChargedJets(selected_, trkJets);
0691     // if(calJets->size()>0)fillUE_with_CaloJets(selected_,calJets);
0692   }
0693 }
0694 
0695 void QcdUeDQM::book1D(DQMStore::IBooker &iBooker,
0696                       std::vector<MonitorElement *> &mes,
0697                       const std::string &name,
0698                       const std::string &title,
0699                       int nx,
0700                       double x1,
0701                       double x2,
0702                       bool sumw2,
0703                       bool sbox) {
0704   // Book 1D histos.
0705   for (size_t i = 0; i < hltTrgUsedNames_.size(); ++i) {
0706     std::string folderName = "Physics/QcdUe/" + hltTrgUsedNames_.at(i);
0707     iBooker.setCurrentFolder(folderName);
0708     MonitorElement *e = iBooker.book1D(Form("%s_%s", name.c_str(), hltTrgUsedNames_.at(i).c_str()),
0709                                        Form("%s: %s", hltTrgUsedNames_.at(i).c_str(), title.c_str()),
0710                                        nx,
0711                                        x1,
0712                                        x2);
0713     TH1 *h1 = e->getTH1();
0714     if (sumw2) {
0715       if (0 == h1->GetSumw2N()) {  // protect against re-summing (would cause
0716                                    // exception)
0717         h1->Sumw2();
0718       }
0719     }
0720     h1->SetStats(sbox);
0721     mes.push_back(e);
0722   }
0723 }
0724 
0725 void QcdUeDQM::bookProfile(DQMStore::IBooker &iBooker,
0726                            std::vector<MonitorElement *> &mes,
0727                            const std::string &name,
0728                            const std::string &title,
0729                            int nx,
0730                            double x1,
0731                            double x2,
0732                            double y1,
0733                            double y2,
0734                            bool sumw2,
0735                            bool sbox) {
0736   // Book Profile histos.
0737 
0738   for (size_t i = 0; i < hltTrgUsedNames_.size(); ++i) {
0739     std::string folderName = "Physics/QcdUe/" + hltTrgUsedNames_.at(i);
0740     iBooker.setCurrentFolder(folderName);
0741     MonitorElement *e = iBooker.bookProfile(Form("%s_%s", name.c_str(), hltTrgUsedNames_.at(i).c_str()),
0742                                             Form("%s: %s", hltTrgUsedNames_.at(i).c_str(), title.c_str()),
0743                                             nx,
0744                                             x1,
0745                                             x2,
0746                                             y1,
0747                                             y2,
0748                                             " ");
0749     mes.push_back(e);
0750   }
0751 }
0752 
0753 void QcdUeDQM::fill1D(std::vector<TH1F *> &hs, double val, double w) {
0754   // Loop over histograms and fill if trigger has fired.
0755 
0756   for (size_t i = 0; i < hs.size(); ++i) {
0757     if (!hltTrgDeci_.at(i))
0758       continue;
0759     hs.at(i)->Fill(val, w);
0760   }
0761 }
0762 
0763 //--------------------------------------------------------------------------------------------------
0764 void QcdUeDQM::fill1D(std::vector<MonitorElement *> &mes, double val, double w) {
0765   // Loop over histograms and fill if trigger has fired.
0766 
0767   for (size_t i = 0; i < mes.size(); ++i) {
0768     if (!hltTrgDeci_.at(i))
0769       continue;
0770     mes.at(i)->Fill(val, w);
0771   }
0772 }
0773 
0774 //--------------------------------------------------------------------------------------------------
0775 void QcdUeDQM::setLabel1D(std::vector<MonitorElement *> &mes) {
0776   // Loop over histograms and fill if trigger has fired.
0777   string cut[5] = {"Nevt", "vtx!=bmspt", "Zvtx<10cm", "pT>1GeV", "trackFromVtx"};
0778   for (size_t i = 0; i < mes.size(); ++i) {
0779     if (!hltTrgDeci_.at(i))
0780       continue;
0781     for (size_t j = 1; j < 6; j++)
0782       mes.at(i)->setBinLabel(j, cut[j - 1], 1);
0783   }
0784 }
0785 
0786 //--------------------------------------------------------------------------------------------------
0787 void QcdUeDQM::fill2D(std::vector<TH2F *> &hs, double valx, double valy, double w) {
0788   // Loop over histograms and fill if trigger has fired.
0789 
0790   for (size_t i = 0; i < hs.size(); ++i) {
0791     if (!hltTrgDeci_.at(i))
0792       continue;
0793     hs.at(i)->Fill(valx, valy, w);
0794   }
0795 }
0796 
0797 //--------------------------------------------------------------------------------------------------
0798 void QcdUeDQM::fill2D(std::vector<MonitorElement *> &mes, double valx, double valy, double w) {
0799   // Loop over histograms and fill if trigger has fired.
0800 
0801   for (size_t i = 0; i < mes.size(); ++i) {
0802     if (!hltTrgDeci_.at(i))
0803       continue;
0804     mes.at(i)->Fill(valx, valy, w);
0805   }
0806 }
0807 //--------------------------------------------------------------------------------------------------
0808 void QcdUeDQM::fillProfile(std::vector<TProfile *> &hs, double valx, double valy, double w) {
0809   // Loop over histograms and fill if trigger has fired.
0810 
0811   for (size_t i = 0; i < hs.size(); ++i) {
0812     if (!hltTrgDeci_.at(i))
0813       continue;
0814     hs.at(i)->Fill(valx, valy, w);
0815   }
0816 }
0817 
0818 //--------------------------------------------------------------------------------------------------
0819 void QcdUeDQM::fillProfile(std::vector<MonitorElement *> &mes, double valx, double valy, double w) {
0820   // Loop over histograms and fill if trigger has fired.
0821 
0822   for (size_t i = 0; i < mes.size(); ++i) {
0823     if (!hltTrgDeci_.at(i))
0824       continue;
0825     const double y = valy * w;
0826     mes.at(i)->Fill(valx, y);
0827   }
0828 }
0829 
0830 //--------------------------------------------------------------------------------------------------
0831 bool QcdUeDQM::trackSelection(const reco::Track &trk, const reco::BeamSpot *bs, const reco::Vertex &vtx, int sizevtx) {
0832   bool goodTrk = false;
0833 
0834   if (sizevtx != 1)
0835     return false;  // selection events with only a vertex
0836 
0837   // Fill basic information of all the tracks
0838   fill1D(hNtrackerLayer_, trk.hitPattern().trackerLayersWithMeasurement());
0839   fill1D(hNtrackerPixelLayer_, trk.hitPattern().pixelLayersWithMeasurement());
0840 
0841   fill1D(
0842       hNtrackerStripPixelLayer_,
0843       (trk.hitPattern().pixelLayersWithMeasurement() + trk.hitPattern().numberOfValidStripLayersWithMonoAndStereo()));
0844 
0845   fill1D(hRatioPtErrorPt_, (trk.ptError() / trk.pt()));
0846   fill1D(hTrkPt_, trk.pt());
0847   fill1D(hTrkEta_, trk.eta());
0848   fill1D(hTrkPhi_, trk.phi());
0849   fill1D(hRatioDxySigmaDxyBS_, (trk.dxy(bs->position()) / trk.dxyError()));
0850   fill1D(hRatioDxySigmaDxyPV_, (trk.dxy(vtx.position()) / trk.dxyError()));
0851   fill1D(hRatioDzSigmaDzBS_, (trk.dz(bs->position()) / trk.dzError()));
0852   fill1D(hRatioDzSigmaDzPV_, (trk.dz(vtx.position()) / trk.dzError()));
0853   fill1D(hTrkChi2_, trk.normalizedChi2());
0854   fill1D(hTrkNdof_, trk.ndof());
0855 
0856   fill1D(hBeamSpot_x_, bs->x0());
0857   fill1D(hBeamSpot_y_, bs->y0());
0858   fill1D(hBeamSpot_z_, bs->z0());
0859 
0860   // number of layers
0861   bool layerMinCutbool = (trk.hitPattern().trackerLayersWithMeasurement() >= minHit_ ||
0862                           (trk.hitPattern().trackerLayersWithMeasurement() == 3 &&
0863                            trk.hitPattern().pixelLayersWithMeasurement() == 3 && allowTriplets_));
0864 
0865   // number of pixel layers
0866   bool pxlLayerMinCutbool = (trk.hitPattern().pixelLayersWithMeasurement() >= pxlLayerMinCut_);
0867 
0868   // cut on the hits in pixel layers
0869   bool hasPIX1 = false;
0870   if (requirePIX1_) {
0871     const reco::HitPattern &p = trk.hitPattern();
0872     for (int i = 0; i < p.numberOfAllHits(reco::HitPattern::TRACK_HITS); i++) {
0873       uint32_t hit = p.getHitPattern(reco::HitPattern::TRACK_HITS, i);
0874       if (reco::HitPattern::validHitFilter(hit) && reco::HitPattern::pixelHitFilter(hit) &&
0875           reco::HitPattern::getLayer(hit) == 1) {
0876         hasPIX1 = true;
0877         break;
0878       }
0879     }
0880   } else {
0881     hasPIX1 = true;
0882   }
0883 
0884   // cut on the pT error
0885   bool ptErrorbool =
0886       (trk.ptError() / trk.pt() < ptErr_pt_ || (trk.hitPattern().trackerLayersWithMeasurement() == 3 &&
0887                                                 trk.hitPattern().pixelLayersWithMeasurement() == 3 && allowTriplets_));
0888 
0889   // quality cut
0890   bool quality_ok = true;
0891   if (!quality_.empty()) {
0892     quality_ok = false;
0893     for (unsigned int i = 0; i < quality_.size(); ++i) {
0894       if (trk.quality(quality_[i])) {
0895         quality_ok = true;
0896         break;
0897       }
0898     }
0899   }
0900   //-----
0901   bool algo_ok = true;
0902   if (!algorithm_.empty()) {
0903     if (std::find(algorithm_.begin(), algorithm_.end(), trk.algo()) == algorithm_.end())
0904       algo_ok = false;
0905   }
0906 
0907   if (bsuse_ == 1) {
0908     if (hasPIX1 && pxlLayerMinCutbool && layerMinCutbool &&
0909         (trk.hitPattern().pixelLayersWithMeasurement() +
0910          trk.hitPattern().numberOfValidStripLayersWithMonoAndStereo()) >= min3DHit_ &&
0911         ptErrorbool && fabs(trk.pt()) >= ptMin_ && trk.eta() >= minRapidity_ && trk.eta() <= maxRapidity_ &&
0912         fabs(trk.dxy(bs->position()) / trk.dxyError()) < tip_ && fabs(trk.dz(bs->position()) / trk.dzError()) < lip_ &&
0913         trk.normalizedChi2() <= maxChi2_ && quality_ok && algo_ok) {
0914       goodTrk = true;
0915     }
0916   }
0917 
0918   if (bsuse_ == 0) {
0919     if (hasPIX1 && pxlLayerMinCutbool && layerMinCutbool &&
0920         (trk.hitPattern().pixelLayersWithMeasurement() +
0921          trk.hitPattern().numberOfValidStripLayersWithMonoAndStereo()) >= min3DHit_ &&
0922         ptErrorbool && fabs(trk.pt()) >= ptMin_ && trk.eta() >= minRapidity_ && trk.eta() <= maxRapidity_ &&
0923         fabs(trk.dxy(vtx.position()) / trk.dxyError()) < tip_ && fabs(trk.dz(vtx.position()) / trk.dzError()) < lip_ &&
0924         trk.normalizedChi2() <= maxChi2_ && quality_ok && algo_ok) {
0925       goodTrk = true;
0926     }
0927   }
0928 
0929   return goodTrk;
0930 }
0931 
0932 //--------------------------------------------------------------------------------------------------
0933 bool QcdUeDQM::fillVtxPlots(const reco::BeamSpot *bs, const edm::Handle<reco::VertexCollection> vtxColl) {
0934   const reco::VertexCollection theVertices = *(vtxColl.product());
0935   bool goodVtx = false;
0936   fill1D(hNvertices_, theVertices.size());
0937   for (reco::VertexCollection::const_iterator vertexIt = theVertices.begin(); vertexIt != theVertices.end();
0938        ++vertexIt) {
0939     fill1D(hVertex_z_, vertexIt->z());
0940     fill1D(hVertex_y_, vertexIt->y());
0941     fill1D(hVertex_x_, vertexIt->x());
0942     fill1D(hVertex_ndof_, vertexIt->ndof());
0943     fill1D(hVertex_rho_, vertexIt->position().rho());
0944     fill1D(hVertex_z_bs_, (vertexIt->z() - bs->z0()));
0945 
0946     if (fabs(vertexIt->z() - bs->z0()) < diffvtxbs_ && vertexIt->ndof() >= 4 && vertexIt->position().rho() <= 2.0) {
0947       goodVtx = true;
0948       vtx1 = (*vertexIt);
0949 
0950       break;
0951     }
0952   }  // Loop over vertcies
0953   return goodVtx;
0954 }
0955 //--------------------------------------------------------------------------------------------------
0956 void QcdUeDQM::fillpTMaxRelated(const std::vector<const reco::Track *> &track) {
0957   fill1D(hNgoodTrk_, track.size());
0958   if (!track.empty()) {
0959     fill1D(hLeadingTrack_pTSpectrum_, track[0]->pt());
0960     fill1D(hLeadingTrack_phiSpectrum_, track[0]->phi());
0961     fill1D(hLeadingTrack_etaSpectrum_, track[0]->eta());
0962   }
0963   for (size_t i = 0; i < track.size(); i++) {
0964     fill1D(hGoodTrkPt500_, track[i]->pt());
0965     fill1D(hGoodTrkEta500_, track[i]->eta());
0966     fill1D(hGoodTrkPhi500_, track[i]->phi());
0967     if (track[i]->pt() > 0.9) {
0968       fill1D(hGoodTrkPt900_, track[i]->pt());
0969       fill1D(hGoodTrkEta900_, track[i]->eta());
0970       fill1D(hGoodTrkPhi900_, track[i]->phi());
0971     }
0972   }
0973 }
0974 
0975 void QcdUeDQM::fillChargedJetSpectra(const edm::Handle<reco::TrackJetCollection> trackJets) {
0976   fill1D(hChargedJetMulti_, trackJets->size());
0977   for (reco::TrackJetCollection::const_iterator f = trackJets->begin(); f != trackJets->end(); f++) {
0978     if (f != trackJets->begin())
0979       continue;
0980     fill1D(hLeadingChargedJet_pTSpectrum_, f->pt());
0981     fill1D(hLeadingChargedJet_etaSpectrum_, f->eta());
0982     fill1D(hLeadingChargedJet_phiSpectrum_, f->phi());
0983   }
0984 }
0985 
0986 /*
0987 void QcdUeDQM::fillCaloJetSpectra(const edm::Handle<reco::CaloJetCollection>
0988 caloJets)
0989 {
0990   fill1D(hCaloJetMulti_,caloJets->size());
0991    for( reco::CaloJetCollection::const_iterator f  = caloJets->begin();  f !=
0992 caloJets->end(); f++)
0993      {
0994        if(f != caloJets->begin())continue;
0995        fill1D(hLeadingCaloJet_pTSpectrum_,f->pt());
0996        fill1D(hLeadingCaloJet_etaSpectrum_,f->eta());
0997        fill1D(hLeadingCaloJet_phiSpectrum_,f->phi());
0998      }
0999 
1000 }
1001 
1002 
1003 */
1004 
1005 /*
1006  weight for transverse/toward/away region = 0.12
1007 
1008 
1009 */
1010 
1011 void QcdUeDQM::fillUE_with_MaxpTtrack(const std::vector<const reco::Track *> &track) {
1012   double w = 0.119;
1013   // double w = 1.;
1014   double nTrk500_TransReg = 0;
1015   double nTrk500_AwayReg = 0;
1016   double nTrk500_TowardReg = 0;
1017 
1018   double pTSum500_TransReg = 0;
1019   double pTSum500_AwayReg = 0;
1020   double pTSum500_TowardReg = 0;
1021 
1022   double nTrk900_TransReg = 0;
1023   double nTrk900_AwayReg = 0;
1024   double nTrk900_TowardReg = 0;
1025 
1026   double pTSum900_TransReg = 0;
1027   double pTSum900_AwayReg = 0;
1028   double pTSum900_TowardReg = 0;
1029   if (!track.empty()) {
1030     if (track[0]->pt() > 1.) {
1031       for (size_t i = 1; i < track.size(); i++) {
1032         double dphi = (180. / PI) * (deltaPhi(track[0]->phi(), track[i]->phi()));
1033         fill1D(hdPhi_maxpTTrack_tracks_, dphi);
1034         if (fabs(dphi) > 60. && fabs(dphi) < 120.) {
1035           pTSum500_TransReg = pTSum500_TransReg + track[i]->pt();
1036           nTrk500_TransReg++;
1037           if (track[i]->pt() > 0.9) {
1038             pTSum900_TransReg = pTSum900_TransReg + track[i]->pt();
1039             nTrk900_TransReg++;
1040           }
1041         }
1042 
1043         if (fabs(dphi) > 120. && fabs(dphi) < 180.) {
1044           pTSum500_AwayReg = pTSum500_AwayReg + track[i]->pt();
1045           nTrk500_AwayReg++;
1046           if (track[i]->pt() > 0.9) {
1047             pTSum900_AwayReg = pTSum900_AwayReg + track[i]->pt();
1048             nTrk900_AwayReg++;
1049           }
1050         }
1051 
1052         if (fabs(dphi) < 60.) {
1053           pTSum500_TowardReg = pTSum500_TowardReg + track[i]->pt();
1054           nTrk500_TowardReg++;
1055           if (track[i]->pt() > 0.9) {
1056             pTSum900_TowardReg = pTSum900_TowardReg + track[i]->pt();
1057             nTrk900_TowardReg++;
1058           }
1059         }
1060       }  // track loop
1061     }    // leading track
1062          // non empty collection
1063     fillProfile(hdNdEtadPhi_pTMax_Toward500_, track[0]->pt(), nTrk500_TowardReg, w);
1064     fillProfile(hdNdEtadPhi_pTMax_Transverse500_, track[0]->pt(), nTrk500_TransReg, w);
1065     fillProfile(hdNdEtadPhi_pTMax_Away500_, track[0]->pt(), nTrk500_AwayReg, w);
1066 
1067     fillProfile(hpTSumdEtadPhi_pTMax_Toward500_, track[0]->pt(), pTSum500_TowardReg, w);
1068     fillProfile(hpTSumdEtadPhi_pTMax_Transverse500_, track[0]->pt(), pTSum500_TransReg, w);
1069     fillProfile(hpTSumdEtadPhi_pTMax_Away500_, track[0]->pt(), pTSum500_AwayReg, w);
1070 
1071     fillProfile(hdNdEtadPhi_pTMax_Toward900_, track[0]->pt(), nTrk900_TowardReg, w);
1072     fillProfile(hdNdEtadPhi_pTMax_Transverse900_, track[0]->pt(), nTrk900_TransReg, w);
1073     fillProfile(hdNdEtadPhi_pTMax_Away900_, track[0]->pt(), nTrk900_AwayReg, w);
1074 
1075     fillProfile(hpTSumdEtadPhi_pTMax_Toward900_, track[0]->pt(), pTSum900_TowardReg, w);
1076     fillProfile(hpTSumdEtadPhi_pTMax_Transverse900_, track[0]->pt(), pTSum900_TransReg, w);
1077     fillProfile(hpTSumdEtadPhi_pTMax_Away900_, track[0]->pt(), pTSum900_AwayReg, w);
1078   }
1079 }
1080 
1081 void QcdUeDQM::fillUE_with_ChargedJets(const std::vector<const reco::Track *> &track,
1082                                        const edm::Handle<reco::TrackJetCollection> &trackJets) {
1083   double w = 0.119;
1084   double nTrk500_TransReg = 0;
1085   double nTrk500_AwayReg = 0;
1086   double nTrk500_TowardReg = 0;
1087 
1088   double pTSum500_TransReg = 0;
1089   double pTSum500_AwayReg = 0;
1090   double pTSum500_TowardReg = 0;
1091 
1092   double nTrk900_TransReg = 0;
1093   double nTrk900_AwayReg = 0;
1094   double nTrk900_TowardReg = 0;
1095 
1096   double pTSum900_TransReg = 0;
1097   double pTSum900_AwayReg = 0;
1098   double pTSum900_TowardReg = 0;
1099 
1100   if (!(trackJets->empty()) && (trackJets->begin())->pt() > 1.) {
1101     double jetPhi = (trackJets->begin())->phi();
1102     for (size_t i = 0; i < track.size(); i++) {
1103       double dphi = (180. / PI) * (deltaPhi(jetPhi, track[i]->phi()));
1104       fill1D(hdPhi_chargedJet_tracks_, dphi);
1105       if (fabs(dphi) > 60. && fabs(dphi) < 120.) {
1106         pTSum500_TransReg = pTSum500_TransReg + track[i]->pt();
1107         nTrk500_TransReg++;
1108         if (track[i]->pt() > 0.9) {
1109           pTSum900_TransReg = pTSum900_TransReg + track[i]->pt();
1110           nTrk900_TransReg++;
1111         }
1112       }
1113 
1114       if (fabs(dphi) > 120. && fabs(dphi) < 180.) {
1115         pTSum500_AwayReg = pTSum500_AwayReg + track[i]->pt();
1116         nTrk500_AwayReg++;
1117         if (track[i]->pt() > 0.9) {
1118           pTSum900_AwayReg = pTSum900_AwayReg + track[i]->pt();
1119           nTrk900_AwayReg++;
1120         }
1121       }
1122       if (fabs(dphi) < 60.) {
1123         pTSum500_TowardReg = pTSum500_TowardReg + track[i]->pt();
1124         nTrk500_TowardReg++;
1125         if (track[i]->pt() > 0.9) {
1126           pTSum900_TowardReg = pTSum900_TowardReg + track[i]->pt();
1127           nTrk900_TowardReg++;
1128         }
1129       }
1130     }  // tracks loop
1131 
1132   }  // leading track jet
1133 
1134   fillProfile(hdNdEtadPhi_trackJet_Toward500_, (trackJets->begin())->pt(), nTrk500_TowardReg, w);
1135   fillProfile(hdNdEtadPhi_trackJet_Transverse500_, (trackJets->begin())->pt(), nTrk500_TransReg, w);
1136   fillProfile(hdNdEtadPhi_trackJet_Away500_, (trackJets->begin())->pt(), nTrk500_AwayReg, w);
1137 
1138   fillProfile(hpTSumdEtadPhi_trackJet_Toward500_, (trackJets->begin())->pt(), pTSum500_TowardReg, w);
1139   fillProfile(hpTSumdEtadPhi_trackJet_Transverse500_, (trackJets->begin())->pt(), pTSum500_TransReg, w);
1140   fillProfile(hpTSumdEtadPhi_trackJet_Away500_, (trackJets->begin())->pt(), pTSum500_AwayReg, w);
1141 
1142   fillProfile(hdNdEtadPhi_trackJet_Toward900_, (trackJets->begin())->pt(), nTrk900_TowardReg, w);
1143   fillProfile(hdNdEtadPhi_trackJet_Transverse900_, (trackJets->begin())->pt(), nTrk900_TransReg, w);
1144   fillProfile(hdNdEtadPhi_trackJet_Away900_, (trackJets->begin())->pt(), nTrk900_AwayReg, w);
1145 
1146   fillProfile(hpTSumdEtadPhi_trackJet_Toward900_, (trackJets->begin())->pt(), pTSum900_TowardReg, w);
1147   fillProfile(hpTSumdEtadPhi_trackJet_Transverse900_, (trackJets->begin())->pt(), pTSum900_TransReg, w);
1148   fillProfile(hpTSumdEtadPhi_trackJet_Away900_, (trackJets->begin())->pt(), pTSum900_AwayReg, w);
1149 }
1150 
1151 /*
1152 void QcdUeDQM:: fillUE_with_CaloJets(const std::vector<const reco::Track *>
1153 &track, const edm::Handle<reco::CaloJetCollection> &caloJets)
1154 {
1155 double w = 0.119;
1156 double nTrk500_TransReg = 0;
1157 double nTrk500_AwayReg = 0;
1158 double nTrk500_TowardReg = 0;
1159 
1160 double pTSum500_TransReg = 0;
1161 double pTSum500_AwayReg = 0;
1162 double pTSum500_TowardReg = 0;
1163 
1164 double nTrk900_TransReg = 0;
1165 double nTrk900_AwayReg = 0;
1166 double nTrk900_TowardReg = 0;
1167 
1168 double pTSum900_TransReg = 0;
1169 double pTSum900_AwayReg = 0;
1170 double pTSum900_TowardReg = 0;
1171     if(!(caloJets->empty()) && (caloJets->begin())->pt() > 1.)
1172           {
1173             double jetPhi = (caloJets->begin())->phi();
1174              for(size_t i = 0; i < track.size();i++)
1175                    {
1176                          double dphi =
1177 (180./PI)*(deltaPhi(jetPhi,track[i]->phi()));
1178                          fill1D(hdPhi_caloJet_tracks_,dphi);
1179                          if(fabs(dphi)>60. && fabs(dphi)<120.)
1180                            {
1181                                  pTSum500_TransReg =  pTSum500_TransReg +
1182 track[i]->pt();
1183                                  nTrk500_TransReg++;
1184                                  if(track[i]->pt() > 0.9)
1185                                  {
1186                                  pTSum900_TransReg =  pTSum900_TransReg +
1187 track[i]->pt();
1188                                  nTrk900_TransReg++;
1189                                  }
1190                            }
1191                          if(fabs(dphi)>120. && fabs(dphi)<180.)
1192                            {
1193                                  pTSum500_AwayReg =  pTSum500_AwayReg +
1194 track[i]->pt();
1195                                  nTrk500_AwayReg++;
1196                                  if(track[i]->pt() > 0.9)
1197                                  {
1198                                  pTSum900_AwayReg =  pTSum900_AwayReg +
1199 track[i]->pt();
1200                                  nTrk900_AwayReg++;
1201                                  }
1202                            }
1203                          if(fabs(dphi)<60.)
1204                            {
1205                                  pTSum500_TowardReg =  pTSum500_TowardReg +
1206 track[i]->pt();
1207                                  nTrk500_TowardReg++;
1208                                  if(track[i]->pt() > 0.9)
1209                                  {
1210                                  pTSum900_TowardReg =  pTSum900_TowardReg +
1211 track[i]->pt();
1212                                  nTrk900_TowardReg++;
1213                                  }
1214                            }
1215                        }// tracks loop
1216 
1217                    }// leading calo jet
1218                  fillProfile(hdNdEtadPhi_caloJet_Toward500_,
1219 (caloJets->begin())->pt(),nTrk500_TowardReg,w);
1220                  fillProfile(hdNdEtadPhi_caloJet_Transverse500_,
1221 (caloJets->begin())->pt(),nTrk500_TransReg,w);
1222                  fillProfile(hdNdEtadPhi_caloJet_Away500_,
1223 (caloJets->begin())->pt(),nTrk500_AwayReg,w);
1224 
1225                  fillProfile(hpTSumdEtadPhi_caloJet_Toward500_,
1226 (caloJets->begin())->pt(),pTSum500_TowardReg,w);
1227                  fillProfile(hpTSumdEtadPhi_caloJet_Transverse500_,
1228 (caloJets->begin())->pt(),pTSum500_TransReg,w);
1229                  fillProfile(hpTSumdEtadPhi_caloJet_Away500_,
1230 (caloJets->begin())->pt(),pTSum500_AwayReg,w);
1231 
1232                  fillProfile(hdNdEtadPhi_caloJet_Toward900_,
1233 (caloJets->begin())->pt(),nTrk900_TowardReg,w);
1234                  fillProfile(hdNdEtadPhi_caloJet_Transverse900_,
1235 (caloJets->begin())->pt(),nTrk900_TransReg,w);
1236                  fillProfile(hdNdEtadPhi_caloJet_Away900_,
1237 (caloJets->begin())->pt(),nTrk900_AwayReg,w);
1238 
1239                  fillProfile(hpTSumdEtadPhi_caloJet_Toward900_,
1240 (caloJets->begin())->pt(),pTSum900_TowardReg,w);
1241                  fillProfile(hpTSumdEtadPhi_caloJet_Transverse900_,
1242 (caloJets->begin())->pt(),pTSum900_TransReg,w);
1243                  fillProfile(hpTSumdEtadPhi_caloJet_Away900_,
1244 (caloJets->begin())->pt(),pTSum900_AwayReg,w);
1245 
1246 
1247 }
1248 
1249 */
1250 void QcdUeDQM::fillHltBits(const Event &iEvent, const EventSetup &iSetup) {
1251   // Fill HLT trigger bits.
1252 
1253   Handle<TriggerResults> triggerResultsHLT;
1254   getProduct(hltUsedResName_, triggerResultsHLT, iEvent);
1255 
1256   /*  const unsigned int ntrigs(triggerResultsHLT.product()->size());
1257     if( ntrigs != 0 ) {
1258     const edm::TriggerNames & triggerNames =
1259   iEvent.triggerNames(*triggerResultsHLT);
1260     for (unsigned int j=0; j!=ntrigs; j++) {
1261   //    if(triggerResultsHLT->accept(j)){
1262   //    unsigned int hlt_prescale = hltConfig_.prescaleValue(iEvent, iSetup,
1263   triggerName);
1264       cout<<"trigger fired"<<triggerNames.triggerName(j)<<endl;
1265       for(unsigned int itrigName = 0; itrigName < hltTrgNames_.size();
1266   itrigName++ ) {
1267        unsigned int hlt_prescale = hltConfig.prescaleValue(iEvent, iSetup,
1268   hltTrgNames_[itrigName]);
1269 
1270        if(triggerNames.triggerIndex(hltTrgNames_[itrigName]) >= (unsigned
1271   int)ntrigs ) continue;
1272   //    if(
1273   triggerResultsHLT->accept(triggerNames.triggerIndex(hltTrgNames_[itrigName]))
1274   )cout<<hltTrgNames_[itrigName]<<endl;
1275 
1276    }
1277       }
1278      }
1279    // }
1280   */
1281   for (size_t i = 0; i < hltTrgBits_.size(); ++i) {
1282     if (hltTrgBits_.at(i) < 0)
1283       continue;  // ignore unknown trigger
1284     size_t tbit = hltTrgBits_.at(i);
1285     if (tbit < triggerResultsHLT->size()) {
1286       hltTrgDeci_[i] = triggerResultsHLT->accept(tbit);
1287     }
1288   }
1289   // fill correlation histogram
1290   for (size_t i = 0; i < hltTrgBits_.size(); ++i) {
1291     if (hltTrgDeci_.at(i))
1292       h2TrigCorr_->Fill(i, i);
1293     for (size_t j = i + 1; j < hltTrgBits_.size(); ++j) {
1294       if (hltTrgDeci_.at(i) && hltTrgDeci_.at(j))
1295         h2TrigCorr_->Fill(i, j);
1296       if (hltTrgDeci_.at(i) && !hltTrgDeci_.at(j))
1297         h2TrigCorr_->Fill(j, i);
1298     }
1299   }
1300 }