Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /** \class HLTEgammaDoubleEtPhiFilter
0002  *
0003  *
0004  *  \author Jonathan Hollar (LLNL)
0005  *
0006  */
0007 
0008 #include "HLTEgammaDoubleEtPhiFilter.h"
0009 
0010 #include "DataFormats/Common/interface/Handle.h"
0011 
0012 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0015 
0016 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
0017 
0018 class EgammaHLTEtSortCriterium {
0019 public:
0020   bool operator()(edm::Ref<reco::RecoEcalCandidateCollection> lhs, edm::Ref<reco::RecoEcalCandidateCollection> rhs) {
0021     return lhs->et() > rhs->et();
0022   }
0023 };
0024 
0025 //
0026 // constructors and destructor
0027 //
0028 HLTEgammaDoubleEtPhiFilter::HLTEgammaDoubleEtPhiFilter(const edm::ParameterSet& iConfig) : HLTFilter(iConfig) {
0029   candTag_ = iConfig.getParameter<edm::InputTag>("candTag");
0030   etcut1_ = iConfig.getParameter<double>("etcut1");
0031   etcut2_ = iConfig.getParameter<double>("etcut2");
0032   min_Acop_ = iConfig.getParameter<double>("MinAcop");
0033   max_Acop_ = iConfig.getParameter<double>("MaxAcop");
0034   min_EtBalance_ = iConfig.getParameter<double>("MinEtBalance");
0035   max_EtBalance_ = iConfig.getParameter<double>("MaxEtBalance");
0036   npaircut_ = iConfig.getParameter<int>("npaircut");
0037   candToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(candTag_);
0038 }
0039 
0040 HLTEgammaDoubleEtPhiFilter::~HLTEgammaDoubleEtPhiFilter() = default;
0041 
0042 void HLTEgammaDoubleEtPhiFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0043   edm::ParameterSetDescription desc;
0044   makeHLTFilterDescription(desc);
0045   desc.add<edm::InputTag>("candTag", edm::InputTag("hltDoubleL1MatchFilter"));
0046   desc.add<double>("etcut1", 6.0);
0047   desc.add<double>("etcut2", 6.0);
0048   desc.add<double>("MinAcop", -0.1);
0049   desc.add<double>("MaxAcop", 0.6);
0050   desc.add<double>("MinEtBalance", -1.0);
0051   desc.add<double>("MaxEtBalance", 10.0);
0052   desc.add<int>("npaircut", 1);
0053   descriptions.add("hltEgammaDoubleEtPhiFilter", desc);
0054 }
0055 
0056 // ------------ method called to produce the data  ------------
0057 bool HLTEgammaDoubleEtPhiFilter::hltFilter(edm::Event& iEvent,
0058                                            const edm::EventSetup& iSetup,
0059                                            trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0060   using namespace trigger;
0061 
0062   // Ref to Candidate object to be recorded in filter object
0063   edm::Ref<reco::RecoEcalCandidateCollection> ref;
0064 
0065   edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
0066   iEvent.getByToken(candToken_, PrevFilterOutput);
0067 
0068   std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > mysortedrecoecalcands;
0069   PrevFilterOutput->getObjects(TriggerCluster, mysortedrecoecalcands);
0070   if (mysortedrecoecalcands.empty())
0071     PrevFilterOutput->getObjects(TriggerPhoton,
0072                                  mysortedrecoecalcands);  //we dont know if its type trigger cluster or trigger photon
0073   // Sort the list
0074   std::sort(mysortedrecoecalcands.begin(), mysortedrecoecalcands.end(), EgammaHLTEtSortCriterium());
0075   edm::Ref<reco::RecoEcalCandidateCollection> ref1, ref2;
0076 
0077   int n(0);
0078   for (unsigned int i = 0; i < mysortedrecoecalcands.size(); i++) {
0079     ref1 = mysortedrecoecalcands[i];
0080     if (ref1->et() >= etcut1_) {
0081       for (unsigned int j = i + 1; j < mysortedrecoecalcands.size(); j++) {
0082         ref2 = mysortedrecoecalcands[j];
0083         if (ref2->et() >= etcut2_) {
0084           // Acoplanarity
0085           double acop = fabs(ref1->phi() - ref2->phi());
0086           if (acop > M_PI)
0087             acop = 2 * M_PI - acop;
0088           acop = M_PI - acop;
0089           LogDebug("HLTEgammaDoubleEtPhiFilter") << " ... 1-2 acop= " << acop;
0090 
0091           if ((acop >= min_Acop_) && (acop <= max_Acop_)) {
0092             // Et balance
0093             double etbalance = fabs(ref1->et() - ref2->et());
0094             if ((etbalance >= min_EtBalance_) && (etbalance <= max_EtBalance_)) {
0095               filterproduct.addObject(TriggerCluster, ref1);
0096               filterproduct.addObject(TriggerCluster, ref2);
0097               n++;
0098             }
0099           }
0100         }
0101       }
0102     }
0103   }
0104 
0105   // filter decision
0106   bool accept(n >= npaircut_);
0107 
0108   return accept;
0109 }
0110 
0111 // declare this class as a framework plugin
0112 #include "FWCore/Framework/interface/MakerMacros.h"
0113 DEFINE_FWK_MODULE(HLTEgammaDoubleEtPhiFilter);