CMS 3D CMS Logo

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