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