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