File indexing completed on 2024-04-06 12:29:22
0001 #include "CalibFormats/CastorObjects/interface/CastorCalibrations.h"
0002 #include "CalibFormats/CastorObjects/interface/CastorCoderDb.h"
0003 #include "CalibFormats/CastorObjects/interface/CastorDbRecord.h"
0004 #include "CalibFormats/CastorObjects/interface/CastorDbService.h"
0005 #include "DataFormats/Common/interface/Handle.h"
0006 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
0007 #include "DataFormats/HcalDigi/interface/HcalTTPDigi.h"
0008 #include "DataFormats/L1GlobalTrigger/interface/L1GtTechnicalTriggerRecord.h"
0009 #include "FWCore/Framework/interface/one/EDProducer.h"
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "FWCore/Framework/interface/EventSetup.h"
0012 #include "FWCore/Framework/interface/ESHandle.h"
0013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 #include "FWCore/Utilities/interface/ESGetToken.h"
0016
0017 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
0018
0019 class CastorTTRecord : public edm::one::EDProducer<> {
0020 public:
0021 explicit CastorTTRecord(const edm::ParameterSet &ps);
0022 ~CastorTTRecord() override;
0023
0024 void produce(edm::Event &e, const edm::EventSetup &c) override;
0025
0026
0027
0028 void getEnergy_fC(double energy[16][14],
0029 edm::Handle<CastorDigiCollection> &CastorDigiColl,
0030 edm::Event &e,
0031 const edm::EventSetup &c);
0032
0033
0034 void getTriggerDecisions(std::vector<bool> &decision, double energy[16][14]) const;
0035
0036
0037
0038 void getTriggerDecisionsPerOctant(std::vector<bool> tdps[16], double energy[16][14]) const;
0039
0040 private:
0041 edm::EDGetTokenT<CastorDigiCollection> CastorDigiColl_;
0042 unsigned int CastorSignalTS_;
0043
0044 std::vector<unsigned int> ttpBits_;
0045 std::vector<std::string> TrigNames_;
0046 std::vector<double> TrigThresholds_;
0047 edm::ESGetToken<CastorDbService, CastorDbRecord> conditionsToken_;
0048
0049 double reweighted_gain;
0050 };
0051
0052 CastorTTRecord::CastorTTRecord(const edm::ParameterSet &ps) {
0053 CastorDigiColl_ = consumes<CastorDigiCollection>(ps.getParameter<edm::InputTag>("CastorDigiCollection"));
0054 CastorSignalTS_ = ps.getParameter<unsigned int>("CastorSignalTS");
0055
0056 ttpBits_ = ps.getParameter<std::vector<unsigned int>>("ttpBits");
0057 TrigNames_ = ps.getParameter<std::vector<std::string>>("TriggerBitNames");
0058 TrigThresholds_ = ps.getParameter<std::vector<double>>("TriggerThresholds");
0059 conditionsToken_ = esConsumes<CastorDbService, CastorDbRecord>();
0060 reweighted_gain = 1.0;
0061
0062 produces<L1GtTechnicalTriggerRecord>();
0063 }
0064
0065 CastorTTRecord::~CastorTTRecord() {}
0066
0067 void CastorTTRecord::produce(edm::Event &e, const edm::EventSetup &eventSetup) {
0068
0069
0070 std::vector<L1GtTechnicalTrigger> vecTT(ttpBits_.size());
0071
0072
0073 edm::Handle<CastorDigiCollection> CastorDigiColl;
0074 e.getByToken(CastorDigiColl_, CastorDigiColl);
0075
0076 if (!CastorDigiColl.failedToGet()) {
0077 double cas_efC[16][14];
0078 getEnergy_fC(cas_efC, CastorDigiColl, e, eventSetup);
0079
0080 std::vector<bool> decision(ttpBits_.size());
0081
0082 getTriggerDecisions(decision, cas_efC);
0083
0084 for (unsigned int i = 0; i < ttpBits_.size(); i++) {
0085
0086
0087
0088
0089 vecTT.at(i) = L1GtTechnicalTrigger(TrigNames_.at(i), ttpBits_.at(i), 0, decision.at(i));
0090 }
0091
0092 } else {
0093 vecTT.clear();
0094 }
0095
0096
0097 std::unique_ptr<L1GtTechnicalTriggerRecord> output(new L1GtTechnicalTriggerRecord());
0098 output->setGtTechnicalTrigger(vecTT);
0099 e.put(std::move(output));
0100 }
0101
0102 void CastorTTRecord::getEnergy_fC(double energy[16][14],
0103 edm::Handle<CastorDigiCollection> &CastorDigiColl,
0104 edm::Event &e,
0105 const edm::EventSetup &eventSetup) {
0106
0107
0108
0109
0110 edm::ESHandle<CastorDbService> conditions = eventSetup.getHandle(conditionsToken_);
0111 const CastorQIEShape *shape = conditions->getCastorShape();
0112
0113 for (int isec = 0; isec < 16; isec++)
0114 for (int imod = 0; imod < 14; imod++)
0115 energy[isec][imod] = 0;
0116
0117
0118 CastorDigiCollection::const_iterator idigi;
0119 for (idigi = CastorDigiColl->begin(); idigi != CastorDigiColl->end(); idigi++) {
0120 const CastorDataFrame &digi = (*idigi);
0121 HcalCastorDetId cell = digi.id();
0122
0123
0124 const CastorQIECoder *channelCoder = conditions->getCastorCoder(cell);
0125 CastorCoderDb coder(*channelCoder, *shape);
0126
0127
0128 const CastorCalibrations &calibrations = conditions->getCastorCalibrations(cell);
0129
0130
0131 CaloSamples tool;
0132 coder.adc2fC(digi, tool);
0133
0134
0135 int capid = digi[CastorSignalTS_].capid();
0136 double fC = tool[CastorSignalTS_] - calibrations.pedestal(capid);
0137
0138
0139 reweighted_gain = calibrations.gain(capid) / 0.015;
0140
0141 energy[digi.id().sector() - 1][digi.id().module() - 1] = fC;
0142 }
0143 }
0144
0145 void CastorTTRecord::getTriggerDecisions(std::vector<bool> &decision, double energy[16][14]) const {
0146
0147
0148
0149
0150 if (decision.size() < 4)
0151 return;
0152
0153 std::vector<bool> tdpo[8];
0154 getTriggerDecisionsPerOctant(tdpo, energy);
0155
0156
0157 decision.at(0) = true;
0158 decision.at(1) = false;
0159 decision.at(2) = false;
0160 decision.at(3) = false;
0161
0162
0163
0164
0165
0166 for (int ioct = 0; ioct < 8; ioct++) {
0167 int next_oct = (ioct + 1) % 8;
0168 int prev_oct = (ioct + 8 - 1) % 8;
0169
0170
0171 if (!tdpo[ioct].at(0))
0172 decision.at(0) = false;
0173 if (!tdpo[ioct].at(1))
0174 decision.at(0) = false;
0175
0176
0177 if (tdpo[ioct].at(2))
0178 decision.at(1) = true;
0179
0180
0181
0182
0183
0184
0185
0186 if (tdpo[ioct].at(5)) {
0187
0188
0189 if (tdpo[ioct].at(0)) {
0190 if (tdpo[prev_oct].at(1) && tdpo[next_oct].at(0) && tdpo[next_oct].at(1))
0191 decision.at(3) = true;
0192 } else if (tdpo[ioct].at(1)) {
0193 if (tdpo[prev_oct].at(0) && tdpo[prev_oct].at(1) && tdpo[next_oct].at(0))
0194 decision.at(3) = true;
0195 }
0196
0197 }
0198
0199
0200 if (tdpo[ioct].at(6))
0201 decision.at(2) = true;
0202 }
0203
0204
0205
0206
0207 }
0208
0209 void CastorTTRecord::getTriggerDecisionsPerOctant(std::vector<bool> tdpo[8], double energy[16][14]) const {
0210
0211
0212
0213
0214 for (int ioct = 0; ioct < 8; ioct++) {
0215
0216
0217
0218
0219
0220
0221
0222
0223
0224 tdpo[ioct].resize(7);
0225
0226 for (int ibit = 0; ibit < 7; ibit++)
0227 tdpo[ioct].at(ibit) = false;
0228
0229
0230 for (int ioctsec = 0; ioctsec < 2; ioctsec++) {
0231
0232 int isec = 2 * ioct + ioctsec;
0233
0234
0235 double fCsum_mod = 0;
0236 double fCsum_em = 0, fCsum_ha = 0;
0237 double fCsum_jet_had = 0;
0238 double fCsum_col[3] = {0, 0, 0};
0239
0240
0241 for (int imod = 0; imod < 14; imod++) {
0242
0243 fCsum_mod += energy[isec][imod];
0244
0245
0246 if (imod < 2)
0247 fCsum_em += energy[isec][imod];
0248 if (imod > 2 && imod < 12)
0249 fCsum_ha += energy[isec][imod];
0250
0251
0252 if (imod < 4)
0253 fCsum_col[0] += energy[isec][imod];
0254 else if (imod < 8)
0255 fCsum_col[1] += energy[isec][imod];
0256 else if (imod < 12)
0257 fCsum_col[2] += energy[isec][imod];
0258
0259
0260 if (imod > 1 && imod < 5)
0261 fCsum_jet_had += energy[isec][imod];
0262 }
0263
0264
0265 if (fCsum_mod < TrigThresholds_.at(0)) {
0266 if (ioctsec == 0)
0267 tdpo[ioct].at(0) = true;
0268 else if (ioctsec == 1)
0269 tdpo[ioct].at(1) = true;
0270 }
0271
0272
0273
0274
0275
0276
0277
0278 if (fCsum_jet_had > TrigThresholds_.at(1) / reweighted_gain)
0279
0280 if (energy[isec][0] > 26000 / reweighted_gain && energy[isec][1] > 26000 / reweighted_gain)
0281 tdpo[ioct].at(2) = true;
0282
0283
0284 if (fCsum_mod > TrigThresholds_.at(5) / reweighted_gain)
0285 tdpo[ioct].at(6) = true;
0286
0287
0288
0289 if (fCsum_em > TrigThresholds_.at(2) / reweighted_gain)
0290 tdpo[ioct].at(3) = true;
0291 if (fCsum_ha > TrigThresholds_.at(3))
0292 tdpo[ioct].at(4) = true;
0293
0294
0295 int countColumns = 0;
0296 for (int icol = 0; icol < 3; icol++)
0297 if (fCsum_col[icol] > TrigThresholds_.at(4))
0298 countColumns++;
0299 if (countColumns >= 2)
0300 tdpo[ioct].at(5) = true;
0301 }
0302 }
0303 }
0304
0305 #include "FWCore/Framework/interface/MakerMacros.h"
0306 #include "FWCore/PluginManager/interface/ModuleDef.h"
0307
0308 DEFINE_FWK_MODULE(CastorTTRecord);