CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FSPFProducer.cc
Go to the documentation of this file.
4 
8 
9 using namespace std;
10 using namespace edm;
11 using namespace reco;
12 
13 
15 
16  labelPFCandidateCollection_ = iConfig.getParameter < edm::InputTag > ("pfCandidates");
17 
18  pfPatchInHF = iConfig.getParameter<bool>("pfPatchInHF");
19  EM_HF_ScaleFactor = iConfig.getParameter< std::vector <double> >("EM_HF_ScaleFactor");
20  HF_Ratio = iConfig.getParameter<double>("HF_Ratio");
21  par1 = iConfig.getParameter<double>("par1");
22  par2 = iConfig.getParameter<double>("par2");
23  barrel_th = iConfig.getParameter<double>("barrel_th");
24  endcap_th = iConfig.getParameter<double>("endcap_th");
25  middle_th = iConfig.getParameter<double>("middle_th");
26  // register products
27  produces<reco::PFCandidateCollection>();
28 
29 }
30 
32 
33 void
35  const EventSetup& iSetup) {
36 
38  iEvent.getByLabel (labelPFCandidateCollection_, pfCandidates);
39 
40  auto_ptr< reco::PFCandidateCollection > pOutputCandidateCollection(new PFCandidateCollection);
41 
42  /*
43  LogDebug("FSPFProducer")<<"START event: "
44  <<iEvent.id().event()
45  <<" in run "<<iEvent.id().run()<<endl;
46  */
47 
48  double theNeutralFraction, px, py, pz, en;
49  int n_hadron_HF = 0;
50  int n_em_HF = 0;
51  int vEta=-1;
52  reco::PFCandidateCollection::const_iterator itCand = pfCandidates->begin();
53  reco::PFCandidateCollection::const_iterator itCandEnd = pfCandidates->end();
54  for( ; itCand != itCandEnd; itCand++) {
55 
56  // First part: create fake neutral hadrons as a function of charged hadrons and add them to the new collection
57  if(itCand->particleId() == reco::PFCandidate::h){
58  theNeutralFraction = par1 - par2*itCand->energy();
59  if(theNeutralFraction > 0.){
60  px = theNeutralFraction*itCand->px();
61  py = theNeutralFraction*itCand->py();
62  pz = theNeutralFraction*itCand->pz();
63  en = sqrt(px*px + py*py + pz*pz);
64  if (en > energy_threshold(itCand->eta())) {
65  // create a PFCandidate and add it to the particles Collection
66  math::XYZTLorentzVector momentum(px,py,pz,en);
67  reco::PFCandidate FakeNeutralHadron(0, momentum, reco::PFCandidate::h0);
68  pOutputCandidateCollection->push_back(FakeNeutralHadron);
69  }
70  }
71  }
72 
73  // Second part: deal with HF, and put every candidate of the old collection in the new collection
74  if(itCand->particleId() == reco::PFCandidate::egamma_HF){
75  if (pfPatchInHF) {
76  n_em_HF++;
77 
78  if(fabs(itCand->eta())< 4.) vEta = 0;
79  else if(fabs(itCand->eta())<= 5.) vEta = 1;
80 
81  if (vEta==0 || vEta==1) {
82  // copy these PFCandidates after the momentum rescaling
83  px = EM_HF_ScaleFactor[vEta]*itCand->px();
84  py = EM_HF_ScaleFactor[vEta]*itCand->py();
85  pz = EM_HF_ScaleFactor[vEta]*itCand->pz();
86  en = sqrt(px*px + py*py + pz*pz);
87  math::XYZTLorentzVector momentum(px,py,pz,en);
88  reco::PFCandidate EMHF(itCand->charge(), momentum, reco::PFCandidate::egamma_HF);
89  if(en>0.) pOutputCandidateCollection->push_back(EMHF);
90  }
91  } else pOutputCandidateCollection->push_back(*itCand);
92  }
93 
94  else if(itCand->particleId() == reco::PFCandidate::h_HF){
95  if (pfPatchInHF) {
96  // copy these PFCandidates to the new particles Collection only if hadron candidates are currently more than EM candidates
97  n_hadron_HF++;
98  if(n_em_HF < (n_hadron_HF*HF_Ratio)) pOutputCandidateCollection->push_back(*itCand);
99  } else pOutputCandidateCollection->push_back(*itCand);
100  }
101 
102  else pOutputCandidateCollection->push_back(*itCand);
103 
104  }
105  iEvent.put(pOutputCandidateCollection);
106 }
107 
109  if (eta<0) eta = -eta;
110  if (eta < 1.6) return barrel_th;
111  else if (eta < 1.8) return middle_th;
112  else return endcap_th;
113 }
114 
115 
116 
T getParameter(std::string const &) const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual void produce(edm::Event &, const edm::EventSetup &) override
Definition: FSPFProducer.cc:34
T eta() const
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
T sqrt(T t)
Definition: SSEVec.h:48
FSPFProducer(const edm::ParameterSet &)
Definition: FSPFProducer.cc:14
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:38
double energy_threshold(double eta)