File indexing completed on 2024-04-06 12:27:39
0001
0002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0003 #include "RecoPPS/Local/interface/CTPPSPixelClusterProducer.h"
0004
0005 CTPPSPixelClusterProducer::CTPPSPixelClusterProducer(const edm::ParameterSet &conf)
0006 : tokenCTPPSPixelDigi_(consumes<edm::DetSetVector<CTPPSPixelDigi> >(conf.getParameter<edm::InputTag>("tag"))),
0007 tokenCTPPSPixelAnalysisMask_(esConsumes()),
0008 tokenGainCalib_(esConsumes()),
0009 verbosity_(conf.getUntrackedParameter<int>("RPixVerbosity")),
0010 clusterizer_(conf) {
0011 produces<edm::DetSetVector<CTPPSPixelCluster> >();
0012 }
0013
0014 CTPPSPixelClusterProducer::~CTPPSPixelClusterProducer() {}
0015
0016 void CTPPSPixelClusterProducer::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0017 edm::ParameterSetDescription desc;
0018 desc.addUntracked<int>("RPixVerbosity", 0);
0019 desc.add<edm::InputTag>("tag", edm::InputTag("ctppsPixelDigis"));
0020 desc.add<int>("SeedADCThreshold", 2);
0021 desc.add<int>("ADCThreshold", 2);
0022 desc.add<double>("ElectronADCGain", 135.0);
0023 desc.add<int>("VCaltoElectronGain", 50);
0024 desc.add<int>("VCaltoElectronOffset", -411);
0025 desc.add<bool>("doSingleCalibration", false);
0026 descriptions.add("ctppsPixelClusters", desc);
0027 }
0028
0029 void CTPPSPixelClusterProducer::produce(edm::Event &iEvent, const edm::EventSetup &iSetup) {
0030
0031 edm::Handle<edm::DetSetVector<CTPPSPixelDigi> > rpd;
0032 iEvent.getByToken(tokenCTPPSPixelDigi_, rpd);
0033
0034 edm::DetSetVector<CTPPSPixelCluster> output;
0035
0036 if (!rpd->empty()) {
0037
0038 const auto &mask = iSetup.getData(tokenCTPPSPixelAnalysisMask_);
0039
0040 const auto &gainCalibrations = iSetup.getData(tokenGainCalib_);
0041
0042 run(*rpd, output, mask, gainCalibrations);
0043 }
0044
0045 iEvent.put(std::make_unique<edm::DetSetVector<CTPPSPixelCluster> >(output));
0046 }
0047
0048 void CTPPSPixelClusterProducer::run(const edm::DetSetVector<CTPPSPixelDigi> &input,
0049 edm::DetSetVector<CTPPSPixelCluster> &output,
0050 const CTPPSPixelAnalysisMask &mask,
0051 const CTPPSPixelGainCalibrations &gainCalibration) {
0052 for (const auto &ds_digi : input) {
0053 edm::DetSet<CTPPSPixelCluster> &ds_cluster = output.find_or_insert(ds_digi.id);
0054 clusterizer_.buildClusters(ds_digi.id, ds_digi.data, ds_cluster.data, &gainCalibration, &mask);
0055
0056 if (verbosity_) {
0057 unsigned int cluN = 0;
0058 for (std::vector<CTPPSPixelCluster>::iterator iit = ds_cluster.data.begin(); iit != ds_cluster.data.end();
0059 iit++) {
0060 edm::LogInfo("CTPPSPixelClusterProducer") << "Cluster " << ++cluN << " avg row " << (*iit).avg_row()
0061 << " avg col " << (*iit).avg_col() << " ADC.size " << (*iit).size();
0062 }
0063 }
0064 }
0065 }
0066
0067 DEFINE_FWK_MODULE(CTPPSPixelClusterProducer);