File indexing completed on 2024-04-06 12:29:32
0001
0002
0003
0004
0005
0006 #include "CalibFormats/CaloObjects/interface/CaloSamples.h"
0007 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
0008 #include "FWCore/Framework/interface/Event.h"
0009 #include "FWCore/Framework/interface/EventSetup.h"
0010 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0013 #include "FWCore/Utilities/interface/InputTag.h"
0014 #include "SimCalorimetry/CaloSimAlgos/interface/CaloHitAnalyzer.h"
0015 #include "SimCalorimetry/CaloSimAlgos/interface/CaloValidationStatistics.h"
0016 #include "SimCalorimetry/HcalSimAlgos/interface/HcalHitFilter.h"
0017 #include "SimCalorimetry/HcalSimAlgos/interface/HcalSimParameterMap.h"
0018 #include "SimCalorimetry/HcalSimAlgos/interface/ZDCHitFilter.h"
0019
0020 #include <iostream>
0021 #include <string>
0022
0023 class HcalDigiStatistics {
0024 public:
0025 HcalDigiStatistics(std::string name,
0026 int maxBin,
0027 float amplitudeThreshold,
0028 float expectedPedestal,
0029 float binPrevToBinMax,
0030 float binNextToBinMax,
0031 CaloHitAnalyzer &litudeAnalyzer)
0032 : maxBin_(maxBin),
0033 amplitudeThreshold_(amplitudeThreshold),
0034 pedestal_(name + " pedestal", expectedPedestal, 0.),
0035 binPrevToBinMax_(name + " binPrevToBinMax", binPrevToBinMax, 0.),
0036 binNextToBinMax_(name + " binNextToBinMax", binNextToBinMax, 0.),
0037 amplitudeAnalyzer_(amplitudeAnalyzer) {}
0038
0039 template <class Digi>
0040 void analyze(const Digi &digi);
0041
0042 private:
0043 int maxBin_;
0044 float amplitudeThreshold_;
0045 CaloValidationStatistics pedestal_;
0046 CaloValidationStatistics binPrevToBinMax_;
0047 CaloValidationStatistics binNextToBinMax_;
0048 CaloHitAnalyzer &litudeAnalyzer_;
0049 };
0050
0051 template <class Digi>
0052 void HcalDigiStatistics::analyze(const Digi &digi) {
0053 pedestal_.addEntry(digi[0].adc());
0054 pedestal_.addEntry(digi[1].adc());
0055
0056 double pedestal_fC = 0.5 * (digi[0].nominal_fC() + digi[1].nominal_fC());
0057
0058 double maxAmplitude = digi[maxBin_].nominal_fC() - pedestal_fC;
0059
0060 if (maxAmplitude > amplitudeThreshold_) {
0061 double binPrevToBinMax = (digi[maxBin_ - 1].nominal_fC() - pedestal_fC) / maxAmplitude;
0062 binPrevToBinMax_.addEntry(binPrevToBinMax);
0063
0064 double binNextToBinMax = (digi[maxBin_ + 1].nominal_fC() - pedestal_fC) / maxAmplitude;
0065 binNextToBinMax_.addEntry(binNextToBinMax);
0066
0067 double amplitude = digi[maxBin_].nominal_fC() + digi[maxBin_ + 1].nominal_fC() - 2 * pedestal_fC;
0068
0069 amplitudeAnalyzer_.analyze(digi.id().rawId(), amplitude);
0070 }
0071 }
0072
0073 class HcalDigiAnalyzer : public edm::one::EDAnalyzer<> {
0074 public:
0075 explicit HcalDigiAnalyzer(edm::ParameterSet const &conf);
0076 void analyze(edm::Event const &e, edm::EventSetup const &c) override;
0077
0078 private:
0079 std::string hitReadoutName_;
0080 HcalSimParameterMap simParameterMap_;
0081 HBHEHitFilter hbheFilter_;
0082 HOHitFilter hoFilter_;
0083 HFHitFilter hfFilter_;
0084 ZDCHitFilter zdcFilter_;
0085 CaloHitAnalyzer hbheHitAnalyzer_;
0086 CaloHitAnalyzer hoHitAnalyzer_;
0087 CaloHitAnalyzer hfHitAnalyzer_;
0088 CaloHitAnalyzer zdcHitAnalyzer_;
0089 HcalDigiStatistics hbheDigiStatistics_;
0090 HcalDigiStatistics hoDigiStatistics_;
0091 HcalDigiStatistics hfDigiStatistics_;
0092 HcalDigiStatistics zdcDigiStatistics_;
0093
0094 const edm::InputTag hbheDigiCollectionTag_;
0095 const edm::InputTag hoDigiCollectionTag_;
0096 const edm::InputTag hfDigiCollectionTag_;
0097 const edm::EDGetTokenT<HBHEDigiCollection> hbheDigiCollectionToken_;
0098 const edm::EDGetTokenT<HODigiCollection> hoDigiCollectionToken_;
0099 const edm::EDGetTokenT<HFDigiCollection> hfDigiCollectionToken_;
0100 const edm::EDGetTokenT<CrossingFrame<PCaloHit>> cfToken_;
0101 const edm::EDGetTokenT<CrossingFrame<PCaloHit>> zdccfToken_;
0102 };
0103
0104 HcalDigiAnalyzer::HcalDigiAnalyzer(edm::ParameterSet const &conf)
0105 : hitReadoutName_("HcalHits"),
0106 simParameterMap_(),
0107 hbheFilter_(),
0108 hoFilter_(),
0109 hfFilter_(),
0110 hbheHitAnalyzer_("HBHEDigi", 1., &simParameterMap_, &hbheFilter_),
0111 hoHitAnalyzer_("HODigi", 1., &simParameterMap_, &hoFilter_),
0112 hfHitAnalyzer_("HFDigi", 1., &simParameterMap_, &hfFilter_),
0113 zdcHitAnalyzer_("ZDCDigi", 1., &simParameterMap_, &zdcFilter_),
0114 hbheDigiStatistics_("HBHEDigi", 4, 10., 6., 0.1, 0.5, hbheHitAnalyzer_),
0115 hoDigiStatistics_("HODigi", 4, 10., 6., 0.1, 0.5, hoHitAnalyzer_),
0116 hfDigiStatistics_("HFDigi", 3, 10., 6., 0.1, 0.5, hfHitAnalyzer_),
0117 zdcDigiStatistics_("ZDCDigi", 3, 10., 6., 0.1, 0.5, zdcHitAnalyzer_),
0118 hbheDigiCollectionTag_(conf.getParameter<edm::InputTag>("hbheDigiCollectionTag")),
0119 hoDigiCollectionTag_(conf.getParameter<edm::InputTag>("hoDigiCollectionTag")),
0120 hfDigiCollectionTag_(conf.getParameter<edm::InputTag>("hfDigiCollectionTag")),
0121 hbheDigiCollectionToken_(consumes(hbheDigiCollectionTag_)),
0122 hoDigiCollectionToken_(consumes(hoDigiCollectionTag_)),
0123 hfDigiCollectionToken_(consumes(hfDigiCollectionTag_)),
0124 cfToken_(consumes(edm::InputTag("mix", "HcalHits"))),
0125 zdccfToken_(consumes(edm::InputTag("mix", "ZDCHits"))) {}
0126
0127 namespace HcalDigiAnalyzerImpl {
0128 template <class Collection>
0129 void analyze(edm::Event const &e, HcalDigiStatistics &statistics, const edm::EDGetTokenT<Collection> &token) {
0130 const edm::Handle<Collection> &digis = e.getHandle(token);
0131 for (unsigned i = 0; i < digis->size(); ++i) {
0132 edm::LogVerbatim("HcalSim") << (*digis)[i];
0133 statistics.analyze((*digis)[i]);
0134 }
0135 }
0136 }
0137
0138 void HcalDigiAnalyzer::analyze(edm::Event const &e, edm::EventSetup const &c) {
0139
0140 const edm::Handle<CrossingFrame<PCaloHit>> &cf = e.getHandle(cfToken_);
0141
0142
0143
0144 std::unique_ptr<MixCollection<PCaloHit>> hits(new MixCollection<PCaloHit>(cf.product()));
0145
0146 hbheHitAnalyzer_.fillHits(*hits);
0147 hoHitAnalyzer_.fillHits(*hits);
0148 hfHitAnalyzer_.fillHits(*hits);
0149
0150 HcalDigiAnalyzerImpl::analyze<HBHEDigiCollection>(e, hbheDigiStatistics_, hbheDigiCollectionToken_);
0151 HcalDigiAnalyzerImpl::analyze<HODigiCollection>(e, hoDigiStatistics_, hoDigiCollectionToken_);
0152 HcalDigiAnalyzerImpl::analyze<HFDigiCollection>(e, hfDigiStatistics_, hfDigiCollectionToken_);
0153
0154 }
0155
0156 #include "FWCore/Framework/interface/MakerMacros.h"
0157
0158 DEFINE_FWK_MODULE(HcalDigiAnalyzer);