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