CMS 3D CMS Logo

HLTElectronPixelMatchFilter.cc
Go to the documentation of this file.
1 
10 
12 
16 
20 
23 
28 
32 
34 
36  candTag_ = iConfig.getParameter< edm::InputTag > ("candTag");
37  l1PixelSeedsTag_ = iConfig.getParameter< edm::InputTag > ("l1PixelSeedsTag");
38  npixelmatchcut_ = iConfig.getParameter< double > ("npixelmatchcut");
39  ncandcut_ = iConfig.getParameter< int > ("ncandcut");
40  l1EGTag_ = iConfig.getParameter< edm::InputTag > ("l1EGCand");
41 
42  candToken_ = consumes<trigger::TriggerFilterObjectWithRefs> (candTag_);
43  l1PixelSeedsToken_ = consumes<reco::ElectronSeedCollection> (l1PixelSeedsTag_);
44 
45  sPhi1B_ = iConfig.getParameter< double >("s_a_phi1B") ;
46  sPhi1I_ = iConfig.getParameter< double >("s_a_phi1I") ;
47  sPhi1F_ = iConfig.getParameter< double >("s_a_phi1F") ;
48  sPhi2B_ = iConfig.getParameter< double >("s_a_phi2B") ;
49  sPhi2I_ = iConfig.getParameter< double >("s_a_phi2I") ;
50  sPhi2F_ = iConfig.getParameter< double >("s_a_phi2F") ;
51  sZ2B_ = iConfig.getParameter< double >("s_a_zB" ) ;
52  sR2I_ = iConfig.getParameter< double >("s_a_rI" ) ;
53  sR2F_ = iConfig.getParameter< double >("s_a_rF" ) ;
54  s2BarrelThres_ = std::pow(std::atanh(iConfig.getParameter< double >("tanhSO10BarrelThres"))*10., 2);
55  s2InterThres_ = std::pow(std::atanh(iConfig.getParameter< double >("tanhSO10InterThres"))*10., 2);
56  s2ForwardThres_ = std::pow(std::atanh(iConfig.getParameter< double >("tanhSO10ForwardThres"))*10., 2);
57 
58  isPixelVeto_ = iConfig.getParameter< bool >("pixelVeto");
59  useS_ = iConfig.getParameter< bool >("useS" );
60 
61 }
62 
64 
65 float HLTElectronPixelMatchFilter::calDPhi1Sq(reco::ElectronSeedCollection::const_iterator seed, int charge)const
66 {
67  const float dPhi1Const = seed->subDet1()==1 ? seed->subDet2()==1 ? sPhi1B_ : sPhi1I_ : sPhi1F_;
68  float dPhi1 = charge<0 ? seed->dPhi1()/dPhi1Const : seed->dPhi1Pos()/dPhi1Const;
69  return dPhi1*dPhi1;
70 }
71 
72 float HLTElectronPixelMatchFilter::calDPhi2Sq(reco::ElectronSeedCollection::const_iterator seed, int charge)const
73 {
74  const float dPhi2Const = seed->subDet1()==1 ? seed->subDet2()==1 ? sPhi2B_ : sPhi2I_ : sPhi2F_;
75  float dPhi2 = charge <0 ? seed->dPhi2()/dPhi2Const : seed->dPhi2Pos()/dPhi2Const;
76  return dPhi2*dPhi2;
77 }
78 
79 
80 float HLTElectronPixelMatchFilter::calDZ2Sq(reco::ElectronSeedCollection::const_iterator seed, int charge)const
81 {
82  const float dRZ2Const = seed->subDet1()==1 ? seed->subDet2()==1 ? sZ2B_ : sR2I_ : sR2F_;
83  float dRZ2 = charge<0 ? seed->dRz2()/dRZ2Const : seed->dRz2Pos()/dRZ2Const;
84  return dRZ2*dRZ2;
85 }
86 
90  desc.add<edm::InputTag>("candTag", edm::InputTag("hltEgammaHcalIsolFilter"));
91  desc.add<edm::InputTag>("l1PixelSeedsTag", edm::InputTag("electronPixelSeeds"));
92  desc.add<double>("npixelmatchcut", 1.0);
93  desc.add<int>("ncandcut", 1);
94  desc.add<edm::InputTag>("l1EGCand", edm::InputTag("hltL1IsoRecoEcalCandidate"));
95  desc.add<double>("s_a_phi1B", 0.0069) ;
96  desc.add<double>("s_a_phi1I", 0.0088) ;
97  desc.add<double>("s_a_phi1F", 0.0076) ;
98  desc.add<double>("s_a_phi2B", 0.00037) ;
99  desc.add<double>("s_a_phi2I", 0.00070) ;
100  desc.add<double>("s_a_phi2F", 0.00906) ;
101  desc.add<double>("s_a_zB" , 0.012) ;
102  desc.add<double>("s_a_rI" , 0.027) ;
103  desc.add<double>("s_a_rF" , 0.040) ;
104  desc.add<double>("s2_threshold", 0);
105  desc.add<double>("tanhSO10BarrelThres", 0.35);
106  desc.add<double>("tanhSO10InterThres", 1);
107  desc.add<double>("tanhSO10ForwardThres", 1);
108  desc.add<bool> ("useS" , false);
109  desc.add<bool> ("pixelVeto", false);
110 
111  descriptions.add("hltElectronPixelMatchFilter",desc);
112 }
113 
115  // The filter object
116  using namespace trigger;
117  if (saveTags()) {
118  filterproduct.addCollectionTag(l1EGTag_);
119  }
120 
121  // Ref to Candidate object to be recorded in filter object
123 
125  iEvent.getByToken(candToken_,PrevFilterOutput);
126 
127  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > recoecalcands;
128  PrevFilterOutput->getObjects(TriggerCluster, recoecalcands);
129  if(recoecalcands.empty()) PrevFilterOutput->getObjects(TriggerPhoton,recoecalcands); //we dont know if its type trigger cluster or trigger photon
130 
131  //get hold of the pixel seed - supercluster association map
133  iEvent.getByToken(l1PixelSeedsToken_, l1PixelSeeds);
134 
135  // look at all egammas, check cuts and add to filter object
136  int n = 0;
137  for (auto & recoecalcand : recoecalcands) {
138 
139  ref = recoecalcand;
140  reco::SuperClusterRef recr2 = ref->superCluster();
141 
142  int nmatch = getNrOfMatches(l1PixelSeeds, recr2);
143 
144  if (!isPixelVeto_) {
145  if ( nmatch >= npixelmatchcut_) {
146  n++;
147  filterproduct.addObject(TriggerCluster, ref);
148  }
149  } else {
150  if ( nmatch == 0) {
151  n++;
152  filterproduct.addObject(TriggerCluster, ref);
153  }
154  }
155 
156  }//end of loop over candidates
157 
158  // filter decision
159  const bool accept(n>=ncandcut_);
160  return accept;
161 }
162 
164  reco::SuperClusterRef& candSCRef)const
165 {
166  int nrMatch=0;
167  for(auto seedIt = eleSeeds->begin(); seedIt != eleSeeds->end(); seedIt++){
168  edm::RefToBase<reco::CaloCluster> caloCluster = seedIt->caloCluster() ;
169  reco::SuperClusterRef scRef = caloCluster.castTo<reco::SuperClusterRef>() ;
170  if(&(*candSCRef) == &(*scRef)){
171  if(useS_){
172  float s2Neg = calDPhi1Sq(seedIt,-1) + calDPhi2Sq(seedIt,-1) + calDZ2Sq(seedIt,-1);
173  float s2Pos = calDPhi1Sq(seedIt,1) + calDPhi2Sq(seedIt,1) + calDZ2Sq(seedIt,1);
174 
175  const float s2Thres = seedIt->subDet1()==1 ? seedIt->subDet2()==1 ? s2BarrelThres_ : s2InterThres_ : s2ForwardThres_;
176  if(s2Neg<s2Thres || s2Pos<s2Thres) nrMatch++ ;
177  }
178  else nrMatch++;
179  }//end sc ref match
180  }//end loop over ele seeds
181  return nrMatch;
182 }
183 
184 // declare this class as a framework plugin
T getParameter(std::string const &) const
float calDPhi1Sq(reco::ElectronSeedCollection::const_iterator seed, int charge) const
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > candToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
HLTElectronPixelMatchFilter(const edm::ParameterSet &)
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
edm::EDGetTokenT< reco::ElectronSeedCollection > l1PixelSeedsToken_
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
~HLTElectronPixelMatchFilter() override
ParameterDescriptionBase * add(U const &iLabel, T const &value)
int getNrOfMatches(edm::Handle< reco::ElectronSeedCollection > &eleSeeds, reco::SuperClusterRef &candSCRef) const
float calDZ2Sq(reco::ElectronSeedCollection::const_iterator seed, int charge) const
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
REF castTo() const
Definition: RefToBase.h:289
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool saveTags() const
Definition: HLTFilter.h:45
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
float calDPhi2Sq(reco::ElectronSeedCollection::const_iterator seed, int charge) const