Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:29:21

0001 /**  Castor digis
0002  Author: Panos Katsas
0003 */
0004 
0005 #include "CalibFormats/CaloObjects/interface/CaloSamples.h"
0006 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
0007 #include "FWCore/Framework/interface/Event.h"
0008 #include "FWCore/Framework/interface/EventSetup.h"
0009 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0012 #include "FWCore/Utilities/interface/InputTag.h"
0013 #include "SimCalorimetry/CaloSimAlgos/interface/CaloHitAnalyzer.h"
0014 #include "SimCalorimetry/CaloSimAlgos/interface/CaloValidationStatistics.h"
0015 #include "SimCalorimetry/CastorSim/interface/CastorHitFilter.h"
0016 #include "SimCalorimetry/CastorSim/interface/CastorSimParameterMap.h"
0017 
0018 #include <iostream>
0019 #include <string>
0020 
0021 class CastorDigiStatistics {
0022 public:
0023   CastorDigiStatistics(std::string name,
0024                        int maxBin,
0025                        float amplitudeThreshold,
0026                        float expectedPedestal,
0027                        float binPrevToBinMax,
0028                        float binNextToBinMax,
0029                        CaloHitAnalyzer &amplitudeAnalyzer)
0030       : maxBin_(maxBin),
0031         amplitudeThreshold_(amplitudeThreshold),
0032         pedestal_(name + " pedestal", expectedPedestal, 0.),
0033         binPrevToBinMax_(name + " binPrevToBinMax", binPrevToBinMax, 0.),
0034         binNextToBinMax_(name + " binNextToBinMax", binNextToBinMax, 0.),
0035         amplitudeAnalyzer_(amplitudeAnalyzer) {}
0036 
0037   template <class Digi>
0038   void analyze(const Digi &digi);
0039 
0040 private:
0041   int maxBin_;
0042   float amplitudeThreshold_;
0043   CaloValidationStatistics pedestal_;
0044   CaloValidationStatistics binPrevToBinMax_;
0045   CaloValidationStatistics binNextToBinMax_;
0046   CaloHitAnalyzer &amplitudeAnalyzer_;
0047 };
0048 
0049 template <class Digi>
0050 void CastorDigiStatistics::analyze(const Digi &digi) {
0051   pedestal_.addEntry(digi[0].adc());
0052   pedestal_.addEntry(digi[1].adc());
0053 
0054   double pedestal_fC = 0.5 * (digi[0].nominal_fC() + digi[1].nominal_fC());
0055 
0056   double maxAmplitude = digi[maxBin_].nominal_fC() - pedestal_fC;
0057 
0058   if (maxAmplitude > amplitudeThreshold_) {
0059     double binPrevToBinMax = (digi[maxBin_ - 1].nominal_fC() - pedestal_fC) / maxAmplitude;
0060     binPrevToBinMax_.addEntry(binPrevToBinMax);
0061 
0062     double binNextToBinMax = (digi[maxBin_ + 1].nominal_fC() - pedestal_fC) / maxAmplitude;
0063     binNextToBinMax_.addEntry(binNextToBinMax);
0064 
0065     double amplitude = digi[maxBin_].nominal_fC() + digi[maxBin_ + 1].nominal_fC() - 2 * pedestal_fC;
0066 
0067     amplitudeAnalyzer_.analyze(digi.id().rawId(), amplitude);
0068   }
0069 }
0070 
0071 class CastorDigiAnalyzer : public edm::one::EDAnalyzer<> {
0072 public:
0073   explicit CastorDigiAnalyzer(edm::ParameterSet const &conf);
0074   void analyze(edm::Event const &e, edm::EventSetup const &c) override;
0075 
0076 private:
0077   std::string hitReadoutName_;
0078   CastorSimParameterMap simParameterMap_;
0079   CastorHitFilter castorFilter_;
0080   CaloHitAnalyzer castorHitAnalyzer_;
0081   CastorDigiStatistics castorDigiStatistics_;
0082   const edm::EDGetTokenT<CastorDigiCollection> castordigiToken_;
0083   const edm::EDGetTokenT<CrossingFrame<PCaloHit>> castorcfToken_;
0084 };
0085 
0086 CastorDigiAnalyzer::CastorDigiAnalyzer(edm::ParameterSet const &conf)
0087     : hitReadoutName_("CastorHits"),
0088       simParameterMap_(),
0089       castorHitAnalyzer_("CASTORDigi", 1., &simParameterMap_, &castorFilter_),
0090       castorDigiStatistics_("CASTORDigi", 3, 10., 6., 0.1, 0.5, castorHitAnalyzer_),
0091       castordigiToken_(consumes<CastorDigiCollection>(conf.getParameter<edm::InputTag>("castorDigiCollectionTag"))),
0092       castorcfToken_(consumes<CrossingFrame<PCaloHit>>(edm::InputTag("mix", "g4SimHitsCastorFI"))) {}
0093 
0094 namespace CastorDigiAnalyzerImpl {
0095   template <class Collection>
0096   void analyze(edm::Event const &e, CastorDigiStatistics &statistics, const edm::EDGetTokenT<Collection> &token) {
0097     const edm::Handle<Collection> &digis = e.getHandle(token);
0098     if (!digis.isValid()) {
0099       edm::LogError("CastorDigiAnalyzer") << "Could not find Castor Digi Container ";
0100     } else {
0101       for (unsigned i = 0; i < digis->size(); ++i) {
0102         statistics.analyze((*digis)[i]);
0103       }
0104     }
0105   }
0106 }  // namespace CastorDigiAnalyzerImpl
0107 
0108 void CastorDigiAnalyzer::analyze(edm::Event const &e, edm::EventSetup const &c) {
0109   //  edm::Handle<edm::PCaloHitContainer> hits;
0110   const edm::Handle<CrossingFrame<PCaloHit>> &castorcf = e.getHandle(castorcfToken_);
0111 
0112   // access to SimHits
0113   std::unique_ptr<MixCollection<PCaloHit>> hits(new MixCollection<PCaloHit>(castorcf.product()));
0114   //  if (hits.isValid()) {
0115   castorHitAnalyzer_.fillHits(*hits);
0116   CastorDigiAnalyzerImpl::analyze<CastorDigiCollection>(e, castorDigiStatistics_, castordigiToken_);
0117 }
0118 
0119 #include "FWCore/Framework/interface/MakerMacros.h"
0120 
0121 DEFINE_FWK_MODULE(CastorDigiAnalyzer);