PPSTimingCalibrationPCLWorker

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
/****************************************************************************
 *
 * This is a part of PPS offline software.
 * Authors:
 *   Edoardo Bossini
 *   Piotr Maciej Cwiklicki
 *   Laurent Forthomme
 *
 ****************************************************************************/

#include "DQMServices/Core/interface/DQMGlobalEDAnalyzer.h"
#include "DQMServices/Core/interface/DQMStore.h"

#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"

#include "Geometry/VeryForwardGeometryBuilder/interface/CTPPSGeometry.h"
#include "Geometry/Records/interface/VeryForwardRealGeometryRecord.h"

#include "DataFormats/Common/interface/DetSetVector.h"
#include "DataFormats/CTPPSDetId/interface/CTPPSDiamondDetId.h"
#include "DataFormats/CTPPSReco/interface/CTPPSDiamondRecHit.h"

#include "CalibPPS/TimingCalibration/interface/TimingCalibrationStruct.h"

//------------------------------------------------------------------------------

class PPSTimingCalibrationPCLWorker : public DQMGlobalEDAnalyzer<TimingCalibrationHistograms> {
public:
  explicit PPSTimingCalibrationPCLWorker(const edm::ParameterSet&);

  void dqmAnalyze(const edm::Event&, const edm::EventSetup&, const TimingCalibrationHistograms&) const override;

  static void fillDescriptions(edm::ConfigurationDescriptions&);

private:
  void bookHistograms(DQMStore::IBooker&,
                      const edm::Run&,
                      const edm::EventSetup&,
                      TimingCalibrationHistograms&) const override;

  template <typename T>
  bool searchForProduct(edm::Event const& iEvent,
                        const std::vector<edm::EDGetTokenT<T>>& tokens,
                        const std::vector<edm::InputTag>& tags,
                        edm::Handle<T>& handle) const;

  const std::vector<edm::InputTag> RecHitTags_;
  std::vector<edm::EDGetTokenT<edm::DetSetVector<CTPPSDiamondRecHit>>> diamondRecHitTokens_;
  const edm::ESGetToken<CTPPSGeometry, VeryForwardRealGeometryRecord> geomEsToken_;

  const std::string dqmDir_;
};

//------------------------------------------------------------------------------

PPSTimingCalibrationPCLWorker::PPSTimingCalibrationPCLWorker(const edm::ParameterSet& iConfig)
    : RecHitTags_(iConfig.getParameter<std::vector<edm::InputTag>>("diamondRecHitTags")),
      geomEsToken_(esConsumes<edm::Transition::BeginRun>()),
      dqmDir_(iConfig.getParameter<std::string>("dqmDir")) {
  for (auto& tag : RecHitTags_)
    diamondRecHitTokens_.push_back(consumes<edm::DetSetVector<CTPPSDiamondRecHit>>(tag));
}

//------------------------------------------------------------------------------

void PPSTimingCalibrationPCLWorker::bookHistograms(DQMStore::IBooker& iBooker,
                                                   const edm::Run& iRun,
                                                   const edm::EventSetup& iSetup,
                                                   TimingCalibrationHistograms& iHists) const {
  iBooker.cd();
  iBooker.setCurrentFolder(dqmDir_);
  std::string ch_name;

  const auto& geom = iSetup.getData(geomEsToken_);
  for (auto it = geom.beginSensor(); it != geom.endSensor(); ++it) {
    if (!CTPPSDiamondDetId::check(it->first))
      continue;
    const CTPPSDiamondDetId detid(it->first);

    detid.channelName(ch_name);
    iHists.leadingTime[detid.rawId()] = iBooker.book1D("t_" + ch_name, ch_name + ";t (ns);Entries", 1200, -60., 60.);
    iHists.toT[detid.rawId()] = iBooker.book1D("tot_" + ch_name, ch_name + ";ToT (ns);Entries", 100, -20., 20.);
    iHists.leadingTimeVsToT[detid.rawId()] =
        iBooker.book2D("tvstot_" + ch_name, ch_name + ";ToT (ns);t (ns)", 240, 0., 60., 450, -20., 25.);
  }
}

//------------------------------------------------------------------------------

void PPSTimingCalibrationPCLWorker::dqmAnalyze(const edm::Event& iEvent,
                                               const edm::EventSetup& iSetup,
                                               const TimingCalibrationHistograms& iHists) const {
  edm::Handle<edm::DetSetVector<CTPPSDiamondRecHit>> dsv_rechits;
  // then extract the rechits information for later processing
  searchForProduct(iEvent, diamondRecHitTokens_, RecHitTags_, dsv_rechits);

  // ensure timing detectors rechits are found in the event content
  if (dsv_rechits->empty()) {
    edm::LogWarning("PPSTimingCalibrationPCLWorker:dqmAnalyze") << "No rechits retrieved from the event content.";
    return;
  }
  for (const auto& ds_rechits : *dsv_rechits) {
    const CTPPSDiamondDetId detid(ds_rechits.detId());
    if (iHists.leadingTimeVsToT.count(detid.rawId()) == 0) {
      edm::LogWarning("PPSTimingCalibrationPCLWorker:dqmAnalyze")
          << "Pad with detId=" << detid << " is not set to be monitored.";
      continue;
    }
    for (const auto& rechit : ds_rechits) {
      // skip invalid rechits
      if (rechit.time() == 0. || rechit.toT() < 0.)
        continue;
      iHists.leadingTime.at(detid.rawId())->Fill(rechit.time());
      iHists.toT.at(detid.rawId())->Fill(rechit.toT());
      iHists.leadingTimeVsToT.at(detid.rawId())->Fill(rechit.toT(), rechit.time());
    }
  }
}

//------------------------------------------------------------------------------

void PPSTimingCalibrationPCLWorker::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
  edm::ParameterSetDescription desc;
  desc.add<std::vector<edm::InputTag>>("diamondRecHitTags", {edm::InputTag("ctppsDiamondRecHits")})
      ->setComment("input tag for the PPS diamond detectors rechits");
  desc.add<std::string>("dqmDir", "AlCaReco/PPSTimingCalibrationPCL")
      ->setComment("output path for the various DQM plots");

  descriptions.addWithDefaultLabel(desc);
}

template <typename T>
bool PPSTimingCalibrationPCLWorker::searchForProduct(edm::Event const& iEvent,
                                                     const std::vector<edm::EDGetTokenT<T>>& tokens,
                                                     const std::vector<edm::InputTag>& tags,
                                                     edm::Handle<T>& handle) const {
  bool foundProduct = false;
  for (unsigned int i = 0; i < tokens.size(); i++)
    if (auto h = iEvent.getHandle(tokens[i])) {
      handle = h;
      foundProduct = true;
      edm::LogInfo("searchForProduct") << "Found a product with " << tags[i];
      break;
    }

  if (!foundProduct)
    throw edm::Exception(edm::errors::ProductNotFound) << "Could not find a product with any of the selected labels.";

  return foundProduct;
}

DEFINE_FWK_MODULE(PPSTimingCalibrationPCLWorker);