CMS 3D CMS Logo

CastorDigiAnalyzer.cc
Go to the documentation of this file.
1 
17 
18 #include <iostream>
19 #include <string>
20 
22 public:
24  int maxBin,
25  float amplitudeThreshold,
26  float expectedPedestal,
27  float binPrevToBinMax,
28  float binNextToBinMax,
29  CaloHitAnalyzer &amplitudeAnalyzer)
30  : maxBin_(maxBin),
32  pedestal_(name + " pedestal", expectedPedestal, 0.),
33  binPrevToBinMax_(name + " binPrevToBinMax", binPrevToBinMax, 0.),
34  binNextToBinMax_(name + " binNextToBinMax", binNextToBinMax, 0.),
35  amplitudeAnalyzer_(amplitudeAnalyzer) {}
36 
37  template <class Digi>
38  void analyze(const Digi &digi);
39 
40 private:
41  int maxBin_;
47 };
48 
49 template <class Digi>
51  pedestal_.addEntry(digi[0].adc());
52  pedestal_.addEntry(digi[1].adc());
53 
54  double pedestal_fC = 0.5 * (digi[0].nominal_fC() + digi[1].nominal_fC());
55 
56  double maxAmplitude = digi[maxBin_].nominal_fC() - pedestal_fC;
57 
58  if (maxAmplitude > amplitudeThreshold_) {
59  double binPrevToBinMax = (digi[maxBin_ - 1].nominal_fC() - pedestal_fC) / maxAmplitude;
60  binPrevToBinMax_.addEntry(binPrevToBinMax);
61 
62  double binNextToBinMax = (digi[maxBin_ + 1].nominal_fC() - pedestal_fC) / maxAmplitude;
63  binNextToBinMax_.addEntry(binNextToBinMax);
64 
65  double amplitude = digi[maxBin_].nominal_fC() + digi[maxBin_ + 1].nominal_fC() - 2 * pedestal_fC;
66 
67  amplitudeAnalyzer_.analyze(digi.id().rawId(), amplitude);
68  }
69 }
70 
72 public:
73  explicit CastorDigiAnalyzer(edm::ParameterSet const &conf);
74  void analyze(edm::Event const &e, edm::EventSetup const &c) override;
75 
76 private:
84 };
85 
87  : hitReadoutName_("CastorHits"),
88  simParameterMap_(),
89  castorHitAnalyzer_("CASTORDigi", 1., &simParameterMap_, &castorFilter_),
90  castorDigiStatistics_("CASTORDigi", 3, 10., 6., 0.1, 0.5, castorHitAnalyzer_),
91  castordigiToken_(consumes<CastorDigiCollection>(conf.getParameter<edm::InputTag>("castorDigiCollectionTag"))),
92  castorcfToken_(consumes<CrossingFrame<PCaloHit>>(edm::InputTag("mix", "g4SimHitsCastorFI"))) {}
93 
95  template <class Collection>
97  const edm::Handle<Collection> &digis = e.getHandle(token);
98  if (!digis.isValid()) {
99  edm::LogError("CastorDigiAnalyzer") << "Could not find Castor Digi Container ";
100  } else {
101  for (unsigned i = 0; i < digis->size(); ++i) {
102  statistics.analyze((*digis)[i]);
103  }
104  }
105  }
106 } // namespace CastorDigiAnalyzerImpl
107 
109  // edm::Handle<edm::PCaloHitContainer> hits;
110  const edm::Handle<CrossingFrame<PCaloHit>> &castorcf = e.getHandle(castorcfToken_);
111 
112  // access to SimHits
113  std::unique_ptr<MixCollection<PCaloHit>> hits(new MixCollection<PCaloHit>(castorcf.product()));
114  // if (hits.isValid()) {
116  CastorDigiAnalyzerImpl::analyze<CastorDigiCollection>(e, castorDigiStatistics_, castordigiToken_);
117 }
118 
120 
CastorHitFilter castorFilter_
void analyze(int detId, double recEnergy)
to be called for each RecHit
CaloHitAnalyzer castorHitAnalyzer_
CaloValidationStatistics binPrevToBinMax_
void analyze(edm::Event const &e, CastorDigiStatistics &statistics, const edm::EDGetTokenT< Collection > &token)
T const * product() const
Definition: Handle.h:70
CastorDigiStatistics castorDigiStatistics_
const edm::EDGetTokenT< CastorDigiCollection > castordigiToken_
constexpr unsigned int maxBin
Log< level::Error, false > LogError
CastorSimParameterMap simParameterMap_
void fillHits(MixCollection< PCaloHit > &hits)
should be called each event
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
CaloHitAnalyzer & amplitudeAnalyzer_
CaloValidationStatistics pedestal_
const edm::EDGetTokenT< CrossingFrame< PCaloHit > > castorcfToken_
void addEntry(float value, float weight=1.)
CastorDigiStatistics(std::string name, int maxBin, float amplitudeThreshold, float expectedPedestal, float binPrevToBinMax, float binNextToBinMax, CaloHitAnalyzer &amplitudeAnalyzer)
bool isValid() const
Definition: HandleBase.h:70
HLT enums.
void analyze(const Digi &digi)
CastorDigiAnalyzer(edm::ParameterSet const &conf)
uint16_t *__restrict__ uint16_t const *__restrict__ adc
CaloValidationStatistics binNextToBinMax_
void analyze(edm::Event const &e, edm::EventSetup const &c) override