File indexing completed on 2024-04-06 12:18:39
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019 #include "HLTDTActivityFilter.h"
0020
0021
0022 #include <vector>
0023 #include <string>
0024 #include <map>
0025 #include <iostream>
0026 #include <memory>
0027
0028
0029 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0030 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0031 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
0032 #include "DataFormats/GeometryVector/interface/Pi.h"
0033
0034
0035 typedef std::map<uint32_t, std::bitset<4> > activityMap;
0036
0037
0038
0039
0040 HLTDTActivityFilter::HLTDTActivityFilter(const edm::ParameterSet& iConfig)
0041 : HLTFilter(iConfig), muonGeometryRecordToken_(esConsumes()) {
0042 using namespace std;
0043
0044 inputTag_[DCC] = iConfig.getParameter<edm::InputTag>("inputDCC");
0045 inputTag_[DDU] = iConfig.getParameter<edm::InputTag>("inputDDU");
0046 inputTag_[RPC] = iConfig.getParameter<edm::InputTag>("inputRPC");
0047 inputTag_[DIGI] = iConfig.getParameter<edm::InputTag>("inputDigis");
0048
0049 process_[DCC] = iConfig.getParameter<bool>("processDCC");
0050 process_[DDU] = iConfig.getParameter<bool>("processDDU");
0051 process_[RPC] = iConfig.getParameter<bool>("processRPC");
0052 process_[DIGI] = iConfig.getParameter<bool>("processDigis");
0053
0054 orTPG_ = iConfig.getParameter<bool>("orTPG");
0055 orRPC_ = iConfig.getParameter<bool>("orRPC");
0056 orDigi_ = iConfig.getParameter<bool>("orDigi");
0057
0058 minBX_[DCC] = iConfig.getParameter<int>("minDCCBX");
0059 maxBX_[DCC] = iConfig.getParameter<int>("maxDCCBX");
0060 minBX_[DDU] = iConfig.getParameter<int>("minDDUBX");
0061 maxBX_[DDU] = iConfig.getParameter<int>("maxDDUBX");
0062 minBX_[RPC] = iConfig.getParameter<int>("minRPCBX");
0063 maxBX_[RPC] = iConfig.getParameter<int>("maxRPCBX");
0064
0065 minQual_ = iConfig.getParameter<int>("minTPGQual");
0066 maxStation_ = iConfig.getParameter<int>("maxStation");
0067 minChambLayers_ = iConfig.getParameter<int>("minChamberLayers");
0068 minActiveChambs_ = iConfig.getParameter<int>("minActiveChambs");
0069
0070 maxDeltaPhi_ = iConfig.getParameter<double>("maxDeltaPhi");
0071 maxDeltaEta_ = iConfig.getParameter<double>("maxDeltaEta");
0072
0073 activeSecs_.reset();
0074 vector<int> aSectors = iConfig.getParameter<vector<int> >("activeSectors");
0075 vector<int>::const_iterator iSec = aSectors.begin();
0076 vector<int>::const_iterator eSec = aSectors.end();
0077 for (; iSec != eSec; ++iSec)
0078 if ((*iSec) > 0 && (*iSec < 15))
0079 activeSecs_.set((*iSec));
0080
0081 inputDCCToken_ = consumes<L1MuDTChambPhContainer>(inputTag_[DCC]);
0082 inputDDUToken_ = consumes<DTLocalTriggerCollection>(inputTag_[DDU]);
0083 inputRPCToken_ = consumes<L1MuGMTReadoutCollection>(inputTag_[RPC]);
0084 inputDigiToken_ = consumes<DTDigiCollection>(inputTag_[DIGI]);
0085 }
0086
0087 HLTDTActivityFilter::~HLTDTActivityFilter() = default;
0088
0089 void HLTDTActivityFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0090 edm::ParameterSetDescription desc;
0091 makeHLTFilterDescription(desc);
0092 desc.add<edm::InputTag>("inputDCC", edm::InputTag("hltDTTFUnpacker"));
0093 desc.add<edm::InputTag>("inputDDU", edm::InputTag("hltMuonDTDigis"));
0094 desc.add<edm::InputTag>("inputRPC", edm::InputTag("hltGtDigis"));
0095 desc.add<edm::InputTag>("inputDigis", edm::InputTag("hltMuonDTDigis"));
0096 desc.add<bool>("processDCC", true);
0097 desc.add<bool>("processDDU", true);
0098 desc.add<bool>("processRPC", true);
0099 desc.add<bool>("processDigis", true);
0100 desc.add<bool>("orTPG", true);
0101 desc.add<bool>("orRPC", true);
0102 desc.add<bool>("orDigi", false)->setComment(" # && of trig & digi info");
0103 desc.add<int>("minDCCBX", -1);
0104 desc.add<int>("maxDCCBX", 1);
0105 desc.add<int>("minDDUBX", 8);
0106 desc.add<int>("maxDDUBX", 13);
0107 desc.add<int>("minRPCBX", -1);
0108 desc.add<int>("maxRPCBX", 1);
0109 desc.add<int>("minTPGQual", 2.)->setComment(" # 0-1=L 2-3=H 4=LL 5=HL 6=HH");
0110 desc.add<int>("maxStation", 3.);
0111 desc.add<int>("minChamberLayers", 5);
0112 desc.add<int>("minActiveChambs", 1);
0113 desc.add<double>("MaxDeltaPhi", 1.0);
0114 desc.add<double>("MaxDeltaEta", 0.3);
0115 std::vector<int> temp;
0116 for (int i = 1; i <= 12; i++)
0117 temp.push_back(i);
0118 desc.add<std::vector<int> >("activeSectors", temp);
0119 descriptions.add("hltDTActivityFilter", desc);
0120 }
0121
0122
0123
0124
0125
0126
0127 bool HLTDTActivityFilter::hltFilter(edm::Event& iEvent,
0128 const edm::EventSetup& iSetup,
0129 trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0130 using namespace edm;
0131 using namespace std;
0132
0133 activityMap actMap;
0134
0135 if (process_[DCC]) {
0136 edm::Handle<L1MuDTChambPhContainer> l1DTTPGPh;
0137 iEvent.getByToken(inputDCCToken_, l1DTTPGPh);
0138 vector<L1MuDTChambPhDigi> const* phTrigs = l1DTTPGPh->getContainer();
0139 auto iph = phTrigs->begin();
0140 auto iphe = phTrigs->end();
0141
0142 for (; iph != iphe; ++iph) {
0143 int qual = iph->code();
0144 int bx = iph->bxNum();
0145 int ch = iph->stNum();
0146 int sec = iph->scNum() + 1;
0147 int wh = iph->whNum();
0148
0149 if (!activeSecs_[sec])
0150 continue;
0151
0152 if (ch <= maxStation_ && bx >= minBX_[DCC] && bx <= maxBX_[DCC] && qual >= minQual_ && qual < 7) {
0153 actMap[DTChamberId(wh, ch, sec).rawId()].set(DCC);
0154 }
0155 }
0156 }
0157
0158 if (process_[DDU]) {
0159 Handle<DTLocalTriggerCollection> trigsDDU;
0160 iEvent.getByToken(inputDDUToken_, trigsDDU);
0161 DTLocalTriggerCollection::DigiRangeIterator detUnitIt;
0162
0163 for (detUnitIt = trigsDDU->begin(); detUnitIt != trigsDDU->end(); ++detUnitIt) {
0164 int ch = (*detUnitIt).first.station();
0165 if (!activeSecs_[(*detUnitIt).first.sector()])
0166 continue;
0167
0168 const DTLocalTriggerCollection::Range& range = (*detUnitIt).second;
0169
0170 for (DTLocalTriggerCollection::const_iterator trigIt = range.first; trigIt != range.second; ++trigIt) {
0171 int bx = trigIt->bx();
0172 int qual = trigIt->quality();
0173 if (ch <= maxStation_ && bx >= minBX_[DDU] && bx <= maxBX_[DDU] && qual >= minQual_ && qual < 7) {
0174 actMap[(*detUnitIt).first.rawId()].set(DDU);
0175 }
0176 }
0177 }
0178 }
0179
0180 if (process_[DIGI]) {
0181 edm::Handle<DTDigiCollection> dtdigis;
0182 iEvent.getByToken(inputDigiToken_, dtdigis);
0183 std::map<uint32_t, int> hitMap;
0184 DTDigiCollection::DigiRangeIterator dtLayerIdIt;
0185
0186 for (dtLayerIdIt = dtdigis->begin(); dtLayerIdIt != dtdigis->end(); dtLayerIdIt++) {
0187 DTChamberId chId = ((*dtLayerIdIt).first).chamberId();
0188 if (!activeSecs_[(*dtLayerIdIt).first.sector()])
0189 continue;
0190 uint32_t rawId = chId.rawId();
0191 int station = chId.station();
0192
0193 if (station <= maxStation_) {
0194 if (hitMap.find(rawId) != hitMap.end()) {
0195 hitMap[rawId]++;
0196 } else {
0197 hitMap[rawId] = 1;
0198 }
0199 if (hitMap[rawId] >= minChambLayers_) {
0200 actMap[chId.rawId()].set(DIGI);
0201 }
0202 }
0203 }
0204 }
0205
0206 if (process_[RPC]) {
0207 auto const& dtGeom = iSetup.getHandle(muonGeometryRecordToken_);
0208
0209 edm::Handle<L1MuGMTReadoutCollection> gmtrc;
0210 iEvent.getByToken(inputRPCToken_, gmtrc);
0211
0212 std::vector<L1MuGMTReadoutRecord> gmtrr = gmtrc->getRecords();
0213 std::vector<L1MuGMTReadoutRecord>::const_iterator recIt = gmtrr.begin();
0214 std::vector<L1MuGMTReadoutRecord>::const_iterator recEnd = gmtrr.end();
0215
0216 for (; recIt != recEnd; ++recIt) {
0217 std::vector<L1MuRegionalCand> rpcCands = (*recIt).getBrlRPCCands();
0218 std::vector<L1MuRegionalCand>::const_iterator candIt = rpcCands.begin();
0219 std::vector<L1MuRegionalCand>::const_iterator candEnd = rpcCands.end();
0220
0221 for (; candIt != candEnd; ++candIt) {
0222 if (candIt->empty())
0223 continue;
0224 int bx = (*candIt).bx();
0225
0226 if (bx >= minBX_[RPC] && bx <= maxBX_[RPC]) {
0227 auto actMapIt = actMap.begin();
0228 auto actMapEnd = actMap.end();
0229 for (; actMapIt != actMapEnd; ++actMapIt)
0230 if (matchChamber(actMapIt->first, *candIt, dtGeom.product()))
0231 actMapIt->second.set(RPC);
0232 }
0233 }
0234 }
0235 }
0236
0237 int nActCh = 0;
0238 activityMap::const_iterator actMapIt = actMap.begin();
0239 activityMap::const_iterator actMapEnd = actMap.end();
0240
0241 for (; actMapIt != actMapEnd; ++actMapIt)
0242 hasActivity((*actMapIt).second) && nActCh++;
0243
0244 bool result = nActCh >= minActiveChambs_;
0245
0246 return result;
0247 }
0248
0249 bool HLTDTActivityFilter::hasActivity(const std::bitset<4>& actWord) const {
0250 bool actTPG = orTPG_ ? actWord[DCC] || actWord[DDU] : actWord[DCC] && actWord[DDU];
0251 bool actTrig = orRPC_ ? actWord[RPC] || actTPG : actWord[RPC] && actTPG;
0252 bool result = orDigi_ ? actWord[DIGI] || actTrig : actWord[DIGI] && actTrig;
0253
0254 return result;
0255 }
0256
0257 bool HLTDTActivityFilter::matchChamber(uint32_t rawId,
0258 L1MuRegionalCand const& rpcTrig,
0259 DTGeometry const* dtGeom) const {
0260 const GlobalPoint chPos = dtGeom->chamber(DTChamberId(rawId))->position();
0261
0262 float fDeltaPhi = fabs(chPos.phi() - rpcTrig.phiValue());
0263 if (fDeltaPhi > Geom::pi())
0264 fDeltaPhi = fabs(fDeltaPhi - 2 * Geom::pi());
0265
0266 float fDeltaEta = fabs(chPos.eta() - rpcTrig.etaValue());
0267
0268 bool result = fDeltaPhi < maxDeltaPhi_ && fDeltaEta < maxDeltaEta_;
0269
0270 return result;
0271 }
0272
0273
0274 #include "FWCore/Framework/interface/MakerMacros.h"
0275 DEFINE_FWK_MODULE(HLTDTActivityFilter);