CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GEDPhotonCoreProducer.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <vector>
3 #include <memory>
11 // Framework
31 #include <Math/VectorUtil.h>
32 #include <vector>
33 #include "TLorentzVector.h"
34 #include "TMath.h"
35 
37  conf_(config)
38 
39 {
40 
41  // use onfiguration file to setup input/output collection names
43  consumes<reco::PFCandidateCollection>(conf_.getParameter<edm::InputTag>("pfEgammaCandidates"));
45  consumes<reco::ElectronSeedCollection>(conf_.getParameter<edm::InputTag>("pixelSeedProducer"));
46 
47  GEDPhotonCoreCollection_ = conf_.getParameter<std::string>("gedPhotonCoreCollection");
48 
49 
50 
51  // Register the product
52  produces<reco::PhotonCoreCollection>(GEDPhotonCoreCollection_);
53 
54 }
55 
57 
58 
59 
60 
61 void GEDPhotonCoreProducer::produce(edm::Event &theEvent, const edm::EventSetup& theEventSetup) {
62 
63 
64  using namespace edm;
65  // nEvt_++;
66 
67  reco::PhotonCoreCollection outputPhotonCoreCollection;
68  std::auto_ptr< reco::PhotonCoreCollection > outputPhotonCoreCollection_p(new reco::PhotonCoreCollection);
69 
70  // Get the PF refined cluster collection
71  Handle<reco::PFCandidateCollection> pfCandidateHandle;
72  theEvent.getByToken(pfEgammaCandidates_,pfCandidateHandle);
73  if (!pfCandidateHandle.isValid()) {
74  edm::LogError("GEDPhotonCoreProducer")
75  << "Error! Can't get the pfEgammaCandidates";
76  }
77 
78 
79  // Get ElectronPixelSeeds
80  validPixelSeeds_=true;
83  theEvent.getByToken(pixelSeedProducer_, pixelSeedHandle);
84  if (!pixelSeedHandle.isValid()) {
85  validPixelSeeds_=false;
86  }
87 
88 
89 
90  // std::cout << " GEDPhotonCoreProducer::produce input PFcandidate size " << pfCandidateHandle->size() << std::endl;
91 
92 
93  // Loop over PF candidates and get only photons
94  reco::ElectronSeedCollection::const_iterator pixelSeedItr;
95  for(unsigned int lCand=0; lCand < pfCandidateHandle->size(); lCand++) {
96  reco::PFCandidateRef candRef (reco::PFCandidateRef(pfCandidateHandle,lCand));
97 
98  // Retrieve stuff from the pfPhoton
99  reco::PFCandidateEGammaExtraRef pfPhoRef = candRef->egammaExtraRef();
100  reco::SuperClusterRef refinedSC= pfPhoRef->superClusterRef();
101  reco::SuperClusterRef boxSC= pfPhoRef->superClusterPFECALRef();
102  const reco::ConversionRefVector & doubleLegConv = pfPhoRef->conversionRef();
103  const reco::ConversionRefVector & singleLegConv = pfPhoRef->singleLegConversionRef();
104  reco::CaloClusterPtr refinedSCPtr= edm::refToPtr(refinedSC);
105 
106  // std::cout << "newCandidate doubleLegConv="<<doubleLegConv.size()<< std::endl;
107  //std::cout << "newCandidate singleLegConv="<< pfPhoRef->singleLegConvTrackRef().size()<< std::endl;
108 
110  reco::PhotonCore newCandidate;
111  newCandidate.setPFlowPhoton(true);
112  newCandidate.setStandardPhoton(false);
113  newCandidate.setSuperCluster(refinedSC);
114  newCandidate.setParentSuperCluster(boxSC);
115  // fill conversion infos
116 
117 
118  for(unsigned int lConv=0; lConv < doubleLegConv.size(); lConv++) {
119  newCandidate.addConversion(doubleLegConv[lConv]);
120  }
121 
122  for(unsigned int lConv=0; lConv < singleLegConv.size(); lConv++) {
123  newCandidate.addOneLegConversion(singleLegConv[lConv]);
124  }
125 
126  // std::cout << "newCandidate pf refined SC energy="<< newCandidate.superCluster()->energy()<<std::endl;
127  //std::cout << "newCandidate pf SC energy="<< newCandidate.parentSuperCluster()->energy()<<std::endl;
128  //std::cout << "newCandidate nconv2leg="<<newCandidate.conversions().size()<< std::endl;
129 
130  if ( validPixelSeeds_) {
131  for( unsigned int icp = 0; icp < pixelSeedHandle->size(); icp++) {
132  reco::ElectronSeedRef cpRef(pixelSeedHandle,icp);
133  if ( boxSC.isNonnull() && boxSC.id() == cpRef->caloCluster().id() && boxSC.key() == cpRef->caloCluster().key() ) {
134  newCandidate.addElectronPixelSeed(cpRef);
135  }
136  }
137  }
138 
139  outputPhotonCoreCollection.push_back(newCandidate);
140  }
141 
142  // put the product in the event
143  // edm::LogInfo("GEDPhotonCoreProducer") << " Put in the event " << iSC << " Photon Candidates \n";
144  outputPhotonCoreCollection_p->assign(outputPhotonCoreCollection.begin(),outputPhotonCoreCollection.end());
145  theEvent.put( outputPhotonCoreCollection_p, GEDPhotonCoreCollection_);
146 
147 
148 
149 
150 }
T getParameter(std::string const &) const
void setSuperCluster(const reco::SuperClusterRef &r)
set reference to SuperCluster
Definition: PhotonCore.h:50
void addOneLegConversion(const reco::ConversionRef &r)
add single ConversionRef to the vector of Refs
Definition: PhotonCore.h:56
void addConversion(const reco::ConversionRef &r)
add single ConversionRef to the vector of Refs
Definition: PhotonCore.h:54
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
Ptr< typename C::value_type > refToPtr(Ref< C, typename C::value_type, refhelper::FindUsingAdvance< C, typename C::value_type > > const &ref)
Definition: RefToPtr.h:18
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
key_type key() const
Accessor for product key.
Definition: Ref.h:264
edm::EDGetTokenT< reco::ElectronSeedCollection > pixelSeedProducer_
ProductID id() const
Accessor for product ID.
Definition: Ref.h:258
void setParentSuperCluster(const reco::SuperClusterRef &r)
set reference to PFlow SuperCluster
Definition: PhotonCore.h:52
void setStandardPhoton(const bool prov)
Definition: PhotonCore.h:61
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
std::vector< ElectronSeed > ElectronSeedCollection
collection of ElectronSeed objects
void setPFlowPhoton(const bool prov)
set the provenance
Definition: PhotonCore.h:60
std::vector< PhotonCore > PhotonCoreCollection
collectin of PhotonCore objects
Definition: PhotonCoreFwd.h:9
GEDPhotonCoreProducer(const edm::ParameterSet &ps)
virtual void produce(edm::Event &evt, const edm::EventSetup &es)
size_type size() const
Size of the RefVector.
Definition: RefVector.h:99
void addElectronPixelSeed(const reco::ElectronSeedRef &r)
set electron pixel seed ref
Definition: PhotonCore.h:58
edm::EDGetTokenT< reco::PFCandidateCollection > pfEgammaCandidates_