Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:18:20

0001 /** \class HLTEgammaDoubleEtDeltaPhiFilter
0002  *  This filter will require only two HLT photons and
0003  *  DeltaPhi between the two photons larger than 2.5
0004  *
0005  *  \author Li Wenbo (PKU)
0006  *
0007  */
0008 
0009 #include "HLTEgammaDoubleEtDeltaPhiFilter.h"
0010 
0011 #include "DataFormats/Common/interface/Handle.h"
0012 
0013 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0016 
0017 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
0018 
0019 //
0020 // constructors and destructor
0021 //
0022 HLTEgammaDoubleEtDeltaPhiFilter::HLTEgammaDoubleEtDeltaPhiFilter(const edm::ParameterSet& iConfig)
0023     : HLTFilter(iConfig) {
0024   inputTag_ = iConfig.getParameter<edm::InputTag>("inputTag");
0025   etcut_ = iConfig.getParameter<double>("etcut");
0026   minDeltaPhi_ = iConfig.getParameter<double>("minDeltaPhi");
0027   l1EGTag_ = iConfig.getParameter<edm::InputTag>("l1EGCand");
0028   inputToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(inputTag_);
0029 }
0030 
0031 HLTEgammaDoubleEtDeltaPhiFilter::~HLTEgammaDoubleEtDeltaPhiFilter() = default;
0032 
0033 void HLTEgammaDoubleEtDeltaPhiFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0034   edm::ParameterSetDescription desc;
0035   makeHLTFilterDescription(desc);
0036   desc.add<edm::InputTag>("inputTag", edm::InputTag("hltDoublePhotonEt5L1MatchFilterRegional"));
0037   desc.add<edm::InputTag>("l1EGCand", edm::InputTag("hltL1IsoRecoEcalCandidate"));
0038   desc.add<double>("etcut", 5.0);
0039   desc.add<double>("minDeltaPhi", 2.5);
0040   descriptions.add("hltEgammaDoubleEtDeltaPhiFilter", desc);
0041 }
0042 
0043 // ------------ method called to produce the data  ------------
0044 bool HLTEgammaDoubleEtDeltaPhiFilter::hltFilter(edm::Event& iEvent,
0045                                                 const edm::EventSetup& iSetup,
0046                                                 trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0047   using namespace trigger;
0048   // The filter object
0049   if (saveTags()) {
0050     filterproduct.addCollectionTag(l1EGTag_);
0051   }
0052 
0053   // get hold of filtered candidates
0054   edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
0055   iEvent.getByToken(inputToken_, PrevFilterOutput);
0056 
0057   std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > recoecalcands;
0058   PrevFilterOutput->getObjects(TriggerCluster, recoecalcands);
0059   if (recoecalcands.empty())
0060     PrevFilterOutput->getObjects(TriggerPhoton,
0061                                  recoecalcands);  //we dont know if its type trigger cluster or trigger photon
0062 
0063   // Refs to the two Candidate objects used to calculate deltaPhi
0064   edm::Ref<reco::RecoEcalCandidateCollection> ref1, ref2;
0065 
0066   // look at all candidates,  check cuts
0067   int n(0);
0068   for (auto& ref : recoecalcands) {
0069     if (ref->et() >= etcut_) {
0070       ++n;
0071       if (n == 1)
0072         ref1 = ref;
0073       if (n == 2)
0074         ref2 = ref;
0075     }
0076   }
0077 
0078   // if there are only two Candidate objects, calculate deltaPhi
0079   double deltaPhi(0.0);
0080   if (n == 2) {
0081     deltaPhi = fabs(ref1->phi() - ref2->phi());
0082     if (deltaPhi > M_PI)
0083       deltaPhi = 2 * M_PI - deltaPhi;
0084 
0085     filterproduct.addObject(TriggerCluster, ref1);
0086     filterproduct.addObject(TriggerCluster, ref2);
0087   }
0088 
0089   // filter decision
0090   bool accept(n == 2 && deltaPhi > minDeltaPhi_);
0091 
0092   return accept;
0093 }
0094 
0095 // define as a framework module
0096 // #include "FWCore/Framework/interface/MakerMacros.h"
0097 // DEFINE_FWK_MODULE(HLTEgammaDoubleEtDeltaPhiFilter);
0098 
0099 // declare this class as a framework plugin
0100 #include "FWCore/Framework/interface/MakerMacros.h"
0101 DEFINE_FWK_MODULE(HLTEgammaDoubleEtDeltaPhiFilter);