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  std::vector<edm::ParameterSet> lAlgos = iConfig.getParameter<std::vector<edm::ParameterSet> >("algos");
22  fNAlgos = lAlgos.size();
23  for (unsigned int i0 = 0; i0 < lAlgos.size(); i0++) {
24  PuppiAlgo pPuppiConfig(lAlgos[i0]);
25  fPuppiAlgo.push_back(pPuppiConfig);
26  }
27 }
28 
29 void PuppiContainer::initialize(const std::vector<RecoObj> &iRecoObjects) {
30  //Clear everything
31  fPFParticles.resize(0);
32  fChargedPV.resize(0);
33  fWeights.resize(0);
34  fVals.resize(0);
35  fRawAlphas.resize(0);
36  fAlphaMed.resize(0);
37  fAlphaRMS.resize(0);
38  //fChargedNoPV.resize(0);
39  //Link to the RecoObjects
40  fPVFrac = 0.;
41  fNPV = 1.;
42  fRecoParticles = &iRecoObjects;
43  for (auto const &rParticle : *fRecoParticles) {
44  PuppiCandidate curPseudoJet;
45  // float nom = sqrt((rParticle.m)*(rParticle.m) + (rParticle.pt)*(rParticle.pt)*(cosh(rParticle.eta))*(cosh(rParticle.eta))) + (rParticle.pt)*sinh(rParticle.eta);//hacked
46  // float denom = sqrt((rParticle.m)*(rParticle.m) + (rParticle.pt)*(rParticle.pt));//hacked
47  // float rapidity = log(nom/denom);//hacked
48  if (edm::isFinite(rParticle.rapidity)) {
49  curPseudoJet.reset_PtYPhiM(rParticle.pt, rParticle.rapidity, rParticle.phi, rParticle.m); //hacked
50  } else {
51  curPseudoJet.reset_PtYPhiM(0, 99., 0, 0); //skipping may have been a better choice
52  }
53  //curPseudoJet.reset_PtYPhiM(rParticle.pt,rParticle.eta,rParticle.phi,rParticle.m);
54  // fill puppi_register
55  curPseudoJet.set_info(rParticle.id);
56  // fill vector of pseudojets for internal references
57  fPFParticles.push_back(curPseudoJet);
58  //Take Charged particles associated to PV
59  if (std::abs(rParticle.id) == 1)
60  fChargedPV.push_back(curPseudoJet);
61  //if(rParticle.id == 3) _chargedNoPV.push_back(curPseudoJet);
62  // if(fNPV < rParticle.vtxId) fNPV = rParticle.vtxId;
63  }
64 }
66 
68  std::vector<PuppiCandidate> const &iParts,
69  int iOpt,
70  const double iRCone) {
71  return var_within_R(iOpt, iParts, iPart, iRCone);
72 }
73 
75  const vector<PuppiCandidate> &particles,
76  const PuppiCandidate &centre,
77  const double R) {
78  if (iId == -1)
79  return 1;
80 
81  //this is a circle in rapidity-phi
82  //it would make more sense to have var definition consistent
83  //fastjet::Selector sel = fastjet::SelectorCircle(R);
84  //sel.set_reference(centre);
85  //the original code used Selector infrastructure: it is too heavy here
86  //logic of SelectorCircle is preserved below
87 
88  vector<double> near_dR2s;
89  near_dR2s.reserve(std::min(50UL, particles.size()));
90  vector<double> near_pts;
91  near_pts.reserve(std::min(50UL, particles.size()));
92  const double r2 = R * R;
93  for (auto const &part : particles) {
94  //squared_distance is in (y,phi) coords: rap() has faster access -> check it first
95  if (std::abs(part.rap() - centre.rap()) < R && part.squared_distance(centre) < r2) {
96  near_dR2s.push_back(reco::deltaR2(part, centre));
97  near_pts.push_back(part.pt());
98  }
99  }
100  double var = 0;
101  //double lSumPt = 0;
102  //if(iId == 1) for(auto pt : near_pts) lSumPt += pt;
103  auto nParts = near_dR2s.size();
104  for (auto i = 0UL; i < nParts; ++i) {
105  auto dr2 = near_dR2s[i];
106  auto pt = near_pts[i];
107  if (dr2 < 0.0001)
108  continue;
109  if (iId == 0)
110  var += (pt / dr2);
111  else if (iId == 1)
112  var += pt;
113  else if (iId == 2)
114  var += (1. / dr2);
115  else if (iId == 3)
116  var += (1. / dr2);
117  else if (iId == 4)
118  var += pt;
119  else if (iId == 5)
120  var += (pt * pt / dr2);
121  }
122  if (iId == 1)
123  var += centre.pt(); //Sum in a cone
124  else if (iId == 0 && var != 0)
125  var = log(var);
126  else if (iId == 3 && var != 0)
127  var = log(var);
128  else if (iId == 5 && var != 0)
129  var = log(var);
130  return var;
131 }
132 //In fact takes the median not the average
134  std::vector<PuppiCandidate> const &iConstits,
135  std::vector<PuppiCandidate> const &iParticles,
136  std::vector<PuppiCandidate> const &iChargedParticles) {
137  for (unsigned int i0 = 0; i0 < iConstits.size(); i0++) {
138  double pVal = -1;
139  //Calculate the Puppi Algo to use
140  int pPupId = getPuppiId(iConstits[i0].pt(), iConstits[i0].eta());
141  if (pPupId == -1 || fPuppiAlgo[pPupId].numAlgos() <= iOpt) {
142  fVals.push_back(-1);
143  continue;
144  }
145  //Get the Puppi Sub Algo (given iteration)
146  int pAlgo = fPuppiAlgo[pPupId].algoId(iOpt);
147  bool pCharged = fPuppiAlgo[pPupId].isCharged(iOpt);
148  double pCone = fPuppiAlgo[pPupId].coneSize(iOpt);
149  //Compute the Puppi Metric
150  if (!pCharged)
151  pVal = goodVar(iConstits[i0], iParticles, pAlgo, pCone);
152  if (pCharged)
153  pVal = goodVar(iConstits[i0], iChargedParticles, pAlgo, pCone);
154  fVals.push_back(pVal);
155  //if(std::isnan(pVal) || std::isinf(pVal)) cerr << "====> Value is Nan " << pVal << " == " << iConstits[i0].pt() << " -- " << iConstits[i0].eta() << endl;
156  if (!edm::isFinite(pVal)) {
157  LogDebug("NotFound") << "====> Value is Nan " << pVal << " == " << iConstits[i0].pt() << " -- "
158  << iConstits[i0].eta() << endl;
159  continue;
160  }
161 
162  // // fPuppiAlgo[pPupId].add(iConstits[i0],pVal,iOpt);
163  //code added by Nhan, now instead for every algorithm give it all the particles
164  for (int i1 = 0; i1 < fNAlgos; i1++) {
165  pAlgo = fPuppiAlgo[i1].algoId(iOpt);
166  pCharged = fPuppiAlgo[i1].isCharged(iOpt);
167  pCone = fPuppiAlgo[i1].coneSize(iOpt);
168  double curVal = -1;
169  if (i1 != pPupId) {
170  if (!pCharged)
171  curVal = goodVar(iConstits[i0], iParticles, pAlgo, pCone);
172  if (pCharged)
173  curVal = goodVar(iConstits[i0], iChargedParticles, pAlgo, pCone);
174  } else { //no need to repeat the computation
175  curVal = pVal;
176  }
177  //std::cout << "i1 = " << i1 << ", curVal = " << curVal << ", eta = " << iConstits[i0].eta() << ", pupID = " << pPupId << std::endl;
178  fPuppiAlgo[i1].add(iConstits[i0], curVal, iOpt);
179  }
180  }
181  for (int i0 = 0; i0 < fNAlgos; i0++)
182  fPuppiAlgo[i0].computeMedRMS(iOpt);
183 }
184 //In fact takes the median not the average
186  std::vector<PuppiCandidate> const &iConstits,
187  std::vector<PuppiCandidate> const &iParticles,
188  std::vector<PuppiCandidate> const &iChargedParticles) {
189  for (int j0 = 0; j0 < fNAlgos; j0++) {
190  for (unsigned int i0 = 0; i0 < iConstits.size(); i0++) {
191  double pVal = -1;
192  //Get the Puppi Sub Algo (given iteration)
193  int pAlgo = fPuppiAlgo[j0].algoId(iOpt);
194  bool pCharged = fPuppiAlgo[j0].isCharged(iOpt);
195  double pCone = fPuppiAlgo[j0].coneSize(iOpt);
196  //Compute the Puppi Metric
197  if (!pCharged)
198  pVal = goodVar(iConstits[i0], iParticles, pAlgo, pCone);
199  if (pCharged)
200  pVal = goodVar(iConstits[i0], iChargedParticles, pAlgo, pCone);
201  fRawAlphas.push_back(pVal);
202  if (!edm::isFinite(pVal)) {
203  LogDebug("NotFound") << "====> Value is Nan " << pVal << " == " << iConstits[i0].pt() << " -- "
204  << iConstits[i0].eta() << endl;
205  continue;
206  }
207  }
208  }
209 }
210 int PuppiContainer::getPuppiId(float iPt, float iEta) {
211  int lId = -1;
212  for (int i0 = 0; i0 < fNAlgos; i0++) {
213  int nEtaBinsPerAlgo = fPuppiAlgo[i0].etaBins();
214  for (int i1 = 0; i1 < nEtaBinsPerAlgo; i1++) {
215  if ((std::abs(iEta) > fPuppiAlgo[i0].etaMin(i1)) && (std::abs(iEta) < fPuppiAlgo[i0].etaMax(i1))) {
216  fPuppiAlgo[i0].fixAlgoEtaBin(i1);
217  if (iPt > fPuppiAlgo[i0].ptMin()) {
218  lId = i0;
219  break;
220  }
221  }
222  }
223  }
224  //if(lId == -1) std::cerr << "Error : Full fiducial range is not defined " << std::endl;
225  return lId;
226 }
227 double PuppiContainer::getChi2FromdZ(double iDZ) {
228  //We need to obtain prob of PU + (1-Prob of LV)
229  // Prob(LV) = Gaus(dZ,sigma) where sigma = 1.5mm (its really more like 1mm)
230  //double lProbLV = ROOT::Math::normal_cdf_c(std::abs(iDZ),0.2)*2.; //*2 is to do it double sided
231  //Take iDZ to be corrected by sigma already
232  double lProbLV = ROOT::Math::normal_cdf_c(std::abs(iDZ), 1.) * 2.; //*2 is to do it double sided
233  double lProbPU = 1 - lProbLV;
234  if (lProbPU <= 0)
235  lProbPU = 1e-16; //Quick Trick to through out infs
236  if (lProbPU >= 0)
237  lProbPU = 1 - 1e-16; //Ditto
238  double lChi2PU = TMath::ChisquareQuantile(lProbPU, 1);
239  lChi2PU *= lChi2PU;
240  return lChi2PU;
241 }
242 std::vector<double> const &PuppiContainer::puppiWeights() {
243  int lNParticles = fRecoParticles->size();
244 
245  fWeights.clear();
246  fWeights.reserve(lNParticles);
247  fVals.clear();
248  fVals.reserve(lNParticles);
249  for (int i0 = 0; i0 < fNAlgos; i0++)
250  fPuppiAlgo[i0].reset();
251 
252  int lNMaxAlgo = 1;
253  for (int i0 = 0; i0 < fNAlgos; i0++)
254  lNMaxAlgo = std::max(fPuppiAlgo[i0].numAlgos(), lNMaxAlgo);
255  //Run through all compute mean and RMS
256  for (int i0 = 0; i0 < lNMaxAlgo; i0++) {
257  getRMSAvg(i0, fPFParticles, fPFParticles, fChargedPV);
258  }
259  if (fPuppiDiagnostics)
260  getRawAlphas(0, fPFParticles, fPFParticles, fChargedPV);
261 
262  std::vector<double> pVals;
263  pVals.reserve(lNParticles);
264  for (int i0 = 0; i0 < lNParticles; i0++) {
265  //Refresh
266  pVals.clear();
267  double pWeight = 1;
268  //Get the Puppi Id and if ill defined move on
269  const auto &rParticle = (*fRecoParticles)[i0];
270  int pPupId = getPuppiId(rParticle.pt, rParticle.eta);
271  if (pPupId == -1) {
272  fWeights.push_back(0);
273  fAlphaMed.push_back(-10);
274  fAlphaRMS.push_back(-10);
275  continue;
276  }
277 
278  // fill the p-values
279  double pChi2 = 0;
280  if (fUseExp) {
281  //Compute an Experimental Puppi Weight with delta Z info (very simple example)
282  pChi2 = getChi2FromdZ(rParticle.dZ);
283  //Now make sure Neutrals are not set
284  if ((std::abs(rParticle.pdgId) == 22) || (std::abs(rParticle.pdgId) == 130))
285  pChi2 = 0;
286  }
287  //Fill and compute the PuppiWeight
288  int lNAlgos = fPuppiAlgo[pPupId].numAlgos();
289  for (int i1 = 0; i1 < lNAlgos; i1++)
290  pVals.push_back(fVals[lNParticles * i1 + i0]);
291 
292  pWeight = fPuppiAlgo[pPupId].compute(pVals, pChi2);
293  //Apply the CHS weights
294  if (rParticle.id == 1 && fApplyCHS)
295  pWeight = 1;
296  if (rParticle.id == 2 && fApplyCHS)
297  pWeight = 0;
298  //Basic Weight Checks
299  if (!edm::isFinite(pWeight)) {
300  pWeight = 0.0;
301  LogDebug("PuppiWeightError") << "====> Weight is nan : " << pWeight << " : pt " << rParticle.pt
302  << " -- eta : " << rParticle.eta << " -- Value" << fVals[i0]
303  << " -- id : " << rParticle.id << " -- NAlgos: " << lNAlgos << std::endl;
304  }
305  //Basic Cuts
306  if (pWeight * fPFParticles[i0].pt() < fPuppiAlgo[pPupId].neutralPt(fNPV) && rParticle.id == 0)
307  pWeight = 0; //threshold cut on the neutral Pt
308  // Protect high pT photons (important for gamma to hadronic recoil balance)
309  if ((fPtMaxPhotons > 0) && (rParticle.pdgId == 22) && (std::abs(fPFParticles[i0].eta()) < fEtaMaxPhotons) &&
310  (fPFParticles[i0].pt() > fPtMaxPhotons))
311  pWeight = 1.;
312  // Protect high pT neutrals
313  else if ((fPtMaxNeutrals > 0) && (rParticle.id == 0))
314  pWeight = std::clamp(fPFParticles[i0].pt() / fPtMaxNeutrals, pWeight, 1.);
315  if (pWeight < fPuppiWeightCut)
316  pWeight = 0; //==> Elminate the low Weight stuff
317  if (fInvert)
318  pWeight = 1. - pWeight;
319  //std::cout << "rParticle.pt = " << rParticle.pt << ", rParticle.charge = " << rParticle.charge << ", rParticle.id = " << rParticle.id << ", weight = " << pWeight << std::endl;
320 
321  fWeights.push_back(pWeight);
322  fAlphaMed.push_back(fPuppiAlgo[pPupId].median());
323  fAlphaRMS.push_back(fPuppiAlgo[pPupId].rms());
324  //Now get rid of the thrown out weights for the particle collection
325 
326  // leave these lines in, in case want to move eventually to having no 1-to-1 correspondence between puppi and pf cands
327  // if( std::abs(pWeight) < std::numeric_limits<double>::denorm_min() ) continue; // this line seems not to work like it's supposed to...
328  // if(std::abs(pWeight) <= 0. ) continue;
329  }
330  return fWeights;
331 }
#define LogDebug(id)
T getParameter(std::string const &) const
std::vector< double > const & puppiWeights()
double var_within_R(int iId, const std::vector< PuppiCandidate > &particles, const PuppiCandidate &centre, const double R)
double getChi2FromdZ(double iDZ)
constexpr bool isFinite(T x)
void set_info(int puppi_register)
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)
PuppiContainer(const edm::ParameterSet &iConfig)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T min(T a, T b)
Definition: MathUtil.h:58
int getPuppiId(float iPt, float iEta)
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
part
Definition: HCALResponse.h:20
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