CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/HLTrigger/Muon/src/HLTMuonL1RegionalFilter.cc

Go to the documentation of this file.
00001 #include "HLTrigger/Muon/interface/HLTMuonL1RegionalFilter.h"
00002 #include "DataFormats/Common/interface/Handle.h"
00003 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
00004 #include "DataFormats/HLTReco/interface/TriggerRefsCollections.h"
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006 #include "DataFormats/L1Trigger/interface/L1MuonParticle.h"
00007 #include "DataFormats/L1Trigger/interface/L1MuonParticleFwd.h"
00008 #include "DataFormats/L1GlobalMuonTrigger/interface/L1MuGMTCand.h"
00009 #include "FWCore/Utilities/interface/EDMException.h"
00010 
00011 HLTMuonL1RegionalFilter::HLTMuonL1RegionalFilter(const edm::ParameterSet& iConfig):
00012   candTag_( iConfig.getParameter<edm::InputTag>("CandTag") ),
00013   previousCandTag_( iConfig.getParameter<edm::InputTag>("PreviousCandTag") ),
00014   minN_( iConfig.getParameter<int>("MinN") ),
00015   saveTag_( iConfig.getUntrackedParameter<bool>("SaveTag",false) ) 
00016 {
00017   using namespace std;
00018   using namespace edm;
00019 
00020   // read in the eta-range dependent parameters  
00021   const vector<ParameterSet> cuts = iConfig.getParameter<vector<ParameterSet> >("Cuts");
00022   size_t ranges = cuts.size();
00023   if(ranges==0){
00024     throw edm::Exception(errors::Configuration) << "Please provide at least one PSet in the Cuts VPSet!";
00025   }
00026   etaBoundaries_.reserve(ranges+1);
00027   minPts_.reserve(ranges);
00028   qualityBitMasks_.reserve(ranges);
00029   for(size_t i=0; i<ranges; i++){
00030     //set the eta range
00031     vector<double> etaRange = cuts[i].getParameter<vector<double> >("EtaRange");
00032     if(etaRange.size() != 2 || etaRange[0] >= etaRange[1]){
00033       throw edm::Exception(errors::Configuration) << "EtaRange must have two non-equal values in increasing order!";
00034     }
00035     if(i==0){
00036       etaBoundaries_.push_back( etaRange[0] );
00037     }else if(etaBoundaries_[i] != etaRange[0]){
00038       throw edm::Exception(errors::Configuration) << "EtaRanges must be disjoint without gaps and listed in increasing eta order!";
00039     }
00040     etaBoundaries_.push_back( etaRange[1] );
00041 
00042     //set the minPt
00043     minPts_.push_back( cuts[i].getParameter<double>("MinPt") );
00044 
00045     //set the quality bit masks
00046     qualityBitMasks_.push_back( 0 );
00047     vector<unsigned int> qualities = cuts[i].getParameter<vector<unsigned int> >("QualityBits");
00048     for(size_t j=0; j<qualities.size(); j++){
00049       if(7U < qualities[j]){ // qualities[j] >= 0, since qualities[j] is unsigned
00050         throw edm::Exception(errors::Configuration) << "QualityBits must be between 0 and 7 !";
00051       }
00052       qualityBitMasks_[i] |= 1 << qualities[j];
00053     }
00054   }
00055 
00056   // dump parameters for debugging
00057   if(edm::isDebugEnabled()){
00058     ostringstream ss;
00059     ss<<"Constructed with parameters:"<<endl;
00060     ss<<"    CandTag = "<<candTag_.encode()<<endl;
00061     ss<<"    PreviousCandTag = "<<previousCandTag_.encode()<<endl;
00062     ss<<"    EtaBoundaries = \t"<<etaBoundaries_[0];
00063     for(size_t i=1; i<etaBoundaries_.size(); i++){
00064       ss<<'\t'<<etaBoundaries_[i];
00065     }
00066     ss<<endl;
00067     ss<<"    MinPts =        \t    "<<minPts_[0];
00068     for(size_t i=1; i<minPts_.size(); i++){
00069       ss<<"\t    "<<minPts_[i];
00070     }
00071     ss<<endl;
00072     ss<<"    QualityBitMasks =  \t    "<<qualityBitMasks_[0];
00073     for(size_t i=1; i<qualityBitMasks_.size(); i++){
00074       ss<<"\t    "<<qualityBitMasks_[i];
00075     }
00076     ss<<endl;
00077     ss<<"    MinN = "<<minN_<<endl;
00078     ss<<"    SaveTag = "<<saveTag_;
00079     LogDebug("HLTMuonL1RegionalFilter")<<ss.str();
00080   }
00081 
00082   //register the product
00083   produces<trigger::TriggerFilterObjectWithRefs>();
00084 }
00085 
00086 HLTMuonL1RegionalFilter::~HLTMuonL1RegionalFilter(){
00087 }
00088 
00089 bool HLTMuonL1RegionalFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup){
00090   using namespace std;
00091   using namespace edm;
00092   using namespace trigger;
00093   using namespace l1extra;
00094 
00095   // All HLT filters must create and fill an HLT filter object,
00096   // recording any reconstructed physics objects satisfying (or not)
00097   // this HLT filter, and place it in the Event.
00098 
00099   // The filter object
00100   auto_ptr<TriggerFilterObjectWithRefs> filterproduct(new TriggerFilterObjectWithRefs(path(), module()));
00101 
00102   // get hold of all muons
00103   Handle<L1MuonParticleCollection> allMuons;
00104   iEvent.getByLabel(candTag_, allMuons);
00105 
00106   // get hold of muons that fired the previous level
00107   Handle<TriggerFilterObjectWithRefs> previousLevelCands;
00108   iEvent.getByLabel(previousCandTag_, previousLevelCands);
00109   vector<L1MuonParticleRef> prevMuons;
00110   previousLevelCands->getObjects(TriggerL1Mu, prevMuons);
00111    
00112   // look at all mucands,  check cuts and add to filter object
00113   int n = 0;
00114   for (size_t i = 0; i < allMuons->size(); i++) {
00115     L1MuonParticleRef muon(allMuons, i);
00116 
00117     //check if triggered by the previous level
00118     if( find(prevMuons.begin(), prevMuons.end(), muon) == prevMuons.end() ) continue;
00119 
00120     //find eta region, continue otherwise
00121     float eta   =  muon->eta();
00122     int region = -1;
00123     for(size_t r=0; r<etaBoundaries_.size()-1; r++){
00124       if(etaBoundaries_[r]<=eta && eta<=etaBoundaries_[r+1]){
00125         region = r;
00126         break;
00127       }
00128     }
00129     if(region == -1) continue;
00130 
00131     //check pT cut
00132     if(muon->pt() < minPts_[region]) continue;
00133 
00134     //check quality cut
00135     if(qualityBitMasks_[region]){
00136       int quality = muon->gmtMuonCand().empty() ? 0 : (1 << muon->gmtMuonCand().quality());
00137       if((quality & qualityBitMasks_[region]) == 0) continue;
00138     }
00139 
00140     //we have a good candidate
00141     filterproduct->addObject(TriggerL1Mu, muon);
00142 
00143     n++;
00144   }
00145 
00146   if(saveTag_) filterproduct->addCollectionTag(candTag_);
00147 
00148   // filter decision
00149   const bool accept (n >= minN_);
00150 
00151   // dump event for debugging
00152   if(edm::isDebugEnabled()){
00153     ostringstream ss;
00154     ss.precision(2);
00155     ss<<"L1mu#"<<'\t'<<"q*pt"<<'\t'<<'\t'<<"eta"<<'\t'<<"phi"<<'\t'<<"quality"<<'\t'<<"isPrev"<<'\t'<<"isFired"<<endl;
00156     ss<<"---------------------------------------------------------------"<<endl;
00157 
00158     vector<L1MuonParticleRef> firedMuons;
00159     filterproduct->getObjects(TriggerL1Mu, firedMuons);
00160     for(size_t i=0; i<allMuons->size(); i++){
00161       L1MuonParticleRef mu(allMuons, i);
00162       int quality = mu->gmtMuonCand().empty() ? 0 : mu->gmtMuonCand().quality();
00163       bool isPrev = find(prevMuons.begin(), prevMuons.end(), mu) != prevMuons.end();
00164       bool isFired = find(firedMuons.begin(), firedMuons.end(), mu) != firedMuons.end();
00165       ss<<i<<'\t'<<scientific<<mu->charge()*mu->pt()<<fixed<<'\t'<<mu->eta()<<'\t'<<mu->phi()<<'\t'<<quality<<'\t'<<isPrev<<'\t'<<isFired<<endl;
00166     }
00167     ss<<"---------------------------------------------------------------"<<endl;
00168     LogDebug("HLTMuonL1RegionalFilter")<<ss.str()<<"Decision of filter is "<<accept<<", number of muons passing = "<<filterproduct->l1muonSize();
00169   }
00170 
00171   // put filter object into the Event
00172   iEvent.put(filterproduct);
00173 
00174   return accept;
00175 }
00176