Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /** \class HLTElectronPixelMatchFilter
0002  *
0003  * $Id: HLTElectronPixelMatchFilter.cc,v 1.15 2012/03/06 10:13:59 sharper Exp $
0004  *
0005  *  \author Aidan Randle-Conde (ULB)
0006  *
0007  */
0008 
0009 #include "HLTElectronPixelMatchFilter.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/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0019 #include "DataFormats/Common/interface/AssociationMap.h"
0020 
0021 #include "DataFormats/CaloRecHit/interface/CaloCluster.h"
0022 #include "DataFormats/CaloRecHit/interface/CaloClusterFwd.h"
0023 
0024 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
0025 #include "DataFormats/EgammaReco/interface/SuperClusterFwd.h"
0026 #include "DataFormats/EgammaReco/interface/ElectronSeed.h"
0027 #include "DataFormats/EgammaReco/interface/ElectronSeedFwd.h"
0028 
0029 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
0030 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidateFwd.h"
0031 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
0032 
0033 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
0034 
0035 HLTElectronPixelMatchFilter::HLTElectronPixelMatchFilter(const edm::ParameterSet& iConfig) : HLTFilter(iConfig) {
0036   candTag_ = iConfig.getParameter<edm::InputTag>("candTag");
0037   l1PixelSeedsTag_ = iConfig.getParameter<edm::InputTag>("l1PixelSeedsTag");
0038   npixelmatchcut_ = iConfig.getParameter<double>("npixelmatchcut");
0039   ncandcut_ = iConfig.getParameter<int>("ncandcut");
0040   l1EGTag_ = iConfig.getParameter<edm::InputTag>("l1EGCand");
0041 
0042   candToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(candTag_);
0043   l1PixelSeedsToken_ = consumes<reco::ElectronSeedCollection>(l1PixelSeedsTag_);
0044 
0045   sPhi1B_ = iConfig.getParameter<double>("s_a_phi1B");
0046   sPhi1I_ = iConfig.getParameter<double>("s_a_phi1I");
0047   sPhi1F_ = iConfig.getParameter<double>("s_a_phi1F");
0048   sPhi2B_ = iConfig.getParameter<double>("s_a_phi2B");
0049   sPhi2I_ = iConfig.getParameter<double>("s_a_phi2I");
0050   sPhi2F_ = iConfig.getParameter<double>("s_a_phi2F");
0051   sZ2B_ = iConfig.getParameter<double>("s_a_zB");
0052   sR2I_ = iConfig.getParameter<double>("s_a_rI");
0053   sR2F_ = iConfig.getParameter<double>("s_a_rF");
0054   s2BarrelThres_ = std::pow(std::atanh(iConfig.getParameter<double>("tanhSO10BarrelThres")) * 10., 2);
0055   s2InterThres_ = std::pow(std::atanh(iConfig.getParameter<double>("tanhSO10InterThres")) * 10., 2);
0056   s2ForwardThres_ = std::pow(std::atanh(iConfig.getParameter<double>("tanhSO10ForwardThres")) * 10., 2);
0057 
0058   isPixelVeto_ = iConfig.getParameter<bool>("pixelVeto");
0059   useS_ = iConfig.getParameter<bool>("useS");
0060 }
0061 
0062 HLTElectronPixelMatchFilter::~HLTElectronPixelMatchFilter() = default;
0063 
0064 float HLTElectronPixelMatchFilter::calDPhi1Sq(reco::ElectronSeedCollection::const_iterator seed, int charge) const {
0065   const float dPhi1Const = seed->subDet(0) == 1 ? seed->subDet(1) == 1 ? sPhi1B_ : sPhi1I_ : sPhi1F_;
0066   float dPhi1 = charge < 0 ? seed->dPhiNeg(0) / dPhi1Const : seed->dPhiPos(0) / dPhi1Const;
0067   return dPhi1 * dPhi1;
0068 }
0069 
0070 float HLTElectronPixelMatchFilter::calDPhi2Sq(reco::ElectronSeedCollection::const_iterator seed, int charge) const {
0071   const float dPhi2Const = seed->subDet(0) == 1 ? seed->subDet(1) == 1 ? sPhi2B_ : sPhi2I_ : sPhi2F_;
0072   float dPhi2 = charge < 0 ? seed->dPhiNeg(1) / dPhi2Const : seed->dPhiPos(1) / dPhi2Const;
0073   return dPhi2 * dPhi2;
0074 }
0075 
0076 float HLTElectronPixelMatchFilter::calDZ2Sq(reco::ElectronSeedCollection::const_iterator seed, int charge) const {
0077   const float dRZ2Const = seed->subDet(0) == 1 ? seed->subDet(1) == 1 ? sZ2B_ : sR2I_ : sR2F_;
0078   float dRZ2 = charge < 0 ? seed->dRZNeg(1) / dRZ2Const : seed->dRZPos(1) / dRZ2Const;
0079   return dRZ2 * dRZ2;
0080 }
0081 
0082 void HLTElectronPixelMatchFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0083   edm::ParameterSetDescription desc;
0084   makeHLTFilterDescription(desc);
0085   desc.add<edm::InputTag>("candTag", edm::InputTag("hltEgammaHcalIsolFilter"));
0086   desc.add<edm::InputTag>("l1PixelSeedsTag", edm::InputTag("electronPixelSeeds"));
0087   desc.add<double>("npixelmatchcut", 1.0);
0088   desc.add<int>("ncandcut", 1);
0089   desc.add<edm::InputTag>("l1EGCand", edm::InputTag("hltL1IsoRecoEcalCandidate"));
0090   desc.add<double>("s_a_phi1B", 0.0069);
0091   desc.add<double>("s_a_phi1I", 0.0088);
0092   desc.add<double>("s_a_phi1F", 0.0076);
0093   desc.add<double>("s_a_phi2B", 0.00037);
0094   desc.add<double>("s_a_phi2I", 0.00070);
0095   desc.add<double>("s_a_phi2F", 0.00906);
0096   desc.add<double>("s_a_zB", 0.012);
0097   desc.add<double>("s_a_rI", 0.027);
0098   desc.add<double>("s_a_rF", 0.040);
0099   desc.add<double>("s2_threshold", 0);
0100   desc.add<double>("tanhSO10BarrelThres", 0.35);
0101   desc.add<double>("tanhSO10InterThres", 1);
0102   desc.add<double>("tanhSO10ForwardThres", 1);
0103   desc.add<bool>("useS", false);
0104   desc.add<bool>("pixelVeto", false);
0105 
0106   descriptions.add("hltElectronPixelMatchFilter", desc);
0107 }
0108 
0109 bool HLTElectronPixelMatchFilter::hltFilter(edm::Event& iEvent,
0110                                             const edm::EventSetup& iSetup,
0111                                             trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0112   // The filter object
0113   using namespace trigger;
0114   if (saveTags()) {
0115     filterproduct.addCollectionTag(l1EGTag_);
0116   }
0117 
0118   // Ref to Candidate object to be recorded in filter object
0119   edm::Ref<reco::RecoEcalCandidateCollection> ref;
0120 
0121   edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
0122   iEvent.getByToken(candToken_, PrevFilterOutput);
0123 
0124   std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > recoecalcands;
0125   PrevFilterOutput->getObjects(TriggerCluster, recoecalcands);
0126   if (recoecalcands.empty())
0127     PrevFilterOutput->getObjects(TriggerPhoton,
0128                                  recoecalcands);  //we dont know if its type trigger cluster or trigger photon
0129 
0130   //get hold of the pixel seed - supercluster association map
0131   edm::Handle<reco::ElectronSeedCollection> l1PixelSeeds;
0132   iEvent.getByToken(l1PixelSeedsToken_, l1PixelSeeds);
0133 
0134   // look at all egammas,  check cuts and add to filter object
0135   int n = 0;
0136   for (auto& recoecalcand : recoecalcands) {
0137     ref = recoecalcand;
0138     reco::SuperClusterRef recr2 = ref->superCluster();
0139 
0140     int nmatch = getNrOfMatches(l1PixelSeeds, recr2);
0141 
0142     if (!isPixelVeto_) {
0143       if (nmatch >= npixelmatchcut_) {
0144         n++;
0145         filterproduct.addObject(TriggerCluster, ref);
0146       }
0147     } else {
0148       if (nmatch == 0) {
0149         n++;
0150         filterproduct.addObject(TriggerCluster, ref);
0151       }
0152     }
0153 
0154   }  //end of loop over candidates
0155 
0156   // filter decision
0157   const bool accept(n >= ncandcut_);
0158   return accept;
0159 }
0160 
0161 int HLTElectronPixelMatchFilter::getNrOfMatches(edm::Handle<reco::ElectronSeedCollection>& eleSeeds,
0162                                                 reco::SuperClusterRef& candSCRef) const {
0163   int nrMatch = 0;
0164   for (auto seedIt = eleSeeds->begin(); seedIt != eleSeeds->end(); seedIt++) {
0165     edm::RefToBase<reco::CaloCluster> caloCluster = seedIt->caloCluster();
0166     reco::SuperClusterRef scRef = caloCluster.castTo<reco::SuperClusterRef>();
0167     if (&(*candSCRef) == &(*scRef)) {
0168       if (useS_) {
0169         float s2Neg = calDPhi1Sq(seedIt, -1) + calDPhi2Sq(seedIt, -1) + calDZ2Sq(seedIt, -1);
0170         float s2Pos = calDPhi1Sq(seedIt, 1) + calDPhi2Sq(seedIt, 1) + calDZ2Sq(seedIt, 1);
0171 
0172         const float s2Thres =
0173             seedIt->subDet(0) == 1 ? seedIt->subDet(1) == 1 ? s2BarrelThres_ : s2InterThres_ : s2ForwardThres_;
0174         if (s2Neg < s2Thres || s2Pos < s2Thres)
0175           nrMatch++;
0176       } else
0177         nrMatch++;
0178     }  //end sc ref match
0179   }    //end loop over ele seeds
0180   return nrMatch;
0181 }
0182 
0183 // declare this class as a framework plugin
0184 #include "FWCore/Framework/interface/MakerMacros.h"
0185 DEFINE_FWK_MODULE(HLTElectronPixelMatchFilter);