test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FFTGenericScaleCalculator.cc
Go to the documentation of this file.
1 #include <cassert>
2 #include <cfloat>
3 
7 
8 #define int_param(varname) m_ ## varname (ps.getParameter< int >( #varname ))
9 
10 #define check_param(varname) if ((m_ ## varname) >= 0) {\
11  if ((m_ ## varname) >= nFactors)\
12  throw cms::Exception("FFTJetBadConfig")\
13  << "In FFTGenericScaleCalculator constructor: "\
14  << "out of range mapping for variable \""\
15  << #varname << "\"" << std::endl;\
16  mask[(m_ ## varname)] = 1;\
17  ++dim;\
18  }
19 
20 
21 static inline double delPhi(const double phi1, const double phi2)
22 {
23  double dphi = phi1 - phi2;
24  if (dphi > M_PI)
25  dphi -= 2.0*M_PI;
26  else if (dphi < -M_PI)
27  dphi += 2.0*M_PI;
28  return dphi;
29 }
30 
31 
33  : m_factors(ps.getParameter<std::vector<double> >("factors")),
34  m_minLog(ps.getUntrackedParameter<double>("minLog", -800.0)),
35  int_param(eta),
36  int_param(phi),
37  int_param(pt),
38  int_param(logPt),
39  int_param(mass),
40  int_param(logMass),
42  int_param(logEnergy),
43  int_param(gamma),
44  int_param(logGamma),
46  int_param(ncells),
47  int_param(etSum),
48  int_param(etaWidth),
49  int_param(phiWidth),
50  int_param(averageWidth),
51  int_param(widthRatio),
52  int_param(etaPhiCorr),
53  int_param(fuzziness),
54  int_param(convergenceDistance),
55  int_param(recoScale),
56  int_param(recoScaleRatio),
57  int_param(membershipFactor),
58  int_param(magnitude),
59  int_param(logMagnitude),
60  int_param(magS1),
61  int_param(LogMagS1),
62  int_param(magS2),
63  int_param(LogMagS2),
64  int_param(driftSpeed),
65  int_param(magSpeed),
66  int_param(lifetime),
67  int_param(splitTime),
68  int_param(mergeTime),
70  int_param(logScale),
71  int_param(nearestNeighborDistance),
72  int_param(clusterRadius),
73  int_param(clusterSeparation),
74  int_param(dRFromJet),
75  int_param(LaplacianS1),
76  int_param(LaplacianS2),
77  int_param(LaplacianS3),
78  int_param(HessianS2),
79  int_param(HessianS4),
80  int_param(HessianS6),
81  int_param(nConstituents),
82  int_param(aveConstituentPt),
83  int_param(logAveConstituentPt),
84  int_param(constituentPtDistribution),
85  int_param(constituentEtaPhiSpread),
91  int_param(HFHadronEnergyFraction),
92  int_param(HFEMEnergyFraction),
98  int_param(HFHadronMultiplicity),
99  int_param(HFEMMultiplicity),
100  int_param(chargedEmEnergyFraction),
101  int_param(chargedMuEnergyFraction),
102  int_param(neutralEmEnergyFraction),
103  int_param(EmEnergyFraction),
104  int_param(chargedMultiplicity),
105  int_param(neutralMultiplicity)
106 {
107  const int nFactors = m_factors.size();
108  std::vector<int> mask(nFactors, 0);
109  int dim = 0;
110 
111  check_param(eta);
112  check_param(phi);
113  check_param(pt);
114  check_param(logPt);
115  check_param(mass);
116  check_param(logMass);
118  check_param(logEnergy);
119  check_param(gamma);
120  check_param(logGamma);
122  check_param(ncells);
123  check_param(etSum);
124  check_param(etaWidth);
125  check_param(phiWidth);
126  check_param(averageWidth);
127  check_param(widthRatio);
128  check_param(etaPhiCorr);
129  check_param(fuzziness);
130  check_param(convergenceDistance);
131  check_param(recoScale);
132  check_param(recoScaleRatio);
133  check_param(membershipFactor);
134  check_param(magnitude);
135  check_param(logMagnitude);
136  check_param(magS1);
137  check_param(LogMagS1);
138  check_param(magS2);
139  check_param(LogMagS2);
140  check_param(driftSpeed);
141  check_param(magSpeed);
142  check_param(lifetime);
143  check_param(splitTime);
144  check_param(mergeTime);
146  check_param(logScale);
147  check_param(nearestNeighborDistance);
148  check_param(clusterRadius);
149  check_param(clusterSeparation);
150  check_param(dRFromJet);
151  check_param(LaplacianS1);
152  check_param(LaplacianS2);
153  check_param(LaplacianS3);
154  check_param(HessianS2);
155  check_param(HessianS4);
156  check_param(HessianS6);
157  check_param(nConstituents);
158  check_param(aveConstituentPt);
159  check_param(logAveConstituentPt);
160  check_param(constituentPtDistribution);
161  check_param(constituentEtaPhiSpread);
167  check_param(HFHadronEnergyFraction);
168  check_param(HFEMEnergyFraction);
174  check_param(HFHadronMultiplicity);
175  check_param(HFEMMultiplicity);
176  check_param(chargedEmEnergyFraction);
177  check_param(chargedMuEnergyFraction);
178  check_param(neutralEmEnergyFraction);
179  check_param(EmEnergyFraction);
180  check_param(chargedMultiplicity);
181  check_param(neutralMultiplicity);
182 
183  if (dim != nFactors)
184  throw cms::Exception("FFTJetBadConfig")
185  << "In FFTGenericScaleCalculator constructor: "
186  << "incompatible number of scaling factors: expected "
187  << dim << ", got " << nFactors << std::endl;
188  for (int i=0; i<nFactors; ++i)
189  if (mask[i] == 0)
190  throw cms::Exception("FFTJetBadConfig")
191  << "In FFTGenericScaleCalculator constructor: "
192  << "variable number " << i << " is not mapped" << std::endl;
193 }
194 
196  const reco::Jet& jet, const reco::FFTJet<float>& fftJet,
197  const math::XYZTLorentzVector& current,
198  double* buf, const unsigned dim) const
199 {
200  // Verify that the input is reasonable
201  if (dim != m_factors.size())
202  throw cms::Exception("FFTJetBadConfig")
203  << "In FFTGenericScaleCalculator::mapFFTJet: "
204  << "incompatible table dimensionality: expected "
205  << m_factors.size() << ", got " << dim << std::endl;
206  if (dim)
207  assert(buf);
208  else
209  return;
210 
211  // Go over all variables and map them as configured.
212  // Variables from the "current" Lorentz vector.
213  if (m_eta >= 0)
214  buf[m_eta] = current.eta();
215 
216  if (m_phi >= 0)
217  buf[m_phi] = current.phi();
218 
219  if (m_pt >= 0)
220  buf[m_pt] = current.pt();
221 
222  if (m_logPt >= 0)
223  buf[m_logPt] = f_safeLog(current.pt());
224 
225  if (m_mass >= 0)
226  buf[m_mass] = current.M();
227 
228  if (m_logMass >= 0)
229  buf[m_mass] = f_safeLog(current.M());
230 
231  if (m_energy >= 0)
232  buf[m_energy] = current.e();
233 
234  if (m_logEnergy >= 0)
235  buf[m_energy] = f_safeLog(current.e());
236 
237  if (m_gamma >= 0)
238  {
239  const double m = current.M();
240  if (m > 0.0)
241  buf[m_gamma] = current.e()/m;
242  else
243  buf[m_gamma] = DBL_MAX;
244  }
245 
246  if (m_logGamma >= 0)
247  {
248  const double m = current.M();
249  if (m > 0.0)
250  buf[m_gamma] = current.e()/m;
251  else
252  buf[m_gamma] = DBL_MAX;
253  buf[m_gamma] = log(buf[m_gamma]);
254  }
255 
256  // Variables from fftJet
257  if (m_pileup >= 0)
258  buf[m_pileup] = fftJet.f_pileup().pt();
259 
260  if (m_ncells >= 0)
261  buf[m_ncells] = fftJet.f_ncells();
262 
263  if (m_etSum >= 0)
264  buf[m_etSum] = fftJet.f_etSum();
265 
266  if (m_etaWidth >= 0)
267  buf[m_etaWidth] = fftJet.f_etaWidth();
268 
269  if (m_phiWidth >= 0)
270  buf[m_phiWidth] = fftJet.f_phiWidth();
271 
272  if (m_averageWidth >= 0)
273  {
274  const double etaw = fftJet.f_etaWidth();
275  const double phiw = fftJet.f_phiWidth();
276  buf[m_averageWidth] = sqrt(etaw*etaw + phiw*phiw);
277  }
278 
279  if (m_widthRatio >= 0)
280  {
281  const double etaw = fftJet.f_etaWidth();
282  const double phiw = fftJet.f_phiWidth();
283  if (phiw > 0.0)
284  buf[m_widthRatio] = etaw/phiw;
285  else
286  buf[m_widthRatio] = DBL_MAX;
287  }
288 
289  if (m_etaPhiCorr >= 0)
290  buf[m_etaPhiCorr] = fftJet.f_etaPhiCorr();
291 
292  if (m_fuzziness >= 0)
293  buf[m_fuzziness] = fftJet.f_fuzziness();
294 
295  if (m_convergenceDistance >= 0)
297 
298  if (m_recoScale >= 0)
299  buf[m_recoScale] = fftJet.f_recoScale();
300 
301  if (m_recoScaleRatio >= 0)
302  buf[m_recoScaleRatio] = fftJet.f_recoScaleRatio();
303 
304  if (m_membershipFactor >= 0)
305  buf[m_membershipFactor] = fftJet.f_membershipFactor();
306 
307  // Get most often used precluster quantities
308  const reco::PattRecoPeak<float>& preclus = fftJet.f_precluster();
309  const double scale = preclus.scale();
310 
311  if (m_magnitude >= 0)
312  buf[m_magnitude] = preclus.magnitude();
313 
314  if (m_logMagnitude >= 0)
315  buf[m_logMagnitude] = f_safeLog(preclus.magnitude());
316 
317  if (m_magS1 >= 0)
318  buf[m_magS1] = preclus.magnitude()*scale;
319 
320  if (m_LogMagS1 >= 0)
321  buf[m_LogMagS1] = f_safeLog(preclus.magnitude()*scale);
322 
323  if (m_magS2 >= 0)
324  buf[m_magS2] = preclus.magnitude()*scale*scale;
325 
326  if (m_LogMagS2 >= 0)
327  buf[m_LogMagS2] = f_safeLog(preclus.magnitude()*scale*scale);
328 
329  if (m_driftSpeed >= 0)
330  buf[m_driftSpeed] = preclus.driftSpeed();
331 
332  if (m_magSpeed >= 0)
333  buf[m_magSpeed] = preclus.magSpeed();
334 
335  if (m_lifetime >= 0)
336  buf[m_lifetime] = preclus.lifetime();
337 
338  if (m_splitTime >= 0)
339  buf[m_splitTime] = preclus.splitTime();
340 
341  if (m_mergeTime >= 0)
342  buf[m_mergeTime] = preclus.mergeTime();
343 
344  if (m_scale >= 0)
345  buf[m_scale] = scale;
346 
347  if (m_logScale >= 0)
348  buf[m_logScale] = f_safeLog(scale);
349 
350  if (m_nearestNeighborDistance >= 0)
352 
353  if (m_clusterRadius >= 0)
354  buf[m_clusterRadius] = preclus.clusterRadius();
355 
356  if (m_clusterSeparation >= 0)
357  buf[m_clusterSeparation] = preclus.clusterSeparation();
358 
359  if (m_dRFromJet >= 0)
360  {
361  const double deta = preclus.eta() - current.eta();
362  const double dphi = delPhi(preclus.phi(), current.phi());
363  buf[m_dRFromJet] = sqrt(deta*deta + dphi*dphi);
364  }
365 
366  if (m_LaplacianS1 >= 0)
367  {
368  double h[3];
369  preclus.hessian(h);
370  buf[m_LaplacianS1] = fabs(h[0] + h[2])*scale;
371  }
372 
373  if (m_LaplacianS2 >= 0)
374  {
375  double h[3];
376  preclus.hessian(h);
377  buf[m_LaplacianS2] = fabs(h[0] + h[2])*scale*scale;
378  }
379 
380  if (m_LaplacianS3 >= 0)
381  {
382  double h[3];
383  preclus.hessian(h);
384  buf[m_LaplacianS3] = fabs(h[0] + h[2])*scale*scale*scale;
385  }
386 
387  if (m_HessianS2 >= 0)
388  {
389  double h[3];
390  preclus.hessian(h);
391  buf[m_HessianS2] = fabs(h[0]*h[2] - h[1]*h[1])*scale*scale;
392  }
393 
394  if (m_HessianS4 >= 0)
395  {
396  double h[3];
397  preclus.hessian(h);
398  buf[m_HessianS4] = fabs(h[0]*h[2] - h[1]*h[1])*pow(scale, 4);
399  }
400 
401  if (m_HessianS6 >= 0)
402  {
403  double h[3];
404  preclus.hessian(h);
405  buf[m_HessianS6] = fabs(h[0]*h[2] - h[1]*h[1])*pow(scale, 6);
406  }
407 
408  // Variables from reco::Jet
409  if (m_nConstituents >= 0)
410  buf[m_nConstituents] = jet.nConstituents();
411 
412  if (m_aveConstituentPt >= 0)
413  buf[m_aveConstituentPt] = current.pt()/jet.nConstituents();
414 
415  if (m_logAveConstituentPt >= 0)
416  buf[m_logAveConstituentPt] = f_safeLog(current.pt()/jet.nConstituents());
417 
420 
421  if (m_constituentEtaPhiSpread >= 0)
423 
424  // Variables from reco::PFJet
425  const reco::PFJet* pfjet = dynamic_cast<const reco::PFJet*>(&jet);
426  if (pfjet)
427  {
428  // Particle flow jet
431 
434 
435  if (m_photonEnergyFraction >= 0)
437 
438  if (m_electronEnergyFraction >= 0)
440 
441  if (m_muonEnergyFraction >= 0)
443 
444  if (m_HFHadronEnergyFraction >= 0)
446 
447  if (m_HFEMEnergyFraction >= 0)
449 
452 
455 
456  if (m_photonMultiplicity >= 0)
458 
459  if (m_electronMultiplicity >= 0)
461 
462  if (m_muonMultiplicity >= 0)
463  buf[m_muonMultiplicity] = pfjet->muonMultiplicity();
464 
465  if (m_HFHadronMultiplicity >= 0)
467 
468  if (m_HFEMMultiplicity >= 0)
469  buf[m_HFEMMultiplicity] = pfjet->HFEMMultiplicity();
470 
471  if (m_chargedEmEnergyFraction >= 0)
473 
474  if (m_chargedMuEnergyFraction >= 0)
476 
477  if (m_neutralEmEnergyFraction >= 0)
479 
480  if (m_EmEnergyFraction >= 0)
482  pfjet->chargedEmEnergyFraction();
483 
484  if (m_chargedMultiplicity >= 0)
486 
487  if (m_neutralMultiplicity >= 0)
489  }
490  else
491  {
492  // Not a particle flow jet
495  m_photonEnergyFraction >= 0 ||
497  m_muonEnergyFraction >= 0 ||
499  m_HFEMEnergyFraction >= 0 ||
502  m_photonMultiplicity >= 0 ||
503  m_electronMultiplicity >= 0 ||
504  m_muonMultiplicity >= 0 ||
505  m_HFHadronMultiplicity >= 0 ||
506  m_HFEMMultiplicity >= 0 ||
510  m_EmEnergyFraction >= 0 ||
511  m_chargedMultiplicity >= 0 ||
513  throw cms::Exception("FFTJetBadConfig")
514  << "In FFTGenericScaleCalculator::mapFFTJet: "
515  << "this configuration is valid for particle flow jets only"
516  << std::endl;
517  }
518 
519  // Apply the scaling factors
520  for (unsigned i=0; i<dim; ++i)
521  buf[i] *= m_factors[i];
522 }
Real scale() const
Definition: PattRecoPeak.h:70
int i
Definition: DBlmapReader.cc:9
int photonMultiplicity() const
photonMultiplicity
Definition: PFJet.h:131
float constituentEtaPhiSpread() const
Real eta() const
Definition: PattRecoPeak.h:62
Real phi() const
Definition: PattRecoPeak.h:63
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
Real f_convergenceDistance() const
Definition: FFTJet.h:72
float chargedHadronEnergyFraction() const
chargedHadronEnergyFraction
Definition: PFJet.h:100
Base class for all types of Jets.
Definition: Jet.h:20
assert(m_qm.get())
Real mergeTime() const
Definition: PattRecoPeak.h:69
Real driftSpeed() const
Definition: PattRecoPeak.h:65
float HFEMEnergyFraction() const
HFEMEnergyFraction.
Definition: PFJet.h:124
int chargedMultiplicity() const
chargedMultiplicity
Definition: PFJet.h:155
float photonEnergyFraction() const
photonEnergyFraction
Definition: PFJet.h:108
Jets made from PFObjects.
Definition: PFJet.h:21
Real f_membershipFactor() const
Definition: FFTJet.h:75
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
static double delPhi(const double phi1, const double phi2)
float neutralEmEnergyFraction() const
neutralEmEnergyFraction
Definition: PFJet.h:152
FFTGenericScaleCalculator(const edm::ParameterSet &ps)
Real f_etSum() const
Definition: FFTJet.h:65
Real f_phiWidth() const
Definition: FFTJet.h:69
Real nearestNeighborDistance() const
Definition: PattRecoPeak.h:71
Real f_recoScale() const
Definition: FFTJet.h:73
Real lifetime() const
Definition: PattRecoPeak.h:67
T sqrt(T t)
Definition: SSEVec.h:18
Real splitTime() const
Definition: PattRecoPeak.h:68
Real f_ncells() const
Definition: FFTJet.h:64
int neutralHadronMultiplicity() const
neutralHadronMultiplicity
Definition: PFJet.h:129
Real clusterSeparation() const
Definition: PattRecoPeak.h:73
float neutralHadronEnergyFraction() const
neutralHadronEnergyFraction
Definition: PFJet.h:104
int neutralMultiplicity() const
neutralMultiplicity
Definition: PFJet.h:157
#define M_PI
Real magSpeed() const
Definition: PattRecoPeak.h:66
const PattRecoPeak< Real > & f_precluster() const
Definition: FFTJet.h:61
Real f_etaWidth() const
Definition: FFTJet.h:68
float chargedEmEnergyFraction() const
chargedEmEnergyFraction
Definition: PFJet.h:144
float HFHadronEnergyFraction() const
HFHadronEnergyFraction.
Definition: PFJet.h:120
float electronEnergyFraction() const
electronEnergyFraction
Definition: PFJet.h:112
#define int_param(varname)
Geom::Phi< T > phi() const
int chargedHadronMultiplicity() const
chargedHadronMultiplicity
Definition: PFJet.h:127
#define check_param(varname)
int muonMultiplicity() const
muonMultiplicity
Definition: PFJet.h:135
int HFEMMultiplicity() const
HFEMMultiplicity.
Definition: PFJet.h:139
float constituentPtDistribution() const
virtual int nConstituents() const
of constituents
Definition: Jet.h:65
Real clusterRadius() const
Definition: PattRecoPeak.h:72
Real f_etaPhiCorr() const
Definition: FFTJet.h:70
float chargedMuEnergyFraction() const
chargedMuEnergyFraction
Definition: PFJet.h:148
const math::XYZTLorentzVector & f_pileup() const
Definition: FFTJet.h:63
double f_safeLog(const double x) const
float muonEnergyFraction() const
muonEnergyFraction
Definition: PFJet.h:116
virtual void mapFFTJet(const reco::Jet &jet, const reco::FFTJet< float > &fftJet, const math::XYZTLorentzVector &current, double *buf, unsigned dim) const
Real magnitude() const
Definition: PattRecoPeak.h:64
Real f_recoScaleRatio() const
Definition: FFTJet.h:74
int HFHadronMultiplicity() const
HFHadronMultiplicity.
Definition: PFJet.h:137
void hessian(double hessianArray[3]) const
Definition: PattRecoPeak.h:74
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
int electronMultiplicity() const
electronMultiplicity
Definition: PFJet.h:133
Real f_fuzziness() const
Definition: FFTJet.h:71