CMS 3D CMS Logo

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