CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EgammaHLTPFChargedIsolationProducer.cc
Go to the documentation of this file.
1 
9 
15 
17 
19 
21  pfCandidateProducer_(consumes<reco::PFCandidateCollection>(config.getParameter<edm::InputTag>("pfCandidatesProducer"))),
22  beamSpotProducer_ (consumes<reco::BeamSpot>(config.getParameter<edm::InputTag>("beamSpotProducer"))),
23  useGsfTrack_ (config.getParameter<bool>("useGsfTrack")),
24  useSCRefs_ (config.getParameter<bool>("useSCRefs")),
25  drMax_ (config.getParameter<double>("drMax")),
26  drVetoBarrel_ (config.getParameter<double>("drVetoBarrel")),
27  drVetoEndcap_ (config.getParameter<double>("drVetoEndcap")),
28  ptMin_ (config.getParameter<double>("ptMin")),
29  dzMax_ (config.getParameter<double>("dzMax")),
30  dxyMax_ (config.getParameter<double>("dxyMax")),
31  pfToUse_ (config.getParameter<int>("pfCandidateType")) {
32 
33  if(useSCRefs_) {
34  recoEcalCandidateProducer_ = consumes<reco::RecoEcalCandidateCollection>(config.getParameter<edm::InputTag>("recoEcalCandidateProducer"));
35  produces < reco::RecoEcalCandidateIsolationMap >();
36  } else {
37  electronProducer_ = consumes<reco::ElectronCollection>(config.getParameter<edm::InputTag>("electronProducer"));
38  produces < reco::ElectronIsolationMap >();
39  }
40 }
41 
44  desc.add<edm::InputTag>("electronProducer", edm::InputTag("hltEle27WP80PixelMatchElectronsL1SeededPF"));
45  desc.add<edm::InputTag>("recoEcalCandidateProducer", edm::InputTag("hltL1SeededRecoEcalCandidatePF"));
46  desc.add<edm::InputTag>("pfCandidatesProducer", edm::InputTag("hltParticleFlowReg"));
47  desc.add<edm::InputTag>("beamSpotProducer", edm::InputTag("hltOnlineBeamSpot"));
48  desc.add<bool>("useGsfTrack", false);
49  desc.add<bool>("useSCRefs", false);
50  desc.add<double>("drMax", 0.3);
51  desc.add<double>("drVetoBarrel", 0.02);
52  desc.add<double>("drVetoEndcap", 0.02);
53  desc.add<double>("ptMin", 0.0);
54  desc.add<double>("dzMax", 0.2);
55  desc.add<double>("dxyMax", 0.1);
56  desc.add<int>("pfCandidateType", 1);
57  descriptions.add(("hltEgammaHLTPFChargedIsolationProducer"), desc);
58 }
59 
61 
65  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
66 
67  iEvent.getByToken(pfCandidateProducer_, pfHandle);
68  const reco::PFCandidateCollection* forIsolation = pfHandle.product();
69 
72 
73  if(useSCRefs_) {
74 
75  iEvent.getByToken(recoEcalCandidateProducer_, recoEcalCandHandle);
76  iEvent.getByToken(beamSpotProducer_, recoBeamSpotHandle);
77  const reco::BeamSpot::Point& beamSpotPosition = recoBeamSpotHandle->position();
78 
79  float dRveto = -1;
80 
81  for(unsigned int iReco=0; iReco<recoEcalCandHandle->size(); iReco++) {
82  reco::RecoEcalCandidateRef candRef(recoEcalCandHandle, iReco);
83 
84  if (fabs(candRef->eta())<1.479)
85  dRveto = drVetoBarrel_;
86  else
87  dRveto = drVetoEndcap_;
88 
89  // Shift the RecoEcalCandidate direction vector according to the vertex
90  math::XYZVector candDirectionWrtVtx(candRef->superCluster()->x() - beamSpotPosition.x(),
91  candRef->superCluster()->y() - beamSpotPosition.y(),
92  candRef->superCluster()->z() - beamSpotPosition.z());
93 
94  float sum = 0;
95 
96  // Loop over the PFCandidates
97  for(unsigned i=0; i<forIsolation->size(); i++) {
98  const reco::PFCandidate& pfc = (*forIsolation)[i];
99 
100  //require that the PFCandidate is a charged hadron
101  if (pfc.particleId() == pfToUse_) {
102 
103  if(pfc.pt() < ptMin_) continue;
104 
105  float dz = fabs(pfc.trackRef()->dz(beamSpotPosition));
106  if(dz > dzMax_) continue;
107 
108  float dxy = fabs(pfc.trackRef()->dxy(beamSpotPosition));
109  if(fabs(dxy) > dxyMax_) continue;
110 
111  float dR = deltaR(candDirectionWrtVtx.Eta(), candDirectionWrtVtx.Phi(), pfc.momentum().Eta(), pfc.momentum().Phi());
112  if(dR > drMax_ || dR < dRveto) continue;
113 
114  sum += pfc.pt();
115  }
116  }
117 
118  recoEcalCandMap.insert(candRef, sum);
119  }
120 
121  } else {
122 
123  iEvent.getByToken(electronProducer_,electronHandle);
124 
125  float dRveto = -1;
126 
127  for(unsigned int iEl=0; iEl<electronHandle->size(); iEl++) {
128  reco::ElectronRef eleRef(electronHandle, iEl);
129  //const reco::Track* eleTrk = useGsfTrack_ ? &*eleRef->gsfTrack() : &*eleRef->track();
130 
131  if (fabs(eleRef->eta())<1.479)
132  dRveto = drVetoBarrel_;
133  else
134  dRveto = drVetoEndcap_;
135 
136  float sum = 0;
137 
138  // Loop over the PFCandidates
139  for(unsigned i=0; i<forIsolation->size(); i++) {
140  const reco::PFCandidate& pfc = (*forIsolation)[i];
141 
142  //require that the PFCandidate is a charged hadron
143  if (pfc.particleId() == pfToUse_) {
144 
145  if(pfc.pt() < ptMin_) continue;
146 
147  float dz = fabs(pfc.trackRef()->dz(eleRef->vertex()));
148  if(dz > dzMax_) continue;
149 
150  float dxy = fabs(pfc.trackRef()->dxy(eleRef->vertex()));
151  if(fabs(dxy) > dxyMax_) continue;
152 
153  float dR = deltaR(eleRef->eta(), eleRef->phi(), pfc.momentum().Eta(), pfc.momentum().Phi());
154  if(dR > drMax_ || dR < dRveto) continue;
155 
156  sum += pfc.pt();
157  }
158  }
159 
160  eleMap.insert(eleRef, sum);
161  }
162 
163  }
164 
165  if(useSCRefs_){
166  std::auto_ptr<reco::RecoEcalCandidateIsolationMap> mapForEvent(new reco::RecoEcalCandidateIsolationMap(recoEcalCandMap));
167  iEvent.put(mapForEvent);
168  }else{
169  std::auto_ptr<reco::ElectronIsolationMap> mapForEvent(new reco::ElectronIsolationMap(eleMap));
170  iEvent.put(mapForEvent);
171  }
172 }
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
const edm::EDGetTokenT< reco::PFCandidateCollection > pfCandidateProducer_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
virtual Vector momentum() const
spatial momentum vector
math::XYZPoint Point
point in the space
Definition: BeamSpot.h:29
virtual double pt() const
transverse momentum
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:433
edm::EDGetTokenT< reco::RecoEcalCandidateCollection > recoEcalCandidateProducer_
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
const edm::EDGetTokenT< reco::BeamSpot > beamSpotProducer_
edm::EDGetTokenT< reco::ElectronCollection > electronProducer_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
T const * product() const
Definition: Handle.h:81
void insert(const key_type &k, const data_type &v)
insert an association
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
virtual ParticleType particleId() const
Definition: PFCandidate.h:373
void produce(edm::Event &, const edm::EventSetup &) override