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