Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023 #include "HLTrigger/special/interface/HLTDTActivityFilter.h"
00024
00025
00026 #include <vector>
00027 #include <string>
00028 #include <map>
00029 #include <iostream>
00030 #include <memory>
00031
00032
00033 #include "DataFormats/L1GlobalMuonTrigger/interface/L1MuGMTReadoutCollection.h"
00034 #include "DataFormats/L1GlobalMuonTrigger/interface/L1MuRegionalCand.h"
00035 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambPhContainer.h"
00036 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambThContainer.h"
00037 #include "DataFormats/DTDigi/interface/DTLocalTriggerCollection.h"
00038 #include "DataFormats/DTDigi/interface/DTDigiCollection.h"
00039
00040 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00041 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00042 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00043 #include "DataFormats/GeometryVector/interface/Pi.h"
00044
00045
00046
00047 typedef std::map<uint32_t,std::bitset<4> > activityMap;
00048
00049
00050
00051
00052
00053 HLTDTActivityFilter::HLTDTActivityFilter(const edm::ParameterSet& iConfig) {
00054
00055 using namespace std;
00056
00057 inputTag_[DCC] = iConfig.getParameter<edm::InputTag>("inputDCC");
00058 inputTag_[DDU] = iConfig.getParameter<edm::InputTag>("inputDDU");
00059 inputTag_[RPC] = iConfig.getParameter<edm::InputTag>("inputRPC");
00060 inputTag_[DIGI] = iConfig.getParameter<edm::InputTag>("inputDigis");
00061
00062 process_[DCC] = iConfig.getParameter<bool>("processDCC");
00063 process_[DDU] = iConfig.getParameter<bool>("processDDU");
00064 process_[RPC] = iConfig.getParameter<bool>("processRPC");
00065 process_[DIGI] = iConfig.getParameter<bool>("processDigis");
00066
00067 orTPG_ = iConfig.getParameter<bool>("orTPG");
00068 orRPC_ = iConfig.getParameter<bool>("orRPC");
00069 orDigi_ = iConfig.getParameter<bool>("orDigi");
00070
00071 minBX_[DCC] = iConfig.getParameter<int>("minDCCBX");
00072 maxBX_[DCC] = iConfig.getParameter<int>("maxDCCBX");
00073 minBX_[DDU] = iConfig.getParameter<int>("minDDUBX");
00074 maxBX_[DDU] = iConfig.getParameter<int>("maxDDUBX");
00075 minBX_[RPC] = iConfig.getParameter<int>("minRPCBX");
00076 maxBX_[RPC] = iConfig.getParameter<int>("maxRPCBX");
00077
00078 minQual_ = iConfig.getParameter<int>("minTPGQual");
00079 maxStation_ = iConfig.getParameter<int>("maxStation");
00080 minChambLayers_ = iConfig.getParameter<int>("minChamberLayers");
00081 minActiveChambs_ = iConfig.getParameter<int>("minActiveChambs");
00082
00083 maxDeltaPhi_ = iConfig.getParameter<double>("maxDeltaPhi");
00084 maxDeltaEta_ = iConfig.getParameter<double>("maxDeltaEta");
00085
00086
00087 activeSecs_.reset();
00088 vector<int> aSectors = iConfig.getParameter<vector<int> >("activeSectors");
00089 vector<int>::const_iterator iSec = aSectors.begin();
00090 vector<int>::const_iterator eSec = aSectors.end();
00091 for (;iSec!=eSec;++iSec)
00092 if ((*iSec)>0 && (*iSec<15)) activeSecs_.set((*iSec));
00093
00094 }
00095
00096
00097 HLTDTActivityFilter::~HLTDTActivityFilter() {
00098
00099 }
00100
00101
00102
00103
00104
00105
00106 bool HLTDTActivityFilter::beginRun(edm::Run& iRun, const edm::EventSetup& iSetup) {
00107
00108 iSetup.get<MuonGeometryRecord>().get(dtGeom_);
00109
00110 return true;
00111
00112 }
00113
00114
00115 bool HLTDTActivityFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) {
00116
00117 using namespace edm;
00118 using namespace std;
00119
00120 activityMap actMap;
00121
00122 if (process_[DCC]) {
00123
00124 edm::Handle<L1MuDTChambPhContainer> l1DTTPGPh;
00125 iEvent.getByLabel(inputTag_[DCC],l1DTTPGPh);
00126 vector<L1MuDTChambPhDigi>* phTrigs = l1DTTPGPh->getContainer();
00127 vector<L1MuDTChambPhDigi>::const_iterator iph = phTrigs->begin();
00128 vector<L1MuDTChambPhDigi>::const_iterator iphe = phTrigs->end();
00129
00130 for(; iph !=iphe ; ++iph) {
00131
00132 int qual = iph->code();
00133 int bx = iph->bxNum();
00134 int ch = iph->stNum();
00135 int sec = iph->scNum() + 1;
00136 int wh = iph->whNum();
00137
00138 if (!activeSecs_[sec]) continue;
00139
00140 if (ch<=maxStation_ && bx>=minBX_[DCC] && bx<=maxBX_[DCC]
00141 && qual>=minQual_ && qual<7) {
00142 actMap[DTChamberId(wh,ch,sec).rawId()].set(DCC);
00143 }
00144
00145 }
00146
00147 }
00148
00149 if (process_[DDU]) {
00150
00151 Handle<DTLocalTriggerCollection> trigsDDU;
00152 iEvent.getByLabel(inputTag_[DDU],trigsDDU);
00153 DTLocalTriggerCollection::DigiRangeIterator detUnitIt;
00154
00155 for (detUnitIt=trigsDDU->begin();detUnitIt!=trigsDDU->end();++detUnitIt){
00156
00157 int ch = (*detUnitIt).first.station();
00158 if (!activeSecs_[(*detUnitIt).first.sector()]) continue;
00159
00160 const DTLocalTriggerCollection::Range& range = (*detUnitIt).second;
00161
00162 for (DTLocalTriggerCollection::const_iterator trigIt = range.first; trigIt!=range.second;++trigIt){
00163 int bx = trigIt->bx();
00164 int qual = trigIt->quality();
00165 if ( ch<=maxStation_ && bx>=minBX_[DDU] && bx<=maxBX_[DDU]
00166 && qual>=minQual_ && qual<7) {
00167 actMap[(*detUnitIt).first.rawId()].set(DDU);
00168 }
00169
00170 }
00171
00172 }
00173
00174 }
00175
00176 if (process_[DIGI]) {
00177
00178 edm::Handle<DTDigiCollection> dtdigis;
00179 iEvent.getByLabel(inputTag_[DIGI], dtdigis);
00180 std::map<uint32_t,int> hitMap;
00181 DTDigiCollection::DigiRangeIterator dtLayerIdIt;
00182
00183 for (dtLayerIdIt=dtdigis->begin(); dtLayerIdIt!=dtdigis->end(); dtLayerIdIt++) {
00184
00185 DTChamberId chId = ((*dtLayerIdIt).first).chamberId();
00186 if (!activeSecs_[(*dtLayerIdIt).first.sector()]) continue;
00187 uint32_t rawId = chId.rawId();
00188 int station = chId.station();
00189
00190 if (station<=maxStation_) {
00191 if (hitMap.find(rawId)!=hitMap.end()) {
00192 hitMap[rawId]++;
00193 } else {
00194 hitMap[rawId]=1;
00195 }
00196 if (hitMap[rawId]>=minChambLayers_) {
00197 actMap[chId.rawId()].set(DIGI);
00198 }
00199 }
00200
00201 }
00202
00203 }
00204
00205 if (process_[RPC]) {
00206
00207 edm::Handle<L1MuGMTReadoutCollection> gmtrc;
00208 iEvent.getByLabel(inputTag_[RPC],gmtrc);
00209
00210 std::vector<L1MuGMTReadoutRecord> gmtrr = gmtrc->getRecords();
00211 std::vector<L1MuGMTReadoutRecord>::const_iterator recIt = gmtrr.begin();
00212 std::vector<L1MuGMTReadoutRecord>::const_iterator recEnd = gmtrr.end();
00213
00214 for(; recIt!=recEnd; ++recIt) {
00215
00216 std::vector<L1MuRegionalCand> rpcCands = (*recIt).getBrlRPCCands();
00217 std::vector<L1MuRegionalCand>::const_iterator candIt = rpcCands.begin();
00218 std::vector<L1MuRegionalCand>::const_iterator candEnd = rpcCands.end();
00219
00220 for(; candIt!=candEnd; ++candIt) {
00221
00222 if (candIt->empty()) continue;
00223 int bx = (*candIt).bx();
00224
00225 if (bx>=minBX_[RPC] && bx<=maxBX_[RPC]) {
00226 activityMap::iterator actMapIt = actMap.begin();
00227 activityMap::iterator actMapEnd = actMap.end();
00228 for (; actMapIt!= actMapEnd; ++ actMapIt)
00229 if (matchChamber((*actMapIt).first,(*candIt)))
00230 (*actMapIt).second.set(RPC);
00231 }
00232 }
00233 }
00234
00235 }
00236
00237 int nActCh = 0;
00238 activityMap::const_iterator actMapIt = actMap.begin();
00239 activityMap::const_iterator actMapEnd = actMap.end();
00240
00241 for (; actMapIt!=actMapEnd; ++actMapIt)
00242 hasActivity((*actMapIt).second) && nActCh++ ;
00243
00244 bool result = nActCh>=minActiveChambs_;
00245
00246 return result;
00247
00248 }
00249
00250
00251 bool HLTDTActivityFilter::hasActivity(const std::bitset<4>& actWord) {
00252
00253 bool actTPG = orTPG_ ? actWord[DCC] || actWord[DDU] : actWord[DCC] && actWord[DDU];
00254 bool actTrig = orRPC_ ? actWord[RPC] || actTPG : actWord[RPC] && actTPG;
00255 bool result = orDigi_ ? actWord[DIGI] || actTrig : actWord[DIGI] && actTrig;
00256
00257 return result;
00258
00259 }
00260
00261 bool HLTDTActivityFilter::matchChamber(const uint32_t& rawId, const L1MuRegionalCand& rpcTrig) {
00262
00263 const GlobalPoint chPos = dtGeom_->chamber(DTChamberId(rawId))->position();
00264
00265 float fDeltaPhi = fabs( chPos.phi() - rpcTrig.phiValue() );
00266 if ( fDeltaPhi>Geom::pi() ) fDeltaPhi = fabs(fDeltaPhi - 2*Geom::pi());
00267
00268 float fDeltaEta = fabs( chPos.eta() - rpcTrig.etaValue() );
00269
00270 bool result = fDeltaPhi<maxDeltaPhi_ && fDeltaEta<maxDeltaEta_;
00271
00272 return result;
00273
00274 }
00275
00276
00277
00278 #include "FWCore/Framework/interface/MakerMacros.h"
00279 DEFINE_FWK_MODULE(HLTDTActivityFilter);