CMS 3D CMS Logo

PFClusterMatchedToPhotonsSelector.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: CommonTools/RecoAlgos
4 // Class: PFClusterMatchedToPhotonsSelector
5 //
13 //
14 // Original Author: RCLSA
15 // Created: Wed, 22 Mar 2017 18:01:40 GMT
16 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 #include <iostream>
23 
24 // user include files
30 
40 
42 bool sortByKey(const EEPSPair& a, const EEPSPair& b) {
43  return a.first < b.first;
44 }
45 
47 public:
50 
51 private:
52  void produce(edm::Event&, const edm::EventSetup&) override;
53 
58  double matchMaxDR2_;
60 
61  double volumeZ_EB_;
63  double volumeZ_EE_;
64 };
65 
67  //now do what ever initialization is needed
68  particleFlowClusterECALToken_ = consumes<reco::PFClusterCollection>(iConfig.getParameter<edm::InputTag>("pfClustersTag"));
69  associationToken_ = consumes<reco::PFCluster::EEtoPSAssociation>(iConfig.getParameter<edm::InputTag>("pfClustersTag"));
70  trackingParticleToken_ = consumes<TrackingParticleCollection>(iConfig.getParameter<edm::InputTag>("trackingParticleTag"));
71  genParticleToken_ = consumes<reco::GenParticleCollection>(iConfig.getParameter<edm::InputTag>("genParticleTag"));
72 
73  matchMaxDR2_ = iConfig.getParameter<double>("maxDR2");
74  matchMaxDEDR2_ = iConfig.getParameter<double>("maxDEDR2");
75  volumeZ_EB_ = iConfig.getParameter<double>("volumeZ_EB");
76  volumeRadius_EB_ = iConfig.getParameter<double>("volumeRadius_EB");
77  volumeZ_EE_ = iConfig.getParameter<double>("volumeZ_EE");
78 
79  produces<reco::PFClusterCollection>();
80  produces<reco::PFCluster::EEtoPSAssociation>();
81  produces<edm::ValueMap<reco::GenParticleRef> >();
82 }
83 
84 
87  desc.add<edm::InputTag>("pfClustersTag", edm::InputTag("particleFlowClusterECAL"));
88  desc.add<edm::InputTag>("trackingParticleTag", edm::InputTag("mix", "MergedTrackTruth"));
89  desc.add<edm::InputTag>("genParticleTag", edm::InputTag("genParticles"));
90  desc.add<double>("maxDR2", 0.1*0.1);
91  desc.add<double>("maxDEDR2", 0.5*0.5);
92  desc.add<double>("volumeZ_EB", 304.5);
93  desc.add<double>("volumeRadius_EB", 123.8);
94  desc.add<double>("volumeZ_EE", 317.0);
95  descriptions.add("pfClusterMatchedToPhotonsSelector", desc);
96 
97 }
98 
100 
101  edm::Handle<reco::PFClusterCollection> particleFlowClusterECALHandle_;
103  edm::Handle<TrackingParticleCollection> trackingParticleHandle_;
104  edm::Handle<reco::GenParticleCollection> genParticleHandle_;
105  iEvent.getByToken(particleFlowClusterECALToken_, particleFlowClusterECALHandle_);
106  iEvent.getByToken(trackingParticleToken_, trackingParticleHandle_);
107  iEvent.getByToken(genParticleToken_, genParticleHandle_);
108  iEvent.getByToken(associationToken_, associationHandle_);
109 
110  std::unique_ptr<reco::PFClusterCollection> out = std::make_unique<reco::PFClusterCollection>();
111  std::unique_ptr<reco::PFCluster::EEtoPSAssociation> association_out = std::make_unique<reco::PFCluster::EEtoPSAssociation>();
112  std::unique_ptr<edm::ValueMap<reco::GenParticleRef> > genmatching_out = std::make_unique<edm::ValueMap<reco::GenParticleRef> >();
113 
114  std::vector<reco::GenParticleRef> genmatching;
115 
116  size_t iN(0);
117  for (size_t iP = 0; iP < particleFlowClusterECALHandle_->size(); iP++) {
118 
119  auto&& pfCluster = particleFlowClusterECALHandle_->at(iP);
120  bool isMatched = false;
122 
123  // Preselect PFclusters
124  if (pfCluster.energy() <= 0) {
125  continue;
126  }
127 
128  for (auto&& trackingParticle : *trackingParticleHandle_) {
129  if (trackingParticle.pdgId() != 22) continue;
130  if (trackingParticle.status() != 1) continue;
131  matchedKey = trackingParticle.genParticles().at(0).key();
132  float dR2 = reco::deltaR2(trackingParticle, pfCluster.position());
133  if (dR2 > matchMaxDR2_) continue;
134  float dE = 1. - trackingParticle.genParticles().at(0)->energy()/pfCluster.energy();
135  if ((dR2 + dE*dE) > matchMaxDEDR2_) continue;
136 
137  bool isConversion = false;
138  for (auto&& vertRef : trackingParticle.decayVertices()) {
139  if (vertRef->position().rho() > volumeRadius_EB_ && std::abs(vertRef->position().z()) < volumeZ_EB_) continue;
140  if (std::abs(vertRef->position().z()) > volumeZ_EE_) continue;
141 
142  for(auto&& tpRef: vertRef->daughterTracks()) {
143  if(std::abs(tpRef->pdgId()) == 11) isConversion = true;
144  break;
145  }
146  if (isConversion) break;
147  }
148  if (isConversion) continue;
149 
150  isMatched = true;
151  break;
152  }
153 
154  if (isMatched) {
155  out->push_back(pfCluster);
156  for (size_t i=0; i<associationHandle_.product()->size(); i++) {
157  if (associationHandle_.product()->at(i).first == iP) {
158  association_out->push_back(std::make_pair(iN, associationHandle_.product()->at(i).second));
159  }
160  }
161  genmatching.push_back(edm::Ref<reco::GenParticleCollection>(genParticleHandle_,matchedKey));
162  }
163  }
164 
165  std::sort(association_out->begin(),association_out->end(),sortByKey);
166  edm::OrphanHandle<reco::PFClusterCollection> pfClusterHandle= iEvent.put(std::move(out));
167  iEvent.put(std::move(association_out));
168 
169  edm::ValueMap<reco::GenParticleRef>::Filler mapFiller(*genmatching_out);
170  mapFiller.insert(pfClusterHandle, genmatching.begin(), genmatching.end());
171  mapFiller.fill();
172  iEvent.put(std::move(genmatching_out));
173 }
174 
175 //define this as a plug-in
T getParameter(std::string const &) const
std::remove_cv< typename std::remove_reference< argument_type >::type >::type key_type
Definition: Ref.h:169
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:136
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< reco::GenParticleCollection > genParticleToken_
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
int iEvent
Definition: GenABIO.cc:230
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< TrackingParticleCollection > trackingParticleToken_
bool isMatched(TrackingRecHit const &hit)
edm::EDGetTokenT< reco::PFClusterCollection > particleFlowClusterECALToken_
static void fillDescriptions(edm::ConfigurationDescriptions &)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
T const * product() const
Definition: Handle.h:81
reco::PFCluster::EEtoPSAssociation::value_type EEPSPair
bool sortByKey(const EEPSPair &a, const EEPSPair &b)
double b
Definition: hdecay.h:120
void add(std::string const &label, ParameterSetDescription const &psetDescription)
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
double a
Definition: hdecay.h:121
void produce(edm::Event &, const edm::EventSetup &) override
edm::EDGetTokenT< reco::PFCluster::EEtoPSAssociation > associationToken_
PFClusterMatchedToPhotonsSelector(const edm::ParameterSet &)
def move(src, dest)
Definition: eostools.py:510