Cache

Pixel

QcdLowPtDQM

Tracklet

Vertex

Macros

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394

#ifndef QcdLowPtDQM_H
#define QcdLowPtDQM_H

#include "DQMServices/Core/interface/DQMOneEDAnalyzer.h"
#include "DataFormats/Common/interface/Handle.h"
#include "DataFormats/GeometryVector/interface/VectorUtil.h"
#include "DataFormats/GeometryVector/interface/GlobalPoint.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
#include "DQMServices/Core/interface/DQMStore.h"
#include <TMath.h>
#include <vector>

class TrackerGeometry;
class TrackerDigiGeometryRecord;
class TrackerTopology;
class TrackerTopologyRcd;
class TH1F;
class TH2F;
class TH3F;

namespace qlpd {
  struct Cache {};
}  // namespace qlpd

class QcdLowPtDQM : public DQMOneEDAnalyzer<edm::LuminosityBlockCache<qlpd::Cache>> {
public:
  class Pixel {
  public:
    Pixel(double x = 0,
          double y = 0,
          double z = 0,
          double eta = 0,
          double phi = 0,
          double adc = 0,
          double sx = 0,
          double sy = 0)
        : x_(x),
          y_(y),
          z_(z),
          rho_(TMath::Sqrt(x_ * x_ + y_ * y_)),
          eta_(eta),
          phi_(phi),
          adc_(adc),
          sizex_(sx),
          sizey_(sy) {}
    Pixel(const GlobalPoint &p, double adc = 0, double sx = 0, double sy = 0)
        : x_(p.x()),
          y_(p.y()),
          z_(p.z()),
          rho_(TMath::Sqrt(x_ * x_ + y_ * y_)),
          eta_(p.eta()),
          phi_(p.phi()),
          adc_(adc),
          sizex_(sx),
          sizey_(sy) {}
    double adc() const { return adc_; }
    double eta() const { return eta_; }
    double rho() const { return rho_; }
    double phi() const { return phi_; }
    double sizex() const { return sizex_; }
    double sizey() const { return sizey_; }
    double x() const { return x_; }
    double y() const { return y_; }
    double z() const { return z_; }

  protected:
    double x_, y_, z_, rho_, eta_, phi_;
    double adc_, sizex_, sizey_;
  };
  class Tracklet {
  public:
    Tracklet() : i1_(-1), i2_(-2), deta_(0), dphi_(0) {}
    Tracklet(const Pixel &p1, const Pixel &p2)
        : p1_(p1), p2_(p2), i1_(-1), i2_(-1), deta_(p1.eta() - p2.eta()), dphi_(Geom::deltaPhi(p1.phi(), p2.phi())) {}
    double deta() const { return deta_; }
    double dphi() const { return dphi_; }
    int i1() const { return i1_; }
    int i2() const { return i2_; }
    double eta() const { return p1_.eta(); }
    const Pixel &p1() const { return p1_; }
    const Pixel &p2() const { return p2_; }
    void seti1(int i1) { i1_ = i1; }
    void seti2(int i2) { i2_ = i2; }

  protected:
    Pixel p1_, p2_;
    int i1_, i2_;
    double deta_, dphi_;
  };
  class Vertex {
  public:
    Vertex(double x = 0, double y = 0, double z = 0, double xs = 0, double ys = 0, double zs = 0, int n = 0)
        : x_(x), y_(y), z_(z), xs_(xs), ys_(ys), zs_(zs), n_(n) {}
    int n() const { return n_; }
    double x() const { return x_; }
    double y() const { return y_; }
    double z() const { return z_; }
    double xs() const { return xs_; }
    double ys() const { return ys_; }
    double zs() const { return zs_; }
    void set(int n, double z, double zs) {
      n_ = n;
      z_ = z;
      zs_ = zs;
    }
    void set(int n, double x, double y, double z, double xs, double ys, double zs) {
      n_ = n;
      x_ = x;
      xs_ = xs;
      y_ = y;
      ys_ = ys;
      z_ = z;
      zs_ = zs;
    }

  protected:
    double x_, y_, z_, xs_, ys_, zs_;
    int n_;
  };

  QcdLowPtDQM(const edm::ParameterSet &parameters);
  ~QcdLowPtDQM() override;
  void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override;
  void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
  void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override;
  std::shared_ptr<qlpd::Cache> globalBeginLuminosityBlock(const edm::LuminosityBlock &,
                                                          const edm::EventSetup &) const override;
  void globalEndLuminosityBlock(const edm::LuminosityBlock &l, const edm::EventSetup &iSetup) override;
  void dqmEndRun(const edm::Run &r, const edm::EventSetup &iSetup) override;

private:
  void book1D(DQMStore::IBooker &,
              std::vector<MonitorElement *> &mes,
              const std::string &name,
              const std::string &title,
              int nx,
              double x1,
              double x2,
              bool sumw2 = true,
              bool sbox = true);
  void book2D(DQMStore::IBooker &,
              std::vector<MonitorElement *> &mes,
              const std::string &name,
              const std::string &title,
              int nx,
              double x1,
              double x2,
              int ny,
              double y1,
              double y2,
              bool sumw2 = true,
              bool sbox = true);
  void create1D(std::vector<TH1F *> &mes,
                const std::string &name,
                const std::string &title,
                int nx,
                double x1,
                double x2,
                bool sumw2 = true,
                bool sbox = true);
  void create2D(std::vector<TH2F *> &mes,
                const std::string &name,
                const std::string &title,
                int nx,
                double x1,
                double x2,
                int ny,
                double y1,
                double y2,
                bool sumw2 = true,
                bool sbox = true);
  void fill1D(std::vector<TH1F *> &hs, double val, double w = 1.);
  void fill1D(std::vector<MonitorElement *> &mes, double val, double w = 1.);
  void fill2D(std::vector<TH2F *> &hs, double valx, double valy, double w = 1.);
  void fill2D(std::vector<MonitorElement *> &mes, double valx, double valy, double w = 1.);
  void fill3D(std::vector<TH3F *> &hs, int gbin, double w = 1.);
  void filldNdeta(const TH3F *AlphaTracklets,
                  const std::vector<TH3F *> &NsigTracklets,
                  const std::vector<TH3F *> &NbkgTracklets,
                  const std::vector<TH1F *> &NEvsPerEta,
                  std::vector<MonitorElement *> &hdNdEtaRawTrkl,
                  std::vector<MonitorElement *> &hdNdEtaSubTrkl,
                  std::vector<MonitorElement *> &hdNdEtaTrklets);
  void fillHltBits(const edm::Event &iEvent);
  void fillPixels(const edm::Event &iEvent, const edm::EventSetup &iSetup);
  void fillPixelClusterInfos(const edm::Event &iEvent, int which = 12);
  void fillPixelClusterInfos(const double vz,
                             const std::vector<Pixel> &pix,
                             std::vector<MonitorElement *> &hClusterYSize,
                             std::vector<MonitorElement *> &hClusterADC);
  void fillTracklets(const edm::Event &iEvent, int which = 12);
  void fillTracklets(std::vector<Tracklet> &tracklets,
                     const std::vector<Pixel> &pix1,
                     const std::vector<Pixel> &pix2,
                     const Vertex &trackletV);
  void fillTracklets(const std::vector<Tracklet> &tracklets,
                     const std::vector<Pixel> &pixels,
                     const Vertex &trackletV,
                     const TH3F *AlphaTracklets,
                     std::vector<TH3F *> &NsigTracklets,
                     std::vector<TH3F *> &NbkgTracklets,
                     std::vector<TH1F *> &eventpereta,
                     std::vector<MonitorElement *> &detaphi,
                     std::vector<MonitorElement *> &deta,
                     std::vector<MonitorElement *> &dphi,
                     std::vector<MonitorElement *> &etavsvtx);
  template <typename TYPE>
  void getProduct(const std::string name, edm::Handle<TYPE> &prod, const edm::Event &event) const;
  template <typename TYPE>
  bool getProductSafe(const std::string name, edm::Handle<TYPE> &prod, const edm::Event &event) const;
  void print(int level, const char *msg);
  void print(int level, const std::string &msg) { print(level, msg.c_str()); }
  void reallyPrint(int level, const char *msg);
  void trackletVertexUnbinned(const edm::Event &iEvent, int which = 12);
  void trackletVertexUnbinned(std::vector<Pixel> &pix1, std::vector<Pixel> &pix2, Vertex &vtx);
  double vertexZFromClusters(const std::vector<Pixel> &pix) const;
  void yieldAlphaHistogram(int which = 12);

  edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> tkGeomToken_;
  edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> tTopoToken_;
  std::string hltResName_;                    // HLT trigger results name
  std::vector<std::string> hltProcNames_;     // HLT process name(s)
  std::vector<std::string> hltTrgNames_;      // HLT trigger name(s)
  std::string pixelName_;                     // pixel reconstructed hits name
  std::string clusterVtxName_;                // cluster vertex name
  double ZVCut_;                              // Z vertex cut for selected events
  double ZVEtaRegion_;                        // Z vertex eta region
  double ZVVtxRegion_;                        // Z vertex vtx region
  double dPhiVc_;                             // dPhi vertex cut for tracklet based vertex
  double dZVc_;                               // dZ vertex cut for tracklet based vertex
  double sigEtaCut_;                          // signal tracklet eta cut
  double sigPhiCut_;                          // signal tracklet phi cut
  double bkgEtaCut_;                          // bkg tracklet eta cut
  double bkgPhiCut_;                          // bgk tracklet phi cut
  int verbose_;                               // verbosity (0=debug,1=warn,2=error,3=throw)
  int pixLayers_;                             // 12 for 12, 13 for 12 and 13, 23 for all
  int clusLayers_;                            // 12 for 12, 13 for 12 and 13, 23 for all
  bool useRecHitQ_;                           // if true use rec hit quality word
  bool usePixelQ_;                            // if true use pixel hit quality word
  std::vector<int> hltTrgBits_;               // HLT trigger bit(s)
  std::vector<bool> hltTrgDeci_;              // HLT trigger descision(s)
  std::vector<std::string> hltTrgUsedNames_;  // HLT used trigger name(s)
  std::string hltUsedResName_;                // used HLT trigger results name
  std::vector<Pixel> bpix1_;                  // barrel pixels layer 1
  std::vector<Pixel> bpix2_;                  // barrel pixels layer 2
  std::vector<Pixel> bpix3_;                  // barrel pixels layer 3
  std::vector<Tracklet> btracklets12_;        // barrel tracklets 12
  std::vector<Tracklet> btracklets13_;        // barrel tracklets 13
  std::vector<Tracklet> btracklets23_;        // barrel tracklets 23
  Vertex trackletV12_;                        // reconstructed tracklet vertex 12
  Vertex trackletV13_;                        // reconstructed tracklet vertex 13
  Vertex trackletV23_;                        // reconstructed tracklet vertex 23
  std::vector<TH3F *> NsigTracklets12_;       // number of signal tracklets 12
  std::vector<TH3F *> NbkgTracklets12_;       // number of background tracklets 12
  std::vector<TH1F *> hEvtCountsPerEta12_;    // event count per tracklet12
                                              // eta-vtx region
  TH3F *AlphaTracklets12_;                    // alpha correction for tracklets 12
  std::vector<TH3F *> NsigTracklets13_;       // number of signal tracklets 13
  std::vector<TH3F *> NbkgTracklets13_;       // number of background tracklets 13
  std::vector<TH1F *> hEvtCountsPerEta13_;    // event count per tracklet13
                                              // eta-vtx region
  TH3F *AlphaTracklets13_;                    // alpha correction for tracklets 13
  std::vector<TH3F *> NsigTracklets23_;       // number of signal tracklets 23
  std::vector<TH3F *> NbkgTracklets23_;       // number of background tracklets 23
  std::vector<TH1F *> hEvtCountsPerEta23_;    // event count per tracklet23
                                              // eta-vtx region
  TH3F *AlphaTracklets23_;                    // alpha correction for tracklets 23
  HLTConfigProvider hltConfig_;
  const TrackerGeometry *tgeo_;                      // tracker geometry
  MonitorElement *repSumMap_;                        // report summary map
  MonitorElement *repSummary_;                       // report summary
  MonitorElement *h2TrigCorr_;                       // trigger correlation plot
  std::vector<MonitorElement *> hNhitsL1_;           // number of hits on layer 1
  std::vector<MonitorElement *> hNhitsL2_;           // number of hits on layer 2
  std::vector<MonitorElement *> hNhitsL3_;           // number of hits on layer 3
  std::vector<MonitorElement *> hNhitsL1z_;          // number of hits on layer 1
                                                     // (zoomed)
  std::vector<MonitorElement *> hNhitsL2z_;          // number of hits on layer 2
                                                     // (zoomed)
  std::vector<MonitorElement *> hNhitsL3z_;          // number of hits on layer 3
                                                     // (zoomed)
  std::vector<MonitorElement *> hdNdEtaHitsL1_;      // dN/dEta of hits on layer 1
  std::vector<MonitorElement *> hdNdEtaHitsL2_;      // dN/dEta of hits on layer 2
  std::vector<MonitorElement *> hdNdEtaHitsL3_;      // dN/dEta of hits on layer 3
  std::vector<MonitorElement *> hdNdPhiHitsL1_;      // dN/dPhi of hits on layer 1
  std::vector<MonitorElement *> hdNdPhiHitsL2_;      // dN/dPhi of hits on layer 2
  std::vector<MonitorElement *> hdNdPhiHitsL3_;      // dN/dPhi of hits on layer 3
  std::vector<MonitorElement *> hTrkVtxZ12_;         // tracklet z vertex 12 histograms
  std::vector<MonitorElement *> hTrkVtxZ13_;         // tracklet z vertex 13 histograms
  std::vector<MonitorElement *> hTrkVtxZ23_;         // tracklet z vertex 23 histograms
  std::vector<MonitorElement *> hRawTrkEtaVtxZ12_;   // tracklet eta vs z vertex
                                                     // 12 histograms
  std::vector<MonitorElement *> hRawTrkEtaVtxZ13_;   // tracklet eta vs z vertex
                                                     // 13 histograms
  std::vector<MonitorElement *> hRawTrkEtaVtxZ23_;   // tracklet eta vs z vertex
                                                     // 23 histograms
  std::vector<MonitorElement *> hTrkRawDetaDphi12_;  // tracklet12 Deta/Dphi
                                                     // distribution
  std::vector<MonitorElement *> hTrkRawDeta12_;      // tracklet12 Deta distribution
  std::vector<MonitorElement *> hTrkRawDphi12_;      // tracklet12 Dphi distribution
  std::vector<MonitorElement *> hTrkRawDetaDphi13_;  // tracklet13 Deta/Dphi
                                                     // distribution
  std::vector<MonitorElement *> hTrkRawDeta13_;      // tracklet13 Deta distribution
  std::vector<MonitorElement *> hTrkRawDphi13_;      // tracklet13 Dphi distribution
  std::vector<MonitorElement *> hTrkRawDetaDphi23_;  // tracklet23 Deta/Dphi
                                                     // distribution
  std::vector<MonitorElement *> hTrkRawDeta23_;      // tracklet23 Deta distribution
  std::vector<MonitorElement *> hTrkRawDphi23_;      // tracklet23 Dphi distribution
  std::vector<MonitorElement *> hdNdEtaRawTrkl12_;   // dN/dEta from raw
                                                     // tracklets 12
  std::vector<MonitorElement *> hdNdEtaSubTrkl12_;   // dN/dEta from beta
                                                     // tracklets 12
  std::vector<MonitorElement *> hdNdEtaTrklets12_;   // dN/dEta corrected by
                                                     // alpha 12
  std::vector<MonitorElement *> hdNdEtaRawTrkl13_;   // dN/dEta from raw
                                                     // tracklets 13
  std::vector<MonitorElement *> hdNdEtaSubTrkl13_;   // dN/dEta from beta
                                                     // tracklets 13
  std::vector<MonitorElement *> hdNdEtaTrklets13_;   // dN/dEta corrected by
                                                     // alpha 13
  std::vector<MonitorElement *> hdNdEtaRawTrkl23_;   // dN/dEta from raw
                                                     // tracklets 23
  std::vector<MonitorElement *> hdNdEtaSubTrkl23_;   // dN/dEta from beta
                                                     // tracklets 23
  std::vector<MonitorElement *> hdNdEtaTrklets23_;   // dN/dEta corrected by
                                                     // alpha 23
  std::vector<MonitorElement *> hClusterVertexZ_;    // cluster z vertex
                                                     // histograms
  std::vector<MonitorElement *> hClusterYSize1_;     // cluster y size histograms
                                                     // on layer 1
  std::vector<MonitorElement *> hClusterYSize2_;     // cluster y size histograms
                                                     // on layer 2
  std::vector<MonitorElement *> hClusterYSize3_;     // cluster y size histograms
                                                     // on layer 3
  std::vector<MonitorElement *> hClusterADC1_;       // cluster adc histograms on
                                                     // layer 1
  std::vector<MonitorElement *> hClusterADC2_;       // cluster adc histograms on
                                                     // layer 2
  std::vector<MonitorElement *> hClusterADC3_;       // cluster adc histograms on
                                                     // layer 3
};

//--------------------------------------------------------------------------------------------------
template <typename TYPE>
inline void QcdLowPtDQM::getProduct(const std::string name, edm::Handle<TYPE> &prod, const edm::Event &event) const {
  // Try to access data collection from EDM file. We check if we really get just
  // one
  // product with the given name. If not we throw an exception.

  event.getByLabel(edm::InputTag(name), prod);
  if (!prod.isValid())
    throw edm::Exception(edm::errors::Configuration, "QcdLowPtDQM::GetProduct()\n")
        << "Collection with label " << name << " is not valid" << std::endl;
}

//--------------------------------------------------------------------------------------------------
template <typename TYPE>
inline bool QcdLowPtDQM::getProductSafe(const std::string name,
                                        edm::Handle<TYPE> &prod,
                                        const edm::Event &event) const {
  // Try to safely access data collection from EDM file. We check if we really
  // get just one
  // product with the given name. If not, we return false.

  if (name.empty())
    return false;

  try {
    event.getByLabel(edm::InputTag(name), prod);
    if (!prod.isValid())
      return false;
  } catch (...) {
    return false;
  }
  return true;
}

//--------------------------------------------------------------------------------------------------
inline void QcdLowPtDQM::print(int level, const char *msg) {
  // Print out message if it

  if (level >= verbose_)
    reallyPrint(level, msg);
}
#endif

/* Local Variables: */
/* show-trailing-whitespace: t */
/* truncate-lines: t */
/* End: */