CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
VoronoiSubtractor.cc
Go to the documentation of this file.
6 using namespace std;
7 
8 bool VoronoiSubtractor::match(const fastjet::PseudoJet cand1, const fastjet::PseudoJet cand2){
9  return (cand1.delta_R(cand2) < rParam_);
10 }
11 
12 
14  PileUpSubtractor(iConfig, std::move(iC)),
15  srcCand_(iC.consumes<reco::CandidateView>(iConfig.getParameter<edm::InputTag>("src"))),
16  srcVor_(iC.consumes<edm::ValueMap<reco::VoronoiBackground> >(iConfig.getParameter<edm::InputTag>("bkg"))),
17  dropZeroTowers_(iConfig.getParameter<bool>("dropZeros")),
18  addNegative_(iConfig.getParameter<bool>("addNegative")),
19  addNegativesFromCone_(iConfig.getParameter<bool>("addNegativesFromCone")),
20  infinitesimalPt_(iConfig.getParameter<double>("infinitesimalPt")),
21  rParam_(iConfig.getParameter<double>("rParam"))
22 {
23 
24 }
25 
26 
27 
29 {
30 
31  LogDebug("VoronoiSubtractor")<<"The subtractor retrieving Voronoi background...\n";
32  geo_ = 0;
33  droppedCandidates_.clear();
34  jetOffset_.clear();
35 
38 
39 }
40 
42 {
43 
44  LogDebug("VoronoiSubtractor")<<"finalizing the output...\n";
45 
46  jetOffset_.reserve(fjJets_->size());
47 
48  for (unsigned int ijet = 0;ijet <fjJets_->size();++ijet) {
49  fastjet::PseudoJet& fjJet = (*fjJets_)[ijet];
50 
51  LogDebug("VoronoiSubtractor")<<"fjJets_ "<<ijet<<" pt : "<<fjJet.pt()
52  <<" --- eta : "<<fjJet.eta()<<" --- phi : "<<fjJet.phi()<<endl;
53 
54  fastjet::PseudoJet subtracted;
55  fastjet::PseudoJet unsubtracted;
56  fastjet::PseudoJet unsubtractedDropped;
57  jetOffset_[ijet] = 0;
58 
59  // Loop over constituents to determine the momentum
60  // before and after subtraction
61 
62  std::vector<fastjet::PseudoJet> fjConstituents = fastjet::sorted_by_pt(fjJet.constituents());
63  for (unsigned int i=0;i<fjConstituents.size();++i) {
64  unsigned int index = fjConstituents[i].user_index();
65 
67  const reco::VoronoiBackground& voronoi = (*backgrounds_)[ref];
68 
69  fastjet::PseudoJet candidate(ref->px(),ref->py(),ref->pz(),ref->energy());
70  double orpt = candidate.perp();
71  unsubtracted += candidate;
72  if(addNegative_ || voronoi.pt() > 0){
73  candidate.reset_PtYPhiM(addNegative_ ? voronoi.pt_subtracted() : voronoi.pt(),ref->rapidity(),ref->phi(),ref->mass());
74  LogDebug("VoronoiSubtractor")<<"candidate "<<index
75  <<" --- original pt : "<<orpt
76  <<" --- voronoi pt : "<<voronoi.pt()
77  <<" --- ref pt : "<<ref->pt()
78  <<" --- constituent "<<index
79  <<" --- pt : "<<fjConstituents[i].perp()
80  <<" --- eta : "<<fjConstituents[i].pseudorapidity()
81  <<" --- phi : "<<fjConstituents[i].phi()
82  <<" --- mass : "<<fjConstituents[i].m()<<endl;
83  subtracted += candidate;
84  }
85  }
86 
87  // Loop over dropped candidates to see whether there is any of them
88  // that would belong to this jet:
90  for(unsigned int i=0; i < droppedCandidates_.size(); ++i){
92  fastjet::PseudoJet dropcand(ref->px(),ref->py(),ref->pz(),ref->energy());
93 
94  if(match(fjJet,dropcand)){
95  unsubtractedDropped += dropcand;
96  unsubtracted += dropcand;
97  }
98  }
99  }
100 
101  fjJet.reset_momentum(subtracted);
102 
103  LogDebug("VoronoiSubtractor")<<"fjJets_ "<<ijet<<" unsubtracted : "<<unsubtracted.pt()<<endl;
104  LogDebug("VoronoiSubtractor")<<"fjJets_ "<<ijet<<" subtracted : "<<subtracted.pt()<<endl;
105  LogDebug("VoronoiSubtractor")<<"fjJets_ "<<ijet<<" dropped : "<<unsubtractedDropped.pt()<<endl;
106  jetOffset_[ijet] = unsubtracted.pt() - fjJet.pt();
107 
108  }
109 
110 
111 
112 
113 
114 
115 
116 
117 
118 
119 }
120 
121 void VoronoiSubtractor::subtractPedestal(vector<fastjet::PseudoJet> & coll)
122 {
123 
124  LogDebug("VoronoiSubtractor")<<"The subtractor subtracting pedestals...\n";
125  vector<fastjet::PseudoJet> newcoll;
126 
127  for (vector<fastjet::PseudoJet>::iterator input_object = coll.begin (),
128  fjInputsEnd = coll.end();
129  input_object != fjInputsEnd; ++input_object) {
130 
131  reco::CandidateViewRef ref(candidates_,input_object->user_index());
132  const reco::VoronoiBackground& voronoi = (*backgrounds_)[ref];
133 
134  double ptold = input_object->pt();
135  double ptnew = voronoi.pt();
136 
137  LogDebug("VoronoiSubtractor")<<"pt old : "<<ptold<<" ; pt new : "<<ptnew
138  <<" E : "<<input_object->e()
139  <<" rap : "<<input_object->rapidity()
140  <<" phi : "<<input_object->phi()
141  <<" MASS : "<<input_object->m()<<endl;
142 
143  // Treatment of candidates with negative pt after subtraction
144  if(ptnew < infinitesimalPt_){
145  if(infinitesimalPt_ > 0){
146  // Low-pt candidate is assigned a very small finite pt
147  // so that the jet clustering includes the candidate
148  // and can associate it to the jet.
149  // The original candidate pt is restored
150  // in the offsetCorrectJets() function.
151  ptnew = infinitesimalPt_;
152  }else{
153  // Low-pt candidate is removed from the input collection,
154  // so that the jet clustering algorithm can function properly.
155  // However, we need to keep track of these candidates
156  // in order to determine how much energy has been subtracted in total.
157  droppedCandidates_.push_back(input_object->user_index());
158  continue;
159  }
160  }
161 
162  int index = input_object->user_index();
163 
164  fastjet::PseudoJet ps(input_object->four_mom());
165  ps.reset_PtYPhiM(ptnew,input_object->rapidity(),input_object->phi(),input_object->m());
166  ps.set_user_index(index);
167 
168  LogDebug("VoronoiSubtractor")<<"New momentum : "<<ps.pt()
169  <<" rap : "<<ps.rap()
170  <<" phi : "<<ps.phi()
171  <<" MASS : "<<ps.m()<<endl;
172 
173  newcoll.push_back(ps);
174  }
175  coll = newcoll;
176 }
177 
178 void VoronoiSubtractor::calculatePedestal( vector<fastjet::PseudoJet> const & coll )
179 {
180  LogDebug("VoronoiSubtractor")<<"do nothing...\n";
181 }
182 
183 
184 void VoronoiSubtractor::calculateOrphanInput(vector<fastjet::PseudoJet> & orphanInput)
185 {
186  LogDebug("VoronoiSubtractor")<<"do nothing...\n";
187 }
188 
189 
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
std::vector< double > jetOffset_
virtual void subtractPedestal(std::vector< fastjet::PseudoJet > &coll)
virtual void offsetCorrectJets()
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
std::vector< fastjet::PseudoJet > * fjJets_
edm::EDGetTokenT< edm::ValueMap< reco::VoronoiBackground > > srcVor_
edm::Handle< reco::VoronoiMap > backgrounds_
virtual void calculateOrphanInput(std::vector< fastjet::PseudoJet > &orphanInput)
std::vector< int > droppedCandidates_
virtual void calculatePedestal(std::vector< fastjet::PseudoJet > const &coll)
bool match(fastjet::PseudoJet, fastjet::PseudoJet)
edm::EDGetTokenT< reco::CandidateView > srcCand_
int iEvent
Definition: GenABIO.cc:230
virtual void setupGeometryMap(edm::Event &iEvent, const edm::EventSetup &iSetup)
JetCorrectorParametersCollection coll
Definition: classes.h:10
double pt_subtracted() const
CaloGeometry const * geo_
edm::Handle< reco::CandidateView > candidates_
VoronoiSubtractor(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)