CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ParametrizedSubtractor.cc
Go to the documentation of this file.
8 
10 #include "string"
11 #include "iostream"
12 using namespace std;
13 
15  for ( std::map<int, double>::iterator iter = esigma_.begin();
16  iter != esigma_.end(); ++iter ){
17  iter->second = s*(iter->second);
18  }
19 }
20 
21 
23  PileUpSubtractor(iConfig),
24  dropZeroTowers_(iConfig.getUntrackedParameter<bool>("dropZeroTowers",true)),
25  cbins_(0)
26 {
27  centTag_ = iConfig.getUntrackedParameter<edm::InputTag>("centTag",edm::InputTag("hiCentrality","","RECO"));
28 
29  interpolate_ = iConfig.getParameter<bool>("interpolate");
30  sumRecHits_ = iConfig.getParameter<bool>("sumRecHits");
31 
32  std::string ifname = "RecoHI/HiJetAlgos/data/PU_DATA.root";
33  TFile* inf = new TFile(edm::FileInPath(ifname).fullPath().data());
34  fPU = (TF1*)inf->Get("fPU");
35  fMean = (TF1*)inf->Get("fMean");
36  fRMS = (TF1*)inf->Get("fRMS");
37  hC = (TH1D*)inf->Get("hC");
38 
39  for(int i = 0; i < 40; ++i){
40  hEta.push_back((TH1D*)inf->Get(Form("hEta_%d",i)));
41  hEtaMean.push_back((TH1D*)inf->Get(Form("hEtaMean_%d",i)));
42  hEtaRMS.push_back((TH1D*)inf->Get(Form("hEtaRMS_%d",i)));
43  }
44 
45 }
46 
47 
49  LogDebug("PileUpSubtractor")<<"The subtractor setting up geometry...\n";
50 
51  // if(!cbins_) getCentralityBinsFromDB(iSetup);
52 
54  iEvent.getByLabel(centTag_,cent);
55 
56  centrality_ = cent->EtHFhitSum();
57  bin_ = 40-hC->FindBin(centrality_);
58  if(bin_ > 39) bin_ = 39;
59  if(bin_ < 0) bin_ = 0;
60 
61  if(!geo_) {
63  iSetup.get<CaloGeometryRecord>().get(pG);
64  geo_ = pG.product();
65  std::vector<DetId> alldid = geo_->getValidDetIds();
66 
67  int ietaold = -10000;
68  ietamax_ = -10000;
69  ietamin_ = 10000;
70  for(std::vector<DetId>::const_iterator did=alldid.begin(); did != alldid.end(); did++){
71  if( (*did).det() == DetId::Hcal ){
72  HcalDetId hid = HcalDetId(*did);
73  if( (hid).depth() == 1 ) {
74  allgeomid_.push_back(*did);
75 
76  if((hid).ieta() != ietaold){
77  ietaold = (hid).ieta();
78  geomtowers_[(hid).ieta()] = 1;
79  if((hid).ieta() > ietamax_) ietamax_ = (hid).ieta();
80  if((hid).ieta() < ietamin_) ietamin_ = (hid).ieta();
81  }
82  else{
83  geomtowers_[(hid).ieta()]++;
84  }
85  }
86  }
87  }
88  }
89 
90  for (int i = ietamin_; i<ietamax_+1; i++) {
91  emean_[i] = 0.;
92  esigma_[i] = 0.;
93  ntowersWithJets_[i] = 0;
94  }
95 }
96 
97 
98 void ParametrizedSubtractor::calculatePedestal( vector<fastjet::PseudoJet> const & coll ){
99  return;
100 }
101 
102 void ParametrizedSubtractor::subtractPedestal(vector<fastjet::PseudoJet> & coll)
103 {
104  if(0){
105  return;
106  }else{
107  LogDebug("PileUpSubtractor")<<"The subtractor subtracting pedestals...\n";
108 
109  int it = -100;
110  vector<fastjet::PseudoJet> newcoll;
111 
112  for (vector<fastjet::PseudoJet>::iterator input_object = coll.begin (),
113  fjInputsEnd = coll.end();
114  input_object != fjInputsEnd; ++input_object) {
115 
116  reco::CandidatePtr const & itow = (*inputs_)[ input_object->user_index() ];
117 
118  it = ieta( itow );
119  iphi( itow );
120 
121  double Original_Et = itow->et();
122  if(sumRecHits_){
123  Original_Et = getEt(itow);
124  }
125 
126  double etnew = Original_Et - getPU(it,1,1);
127  float mScale = etnew/input_object->Et();
128  if(etnew < 0.) mScale = 0.;
129 
130  math::XYZTLorentzVectorD towP4(input_object->px()*mScale, input_object->py()*mScale,
131  input_object->pz()*mScale, input_object->e()*mScale);
132 
133  int index = input_object->user_index();
134  input_object->reset ( towP4.px(),
135  towP4.py(),
136  towP4.pz(),
137  towP4.energy() );
138  input_object->set_user_index(index);
139  if(etnew > 0. && dropZeroTowers_) newcoll.push_back(*input_object);
140  }
141  if(dropZeroTowers_) coll = newcoll;
142 
143  }
144 }
145 
146 
147 void ParametrizedSubtractor::calculateOrphanInput(vector<fastjet::PseudoJet> & orphanInput)
148 {
149  orphanInput = *fjInputs_;
150 }
151 
152 
153 
155 {
156 
157  LogDebug("PileUpSubtractor")<<"The subtractor correcting jets...\n";
158  jetOffset_.clear();
159 
160  using namespace reco;
161 
162  (*fjInputs_) = fjOriginalInputs_;
165 
166  if(0){
167  const fastjet::JetDefinition& def = fjClusterSeq_->jet_def();
168  if ( !doAreaFastjet_ && !doRhoFastjet_) {
169  fastjet::ClusterSequence newseq( *fjInputs_, def );
170  (*fjClusterSeq_) = newseq;
171  } else {
172  fastjet::ClusterSequenceArea newseq( *fjInputs_, def , *fjActiveArea_ );
173  (*fjClusterSeq_) = newseq;
174  }
175 
176  (*fjJets_) = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
177  }
178 
179  jetOffset_.reserve(fjJets_->size());
180 
181  vector<fastjet::PseudoJet>::iterator pseudojetTMP = fjJets_->begin (),
182  jetsEnd = fjJets_->end();
183  for (; pseudojetTMP != jetsEnd; ++pseudojetTMP) {
184 
185  int ijet = pseudojetTMP - fjJets_->begin();
186  jetOffset_[ijet] = 0;
187 
188  std::vector<fastjet::PseudoJet> towers =
189  sorted_by_pt(fjClusterSeq_->constituents(*pseudojetTMP));
190 
191  double newjetet = 0.;
192  for(vector<fastjet::PseudoJet>::const_iterator ito = towers.begin(),
193  towEnd = towers.end();
194  ito != towEnd;
195  ++ito)
196  {
197  const reco::CandidatePtr& originalTower = (*inputs_)[ito->user_index()];
198  int it = ieta( originalTower );
199  double Original_Et = originalTower->et();
200 
201  if(sumRecHits_){
202  Original_Et = getEt(originalTower);
203  }
204 
205  double etnew = Original_Et - getPU(it,1,1);
206  if(etnew < 0.) etnew = 0;
207  newjetet = newjetet + etnew;
208  jetOffset_[ijet] += Original_Et - etnew;
209  }
210 
211  if(sumRecHits_){
212  double mScale = newjetet/pseudojetTMP->Et();
213  int cshist = pseudojetTMP->cluster_hist_index();
214  pseudojetTMP->reset(pseudojetTMP->px()*mScale, pseudojetTMP->py()*mScale,
215  pseudojetTMP->pz()*mScale, pseudojetTMP->e()*mScale);
216  pseudojetTMP->set_cluster_hist_index(cshist);
217  }
218 
219  }
220 }
221 
222 
224  const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
225  const GlobalPoint& pos=geo_->getPosition(ctc->id());
226  double energy = ctc->emEnergy() + ctc->hadEnergy();
227 
228  if(0){
229  energy = 0;
230  const std::vector<DetId>& hitids = ctc->constituents();
231  for(unsigned int i = 0; i< hitids.size(); ++i){
232 
233  }
234  }
235 
236  double et = energy*sin(pos.theta());
237  return et;
238 }
239 
241  const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
242  const GlobalPoint& pos=geo_->getPosition(ctc->id());
243  double eta = pos.eta();
244  return eta;
245 }
246 
248  int it = ieta(in);
249  return getPU(it,1,0);
250 }
251 
253  int it = ieta(in);
254  return getPU(it,0,1);
255 }
256 
258  int it = ieta(in);
259  return getPU(it,1,1);
260 }
261 
262 double ParametrizedSubtractor::getPU(int ieta,bool addMean, bool addSigma) const {
263 
264  //double e = hEta[bin_]->GetBinContent(hEta[bin_]->FindBin(ieta));
265  //double c = fPU->Eval(centrality_);
266 
267  double em = hEtaMean[bin_]->GetBinContent(hEtaMean[bin_]->FindBin(ieta));
268  double cm = fMean->Eval(centrality_);
269 
270  double er = hEtaRMS[bin_]->GetBinContent(hEtaRMS[bin_]->FindBin(ieta));
271  double cr = fRMS->Eval(centrality_);
272 
273  if(interpolate_){
274  double n = 0;
275  int hbin = 40-bin_;
276  double centerweight = (centrality_ - hC->GetBinCenter(hbin));
277  double lowerweight = (centrality_ - hC->GetBinLowEdge(hbin));
278  double upperweight = (centrality_ - hC->GetBinLowEdge(hbin+1));
279 
280  em *= lowerweight*upperweight;
281  er *= lowerweight*upperweight;
282  n += lowerweight*upperweight;
283 
284  if(bin_ > 0){
285  em += upperweight*centerweight*hEtaMean[bin_]->GetBinContent(hEtaMean[bin_-1]->FindBin(ieta));
286  er += upperweight*centerweight*hEtaRMS[bin_]->GetBinContent(hEtaRMS[bin_-1]->FindBin(ieta));
287  n += upperweight*centerweight;
288  }
289 
290  if(bin_ < 39){
291  em += lowerweight*centerweight*hEtaMean[bin_]->GetBinContent(hEtaMean[bin_+1]->FindBin(ieta));
292  er += lowerweight*centerweight*hEtaRMS[bin_]->GetBinContent(hEtaRMS[bin_+1]->FindBin(ieta));
293  n += lowerweight*centerweight;
294  }
295  em /= n;
296  er /= n;
297  }
298 
299  // return e*c;
300  return addMean*em*cm + addSigma*nSigmaPU_*er*cr;
301 }
302 
303 
#define LogDebug(id)
double getEt(const reco::CandidatePtr &in) const
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
std::vector< double > jetOffset_
virtual void calculatePedestal(std::vector< fastjet::PseudoJet > const &coll)
double getEta(const reco::CandidatePtr &in) const
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:15
std::vector< fastjet::PseudoJet > * fjJets_
std::map< int, double > esigma_
std::vector< TH1D * > hEtaRMS
std::map< int, int > geomtowers_
std::vector< fastjet::PseudoJet > fjOriginalInputs_
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
virtual double getSigmaAtTower(const reco::CandidatePtr &in) const
virtual double getPileUpAtTower(const reco::CandidatePtr &in) const
std::vector< TH1D * > hEtaMean
T eta() const
ParametrizedSubtractor(const edm::ParameterSet &iConfig)
int ieta(const reco::CandidatePtr &in) const
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
std::vector< fastjet::PseudoJet > * fjInputs_
std::map< int, int > ntowersWithJets_
int iEvent
Definition: GenABIO.cc:243
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:143
double emEnergy() const
Definition: CaloTower.h:79
std::vector< TH1D * > hEta
string inf
Definition: EcalCondDB.py:94
const std::vector< DetId > & constituents() const
Definition: CaloTower.h:73
const GlobalPoint & getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:68
virtual void subtractPedestal(std::vector< fastjet::PseudoJet > &coll)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
virtual double getMeanAtTower(const reco::CandidatePtr &in) const
double getPU(int ieta, bool addMean, bool addSigma) const
double hadEnergy() const
Definition: CaloTower.h:80
int iphi(const reco::CandidatePtr &in) const
CaloTowerDetId id() const
Definition: CaloTower.h:72
JetCorrectorParametersCollection coll
Definition: classes.h:16
ActiveAreaSpecPtr fjActiveArea_
ClusterSequencePtr fjClusterSeq_
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
std::vector< DetId > getValidDetIds() const
Get the list of all valid detector ids.
Definition: CaloGeometry.cc:90
T eta() const
Definition: PV3DBase.h:76
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
CaloGeometry const * geo_
std::map< int, double > emean_
virtual void calculateOrphanInput(std::vector< fastjet::PseudoJet > &orphanInput)
JetCorrectorParameters::Definitions def
Definition: classes.h:12
std::vector< HcalDetId > allgeomid_
virtual void setupGeometryMap(edm::Event &iEvent, const edm::EventSetup &iSetup)