CMS 3D CMS Logo

PuppiContainer.cc
Go to the documentation of this file.
3 #include "Math/ProbFunc.h"
4 #include "TMath.h"
5 #include <iostream>
6 #include <cmath>
9 
10 using namespace std;
11 
13  fPuppiDiagnostics = iConfig.getParameter<bool>("puppiDiagnostics");
14  fApplyCHS = iConfig.getParameter<bool>("applyCHS");
15  fInvert = iConfig.getParameter<bool>("invertPuppi");
16  fUseExp = iConfig.getParameter<bool>("useExp");
17  fPuppiWeightCut = iConfig.getParameter<double>("MinPuppiWeight");
18  fPtMaxPhotons = iConfig.getParameter<double>("PtMaxPhotons");
19  fEtaMaxPhotons = iConfig.getParameter<double>("EtaMaxPhotons");
20  fPtMaxNeutrals = iConfig.getParameter<double>("PtMaxNeutrals");
21  fPtMaxNeutralsStartSlope = iConfig.getParameter<double>("PtMaxNeutralsStartSlope");
22  std::vector<edm::ParameterSet> lAlgos = iConfig.getParameter<std::vector<edm::ParameterSet> >("algos");
23  fNAlgos = lAlgos.size();
24  for (unsigned int i0 = 0; i0 < lAlgos.size(); i0++) {
25  PuppiAlgo pPuppiConfig(lAlgos[i0]);
26  fPuppiAlgo.push_back(pPuppiConfig);
27  }
28 }
29 
30 void PuppiContainer::initialize(const std::vector<RecoObj> &iRecoObjects) {
31  //Clear everything
32  fPFParticles.resize(0);
33  fPFParticlesForVar.resize(0);
34  fPFParticlesForVarChargedPV.resize(0);
35  fWeights.resize(0);
36  fVals.resize(0);
37  fRawAlphas.resize(0);
38  fAlphaMed.resize(0);
39  fAlphaRMS.resize(0);
40  fPUProxy = 1.;
41  //Link to the RecoObjects
42  fRecoParticles = &iRecoObjects;
43  fPFParticles.reserve(iRecoObjects.size());
44  fPFParticlesForVar.reserve(iRecoObjects.size());
45  fPFParticlesForVarChargedPV.reserve(iRecoObjects.size());
46  for (auto const &rParticle : *fRecoParticles) {
47  PuppiCandidate pCand;
48  pCand.id = rParticle.id;
49  if (edm::isFinite(rParticle.rapidity)) {
50  pCand.pt = rParticle.pt;
51  pCand.eta = rParticle.eta;
52  pCand.rapidity = rParticle.rapidity;
53  pCand.phi = rParticle.phi;
54  pCand.m = rParticle.m;
55  } else {
56  pCand.pt = 0.;
57  pCand.eta = 99.;
58  pCand.rapidity = 99.;
59  pCand.phi = 0.;
60  pCand.m = 0.;
61  }
62 
63  fPFParticles.push_back(pCand);
64 
65  // skip candidates to be ignored in the computation
66  // of PUPPI's alphas (e.g. electrons and muons if puppiNoLep=True)
67  if (std::abs(rParticle.id) == 3)
68  continue;
69 
70  fPFParticlesForVar.push_back(pCand);
71  // charged candidates assigned to LV
72  if (std::abs(rParticle.id) == 1)
73  fPFParticlesForVarChargedPV.push_back(pCand);
74  }
75 }
76 
78 
80  std::vector<PuppiCandidate> const &iParts,
81  int iOpt,
82  const double iRCone) {
83  return var_within_R(iOpt, iParts, iPart, iRCone);
84 }
85 
87  const vector<PuppiCandidate> &particles,
88  const PuppiCandidate &centre,
89  const double R) {
90  if (iId == -1)
91  return 1.;
92 
93  double const r2 = R * R;
94  double var = 0.;
95 
96  for (auto const &cand : particles) {
97  if (std::abs(cand.rapidity - centre.rapidity) < R) {
98  auto const dr2y = reco::deltaR2(cand.rapidity, cand.phi, centre.rapidity, centre.phi);
99  if (dr2y < r2) {
100  auto const dr2 = reco::deltaR2(cand.eta, cand.phi, centre.eta, centre.phi);
101  if (dr2 < 0.0001)
102  continue;
103  auto const pt = cand.pt;
104  switch (iId) {
105  case 5:
106  var += (pt * pt / dr2);
107  break;
108  case 4:
109  var += pt;
110  break;
111  case 3:
112  var += (1. / dr2);
113  break;
114  case 2:
115  var += (1. / dr2);
116  break;
117  case 1:
118  var += pt;
119  break;
120  case 0:
121  var += (pt / dr2);
122  break;
123  }
124  }
125  }
126  }
127 
128  if ((var != 0.) and ((iId == 0) or (iId == 3) or (iId == 5)))
129  var = log(var);
130  else if (iId == 1)
131  var += centre.pt;
132 
133  return var;
134 }
135 
136 //In fact takes the median not the average
138  std::vector<PuppiCandidate> const &iConstits,
139  std::vector<PuppiCandidate> const &iParticles,
140  std::vector<PuppiCandidate> const &iChargedParticles) {
141  for (unsigned int i0 = 0; i0 < iConstits.size(); i0++) {
142  //Calculate the Puppi Algo to use
143  int pPupId = getPuppiId(iConstits[i0].pt, iConstits[i0].eta);
144  if (pPupId == -1 || fPuppiAlgo[pPupId].numAlgos() <= iOpt) {
145  fVals.push_back(-1);
146  continue;
147  }
148  //Get the Puppi Sub Algo (given iteration)
149  int pAlgo = fPuppiAlgo[pPupId].algoId(iOpt);
150  bool pCharged = fPuppiAlgo[pPupId].isCharged(iOpt);
151  double pCone = fPuppiAlgo[pPupId].coneSize(iOpt);
152  // compute the Puppi metric:
153  // - calculate goodVar only for candidates that (1) will not be assigned a predefined weight (e.g 0, 1),
154  // or (2) are required for computations inside puppi-algos (see call to PuppiAlgo::add below)
155  double pVal = -1;
156  bool const getsDefaultWgtIfApplyCHS = iConstits[i0].id == 1 or iConstits[i0].id == 2;
157  if (not((fApplyCHS and getsDefaultWgtIfApplyCHS) or iConstits[i0].id == 3) or
158  (std::abs(iConstits[i0].eta) < fPuppiAlgo[pPupId].etaMaxExtrap() and getsDefaultWgtIfApplyCHS)) {
159  pVal = goodVar(iConstits[i0], pCharged ? iChargedParticles : iParticles, pAlgo, pCone);
160  }
161  fVals.push_back(pVal);
162 
163  if (!edm::isFinite(pVal)) {
164  LogDebug("NotFound") << "====> Value is Nan " << pVal << " == " << iConstits[i0].pt << " -- " << iConstits[i0].eta
165  << endl;
166  continue;
167  }
168 
169  // code added by Nhan: now instead for every algorithm give it all the particles
170  for (int i1 = 0; i1 < fNAlgos; i1++) {
171  // skip cands outside of algo's etaMaxExtrap, as they would anyway be ignored inside PuppiAlgo::add (see end of the block)
172  if (not(std::abs(iConstits[i0].eta) < fPuppiAlgo[i1].etaMaxExtrap() and getsDefaultWgtIfApplyCHS))
173  continue;
174 
175  auto curVal = pVal;
176  // recompute goodVar if algo has changed
177  if (i1 != pPupId) {
178  pAlgo = fPuppiAlgo[i1].algoId(iOpt);
179  pCharged = fPuppiAlgo[i1].isCharged(iOpt);
180  pCone = fPuppiAlgo[i1].coneSize(iOpt);
181  curVal = goodVar(iConstits[i0], pCharged ? iChargedParticles : iParticles, pAlgo, pCone);
182  }
183 
184  fPuppiAlgo[i1].add(iConstits[i0], curVal, iOpt);
185  }
186  }
187 
188  for (int i0 = 0; i0 < fNAlgos; i0++)
189  fPuppiAlgo[i0].computeMedRMS(iOpt);
190 }
191 
192 //In fact takes the median not the average
194  std::vector<PuppiCandidate> const &iConstits,
195  std::vector<PuppiCandidate> const &iParticles,
196  std::vector<PuppiCandidate> const &iChargedParticles) {
197  for (int j0 = 0; j0 < fNAlgos; j0++) {
198  for (unsigned int i0 = 0; i0 < iConstits.size(); i0++) {
199  //Get the Puppi Sub Algo (given iteration)
200  int pAlgo = fPuppiAlgo[j0].algoId(iOpt);
201  bool pCharged = fPuppiAlgo[j0].isCharged(iOpt);
202  double pCone = fPuppiAlgo[j0].coneSize(iOpt);
203  //Compute the Puppi Metric
204  double const pVal = goodVar(iConstits[i0], pCharged ? iChargedParticles : iParticles, pAlgo, pCone);
205  fRawAlphas.push_back(pVal);
206  if (!edm::isFinite(pVal)) {
207  LogDebug("NotFound") << "====> Value is Nan " << pVal << " == " << iConstits[i0].pt << " -- "
208  << iConstits[i0].eta << endl;
209  continue;
210  }
211  }
212  }
213 }
214 
215 int PuppiContainer::getPuppiId(float iPt, float iEta) {
216  int lId = -1;
217  for (int i0 = 0; i0 < fNAlgos; i0++) {
218  int nEtaBinsPerAlgo = fPuppiAlgo[i0].etaBins();
219  for (int i1 = 0; i1 < nEtaBinsPerAlgo; i1++) {
220  if ((std::abs(iEta) >= fPuppiAlgo[i0].etaMin(i1)) && (std::abs(iEta) < fPuppiAlgo[i0].etaMax(i1))) {
221  fPuppiAlgo[i0].fixAlgoEtaBin(i1);
222  if (iPt > fPuppiAlgo[i0].ptMin()) {
223  lId = i0;
224  break;
225  }
226  }
227  }
228  }
229  //if(lId == -1) std::cerr << "Error : Full fiducial range is not defined " << std::endl;
230  return lId;
231 }
232 double PuppiContainer::getChi2FromdZ(double iDZ) {
233  //We need to obtain prob of PU + (1-Prob of LV)
234  // Prob(LV) = Gaus(dZ,sigma) where sigma = 1.5mm (its really more like 1mm)
235  //double lProbLV = ROOT::Math::normal_cdf_c(std::abs(iDZ),0.2)*2.; //*2 is to do it double sided
236  //Take iDZ to be corrected by sigma already
237  double lProbLV = ROOT::Math::normal_cdf_c(std::abs(iDZ), 1.) * 2.; //*2 is to do it double sided
238  double lProbPU = 1 - lProbLV;
239  if (lProbPU <= 0)
240  lProbPU = 1e-16; //Quick Trick to through out infs
241  if (lProbPU >= 0)
242  lProbPU = 1 - 1e-16; //Ditto
243  double lChi2PU = TMath::ChisquareQuantile(lProbPU, 1);
244  lChi2PU *= lChi2PU;
245  return lChi2PU;
246 }
247 std::vector<double> const &PuppiContainer::puppiWeights() {
248  int lNParticles = fRecoParticles->size();
249 
250  fWeights.clear();
251  fWeights.reserve(lNParticles);
252  fVals.clear();
253  fVals.reserve(lNParticles);
254  for (int i0 = 0; i0 < fNAlgos; i0++)
255  fPuppiAlgo[i0].reset();
256 
257  int lNMaxAlgo = 1;
258  for (int i0 = 0; i0 < fNAlgos; i0++)
259  lNMaxAlgo = std::max(fPuppiAlgo[i0].numAlgos(), lNMaxAlgo);
260  //Run through all compute mean and RMS
261  for (int i0 = 0; i0 < lNMaxAlgo; i0++) {
262  getRMSAvg(i0, fPFParticles, fPFParticlesForVar, fPFParticlesForVarChargedPV);
263  }
264  if (fPuppiDiagnostics)
265  getRawAlphas(0, fPFParticles, fPFParticlesForVar, fPFParticlesForVarChargedPV);
266 
267  std::vector<double> pVals;
268  pVals.reserve(lNParticles);
269  for (int i0 = 0; i0 < lNParticles; i0++) {
270  //Refresh
271  pVals.clear();
272  double pWeight = 1;
273  //Get the Puppi Id and if ill defined move on
274  const auto &rParticle = (*fRecoParticles)[i0];
275  int pPupId = getPuppiId(rParticle.pt, rParticle.eta);
276  if (pPupId == -1) {
277  fWeights.push_back(0);
278  fAlphaMed.push_back(-10);
279  fAlphaRMS.push_back(-10);
280  continue;
281  }
282 
283  // fill the p-values
284  double pChi2 = 0;
285  if (fUseExp) {
286  //Compute an Experimental Puppi Weight with delta Z info (very simple example)
287  pChi2 = getChi2FromdZ(rParticle.dZ);
288  //Now make sure Neutrals are not set
289  if ((std::abs(rParticle.pdgId) == 22) || (std::abs(rParticle.pdgId) == 130))
290  pChi2 = 0;
291  }
292  //Fill and compute the PuppiWeight
293  int lNAlgos = fPuppiAlgo[pPupId].numAlgos();
294  for (int i1 = 0; i1 < lNAlgos; i1++)
295  pVals.push_back(fVals[lNParticles * i1 + i0]);
296 
297  pWeight = fPuppiAlgo[pPupId].compute(pVals, pChi2);
298  //Apply the CHS weights
299  if (rParticle.id == 1 && fApplyCHS)
300  pWeight = 1;
301  if (rParticle.id == 2 && fApplyCHS)
302  pWeight = 0;
303  //Apply weight of 1 for leptons if puppiNoLep
304  if (rParticle.id == 3)
305  pWeight = 1;
306  //Basic Weight Checks
307  if (!edm::isFinite(pWeight)) {
308  pWeight = 0.0;
309  LogDebug("PuppiWeightError") << "====> Weight is nan : " << pWeight << " : pt " << rParticle.pt
310  << " -- eta : " << rParticle.eta << " -- Value" << fVals[i0]
311  << " -- id : " << rParticle.id << " -- NAlgos: " << lNAlgos << std::endl;
312  }
313  //Basic Cuts
314  if (pWeight * fPFParticles[i0].pt < fPuppiAlgo[pPupId].neutralPt(fPUProxy) && rParticle.id == 0)
315  pWeight = 0; //threshold cut on the neutral Pt
316  // Protect high pT photons (important for gamma to hadronic recoil balance)
317  if (fPtMaxPhotons > 0 && rParticle.pdgId == 22 && std::abs(fPFParticles[i0].eta) < fEtaMaxPhotons &&
318  fPFParticles[i0].pt > fPtMaxPhotons)
319  pWeight = 1.;
320  // Protect high pT neutrals
321  else if ((fPtMaxNeutrals > 0) && (rParticle.id == 0))
322  pWeight = std::clamp(
323  (fPFParticles[i0].pt - fPtMaxNeutralsStartSlope) / (fPtMaxNeutrals - fPtMaxNeutralsStartSlope), pWeight, 1.);
324  if (pWeight < fPuppiWeightCut)
325  pWeight = 0; //==> Elminate the low Weight stuff
326  if (fInvert)
327  pWeight = 1. - pWeight;
328  //std::cout << "rParticle.pt = " << rParticle.pt << ", rParticle.charge = " << rParticle.charge << ", rParticle.id = " << rParticle.id << ", weight = " << pWeight << std::endl;
329 
330  fWeights.push_back(pWeight);
331  fAlphaMed.push_back(fPuppiAlgo[pPupId].median());
332  fAlphaRMS.push_back(fPuppiAlgo[pPupId].rms());
333  //Now get rid of the thrown out weights for the particle collection
334 
335  // leave these lines in, in case want to move eventually to having no 1-to-1 correspondence between puppi and pf cands
336  // if( std::abs(pWeight) < std::numeric_limits<double>::denorm_min() ) continue; // this line seems not to work like it's supposed to...
337  // if(std::abs(pWeight) <= 0. ) continue;
338  }
339  return fWeights;
340 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
double var_within_R(int iId, const std::vector< PuppiCandidate > &particles, const PuppiCandidate &centre, const double R)
double getChi2FromdZ(double iDZ)
constexpr float ptMin
constexpr bool isFinite(T x)
void getRMSAvg(int iOpt, std::vector< PuppiCandidate > const &iConstits, std::vector< PuppiCandidate > const &iParticles, std::vector< PuppiCandidate > const &iChargeParticles)
void getRawAlphas(int iOpt, std::vector< PuppiCandidate > const &iConstits, std::vector< PuppiCandidate > const &iParticles, std::vector< PuppiCandidate > const &iChargeParticles)
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
PuppiContainer(const edm::ParameterSet &iConfig)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< double > const & puppiWeights()
int getPuppiId(float iPt, float iEta)
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
T median(std::vector< T > values)
Definition: median.h:16
double goodVar(PuppiCandidate const &iPart, std::vector< PuppiCandidate > const &iParts, int iOpt, const double iRCone)
void initialize(const std::vector< RecoObj > &iRecoObjects)
void reset(double vett[256])
Definition: TPedValues.cc:11
#define LogDebug(id)