File indexing completed on 2024-04-06 12:07:32
0001
0002 #include "DQM/HcalTasks/interface/FCDTask.h"
0003
0004 bool operator==(const FCDTask::FCDChannel& lhs, const FCDTask::FCDChannel& rhs) {
0005 return ((lhs.crate == rhs.crate) && (lhs.slot == rhs.slot) && (lhs.fiber == rhs.fiber) &&
0006 (lhs.fiberChannel == rhs.fiberChannel));
0007 }
0008
0009 FCDTask::FCDTask(edm::ParameterSet const& ps)
0010 : hcalDbServiceToken_(esConsumes<HcalDbService, HcalDbRecord, edm::Transition::BeginRun>()) {
0011
0012 _tagQIE10 = ps.getUntrackedParameter<edm::InputTag>("tagQIE10", edm::InputTag("hcalDigis", "ZDC"));
0013 _tokQIE10 = consumes<QIE10DigiCollection>(_tagQIE10);
0014
0015
0016 edm::ParameterSet channelPSet = ps.getParameter<edm::ParameterSet>("fcdChannels");
0017 std::vector<int32_t> crates = channelPSet.getUntrackedParameter<std::vector<int32_t> >("crate");
0018 std::vector<int32_t> slots = channelPSet.getUntrackedParameter<std::vector<int32_t> >("slot");
0019 std::vector<int32_t> fibers = channelPSet.getUntrackedParameter<std::vector<int32_t> >("fiber");
0020 std::vector<int32_t> fiberChannels = channelPSet.getUntrackedParameter<std::vector<int32_t> >("fiber_channel");
0021 for (unsigned int i = 0; i < crates.size(); ++i) {
0022 _channels.push_back({crates[i], slots[i], fibers[i], fiberChannels[i]});
0023 }
0024 }
0025
0026 void FCDTask::bookHistograms(DQMStore::IBooker& ib, edm::Run const& r, edm::EventSetup const& es) {
0027 edm::ESHandle<HcalDbService> dbService = es.getHandle(hcalDbServiceToken_);
0028 _emap = dbService->getHcalMapping();
0029 _ehashmap.initialize(_emap, hcaldqm::electronicsmap::fD2EHashMap);
0030
0031 ib.cd();
0032
0033
0034 std::string histoname;
0035 std::vector<HcalGenericDetId> gids = _emap->allPrecisionId();
0036 for (auto& it_gid : gids) {
0037 if (it_gid.genericSubdet() != HcalGenericDetId::HcalGenZDC) {
0038 continue;
0039 }
0040 HcalElectronicsId eid = HcalElectronicsId(_ehashmap.lookup(it_gid));
0041 for (auto& it_channel : _channels) {
0042 if ((eid.crateId() == it_channel.crate) && (eid.slot() == it_channel.slot) &&
0043 (eid.fiberIndex() == it_channel.fiber) && (eid.fiberChanId() == it_channel.fiberChannel)) {
0044 _fcd_eids.push_back(eid);
0045 }
0046 }
0047 }
0048 for (auto& it_eid : _fcd_eids) {
0049
0050 histoname = std::to_string(it_eid.crateId()) + "-" + std::to_string(it_eid.slot()) + "-" +
0051 std::to_string(it_eid.fiberIndex()) + "-" + std::to_string(it_eid.fiberChanId());
0052 ib.setCurrentFolder("Hcal/FCDTask/ADC");
0053 _cADC[it_eid] = ib.book1DD(histoname.c_str(), histoname.c_str(), 256, 0, 256);
0054 _cADC[it_eid]->setAxisTitle("ADC", 1);
0055 _cADC[it_eid]->setAxisTitle("N", 2);
0056
0057 ib.setCurrentFolder("Hcal/FCDTask/ADC_vs_TS"),
0058 _cADC_vs_TS[it_eid] = ib.book2D(histoname.c_str(), histoname.c_str(), 10, 0, 10, 64, 0, 256);
0059 _cADC_vs_TS[it_eid]->setAxisTitle("TS", 1);
0060 _cADC_vs_TS[it_eid]->setAxisTitle("ADC", 2);
0061
0062 ib.setCurrentFolder("Hcal/FCDTask/TDCTime");
0063 _cTDCTime[it_eid] = ib.book1DD(histoname.c_str(), histoname.c_str(), 500, 0., 250.);
0064 _cTDCTime[it_eid]->setAxisTitle("TDC time [ns]", 1);
0065
0066 ib.setCurrentFolder("Hcal/FCDTask/TDC");
0067 _cTDC[it_eid] = ib.book1DD(histoname.c_str(), histoname.c_str(), 64, -0.5, 63.5);
0068 _cTDC[it_eid]->setAxisTitle("TDC", 1);
0069 }
0070 }
0071
0072 void FCDTask::analyze(edm::Event const& e, edm::EventSetup const&) {
0073 edm::Handle<QIE10DigiCollection> digis;
0074 if (!e.getByToken(_tokQIE10, digis))
0075 edm::LogError("Collection QIE10DigiCollection for ZDC isn't available" + _tagQIE10.label() + " " +
0076 _tagQIE10.instance());
0077
0078 for (auto it = digis->begin(); it != digis->end(); it++) {
0079 const QIE10DataFrame digi = static_cast<const QIE10DataFrame>(*it);
0080 HcalGenericDetId const& gdid = digi.detid();
0081 HcalElectronicsId eid = HcalElectronicsId(_ehashmap.lookup(gdid));
0082 if (std::find(_fcd_eids.begin(), _fcd_eids.end(), eid) == _fcd_eids.end()) {
0083 continue;
0084 }
0085
0086 for (int i = 0; i < digi.samples(); i++) {
0087
0088 _cADC[eid]->Fill(digi[i].adc());
0089 _cADC_vs_TS[eid]->Fill(i, digi[i].adc());
0090 _cTDC[eid]->Fill(digi[i].le_tdc());
0091 if (digi[i].le_tdc() <= 50.) {
0092 double tdctime = 25. * i + 0.5 * digi[i].le_tdc();
0093 _cTDCTime[eid]->Fill(tdctime);
0094 }
0095 }
0096 }
0097 }
0098
0099 DEFINE_FWK_MODULE(FCDTask);