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  fApplyCHS = iConfig.getParameter<bool>("applyCHS");
16  fUseExp = iConfig.getParameter<bool>("useExp");
17  fPuppiWeightCut = iConfig.getParameter<double>("MinPuppiWeight");
18  std::vector<edm::ParameterSet> lAlgos = iConfig.getParameter<std::vector<edm::ParameterSet> >("algos");
19  fNAlgos = lAlgos.size();
20  for(unsigned int i0 = 0; i0 < lAlgos.size(); i0++) {
21  PuppiAlgo pPuppiConfig(lAlgos[i0]);
22  fPuppiAlgo.push_back(pPuppiConfig);
23  }
24 }
25 
26 void PuppiContainer::initialize(const std::vector<RecoObj> &iRecoObjects) {
27  //Clear everything
28  fRecoParticles.resize(0);
29  fPFParticles .resize(0);
30  fChargedPV .resize(0);
31  fPupParticles .resize(0);
32  fWeights .resize(0);
33  fVals.resize(0);
34  //fChargedNoPV.resize(0);
35  //Link to the RecoObjects
36  fPVFrac = 0.;
37  fNPV = 1.;
38  fRecoParticles = iRecoObjects;
39  for (unsigned int i = 0; i < fRecoParticles.size(); i++){
40  fastjet::PseudoJet curPseudoJet;
41  auto fRecoParticle = fRecoParticles[i];
42  curPseudoJet.reset_PtYPhiM(fRecoParticle.pt,fRecoParticle.eta,fRecoParticle.phi,fRecoParticle.m);
43  int puppi_register = 0;
44  if(fRecoParticle.id == 0 or fRecoParticle.charge == 0) puppi_register = 0; // zero is neutral hadron
45  if(fRecoParticle.id == 1 and fRecoParticle.charge != 0) puppi_register = fRecoParticle.charge; // from PV use the
46  if(fRecoParticle.id == 2 and fRecoParticle.charge != 0) puppi_register = fRecoParticle.charge+5; // from NPV use the charge as key +5 as key
47  curPseudoJet.set_user_info( new PuppiUserInfo( puppi_register ) );
48  // fill vector of pseudojets for internal references
49  fPFParticles.push_back(curPseudoJet);
50  //Take Charged particles associated to PV
51  if(std::abs(fRecoParticle.id) == 1) fChargedPV.push_back(curPseudoJet);
52  if(std::abs(fRecoParticle.id) >= 1 ) fPVFrac+=1.;
53  //if((fRecoParticle.id == 0) && (inParticles[i].id == 2)) _genParticles.push_back( curPseudoJet);
54  //if(fRecoParticle.id <= 2 && !(inParticles[i].pt < fNeutralMinE && fRecoParticle.id < 2)) _pfchsParticles.push_back(curPseudoJet);
55  //if(fRecoParticle.id == 3) _chargedNoPV.push_back(curPseudoJet);
56  if(fNPV < fRecoParticle.vtxId) fNPV = fRecoParticle.vtxId;
57  }
58  fPVFrac = double(fChargedPV.size())/fPVFrac;
59 }
61 
62 double PuppiContainer::goodVar(PseudoJet const &iPart,std::vector<PseudoJet> const &iParts, int iOpt,double iRCone) {
63  double lPup = 0;
64  lPup = var_within_R(iOpt,iParts,iPart,iRCone);
65  return lPup;
66 }
67 double PuppiContainer::var_within_R(int iId, const vector<PseudoJet> & particles, const PseudoJet& centre, double R){
68  if(iId == -1) return 1;
69  fastjet::Selector sel = fastjet::SelectorCircle(R);
70  sel.set_reference(centre);
71  vector<PseudoJet> near_particles = sel(particles);
72  double var = 0;
73  //double lSumPt = 0;
74  //if(iId == 1) for(unsigned int i=0; i<near_particles.size(); i++) lSumPt += near_particles[i].pt();
75  for(unsigned int i=0; i<near_particles.size(); i++){
76  double pDEta = near_particles[i].eta()-centre.eta();
77  double pDPhi = std::abs(near_particles[i].phi()-centre.phi());
78  if(pDPhi > 2.*M_PI-pDPhi) pDPhi = 2.*M_PI-pDPhi;
79  double pDR2 = pDEta*pDEta+pDPhi*pDPhi;
80  if(std::abs(pDR2) < 0.0001) continue;
81  if(iId == 0) var += (near_particles[i].pt()/pDR2);
82  if(iId == 1) var += near_particles[i].pt();
83  if(iId == 2) var += (1./pDR2);
84  if(iId == 3) var += (1./pDR2);
85  if(iId == 4) var += near_particles[i].pt();
86  if(iId == 5) var += (near_particles[i].pt() * near_particles[i].pt()/pDR2);
87  }
88  if(iId == 1) var += centre.pt(); //Sum in a cone
89  if(iId == 0 && var != 0) var = log(var);
90  if(iId == 3 && var != 0) var = log(var);
91  if(iId == 5 && var != 0) var = log(var);
92  return var;
93 }
94 //In fact takes the median not the average
95 void PuppiContainer::getRMSAvg(int iOpt,std::vector<fastjet::PseudoJet> const &iConstits,std::vector<fastjet::PseudoJet> const &iParticles,std::vector<fastjet::PseudoJet> const &iChargedParticles) {
96  for(unsigned int i0 = 0; i0 < iConstits.size(); i0++ ) {
97  double pVal = -1;
98  //Calculate the Puppi Algo to use
99  int pPupId = getPuppiId(iConstits[i0].pt(),iConstits[i0].eta());
100  if(fPuppiAlgo[pPupId].numAlgos() <= iOpt) pPupId = -1;
101  if(pPupId == -1) {fVals.push_back(-1); continue;}
102  //Get the Puppi Sub Algo (given iteration)
103  int pAlgo = fPuppiAlgo[pPupId].algoId (iOpt);
104  bool pCharged = fPuppiAlgo[pPupId].isCharged(iOpt);
105  double pCone = fPuppiAlgo[pPupId].coneSize (iOpt);
106  //Compute the Puppi Metric
107  if(!pCharged) pVal = goodVar(iConstits[i0],iParticles ,pAlgo,pCone);
108  if( pCharged) pVal = goodVar(iConstits[i0],iChargedParticles,pAlgo,pCone);
109  fVals.push_back(pVal);
110  //if(std::isnan(pVal) || std::isinf(pVal)) cerr << "====> Value is Nan " << pVal << " == " << iConstits[i0].pt() << " -- " << iConstits[i0].eta() << endl;
111  if( ! edm::isFinite(pVal)) {
112  LogDebug( "NotFound" ) << "====> Value is Nan " << pVal << " == " << iConstits[i0].pt() << " -- " << iConstits[i0].eta() << endl;
113  continue;
114  }
115  fPuppiAlgo[pPupId].add(iConstits[i0],pVal,iOpt);
116  }
117  for(int i0 = 0; i0 < fNAlgos; i0++) fPuppiAlgo[i0].computeMedRMS(iOpt,fPVFrac);
118 }
119 int PuppiContainer::getPuppiId( float iPt, float iEta) {
120  int lId = -1;
121  for(int i0 = 0; i0 < fNAlgos; i0++) {
122  if(std::abs(iEta) < fPuppiAlgo[i0].etaMin()) continue;
123  if(std::abs(iEta) > fPuppiAlgo[i0].etaMax()) continue;
124  if(iPt < fPuppiAlgo[i0].ptMin()) continue;
125  lId = i0;
126  break;
127  }
128  //if(lId == -1) std::cerr << "Error : Full fiducial range is not defined " << std::endl;
129  return lId;
130 }
131 double PuppiContainer::getChi2FromdZ(double iDZ) {
132  //We need to obtain prob of PU + (1-Prob of LV)
133  // Prob(LV) = Gaus(dZ,sigma) where sigma = 1.5mm (its really more like 1mm)
134  //double lProbLV = ROOT::Math::normal_cdf_c(std::abs(iDZ),0.2)*2.; //*2 is to do it double sided
135  //Take iDZ to be corrected by sigma already
136  double lProbLV = ROOT::Math::normal_cdf_c(std::abs(iDZ),1.)*2.; //*2 is to do it double sided
137  double lProbPU = 1-lProbLV;
138  if(lProbPU <= 0) lProbPU = 1e-16; //Quick Trick to through out infs
139  if(lProbPU >= 0) lProbPU = 1-1e-16; //Ditto
140  double lChi2PU = TMath::ChisquareQuantile(lProbPU,1);
141  lChi2PU*=lChi2PU;
142  return lChi2PU;
143 }
144 std::vector<double> const & PuppiContainer::puppiWeights() {
145  fPupParticles .resize(0);
146  fWeights .resize(0);
147  fVals .resize(0);
148  for(int i0 = 0; i0 < fNAlgos; i0++) fPuppiAlgo[i0].reset();
149 
150  int lNMaxAlgo = 1;
151  for(int i0 = 0; i0 < fNAlgos; i0++) lNMaxAlgo = std::max(fPuppiAlgo[i0].numAlgos(),lNMaxAlgo);
152  //Run through all compute mean and RMS
153  int lNParticles = fRecoParticles.size();
154  for(int i0 = 0; i0 < lNMaxAlgo; i0++) {
155  getRMSAvg(i0,fPFParticles,fPFParticles,fChargedPV);
156  }
157  std::vector<double> pVals;
158  for(int i0 = 0; i0 < lNParticles; i0++) {
159  //Refresh
160  pVals.clear();
161  double pWeight = 1;
162  //Get the Puppi Id and if ill defined move on
163  int pPupId = getPuppiId(fRecoParticles[i0].pt,fRecoParticles[i0].eta);
164  if(pPupId == -1) {
165  fWeights .push_back(pWeight);
166  continue;
167  }
168  // fill the p-values
169  double pChi2 = 0;
170  if(fUseExp){
171  //Compute an Experimental Puppi Weight with delta Z info (very simple example)
172  pChi2 = getChi2FromdZ(fRecoParticles[i0].dZ);
173  //Now make sure Neutrals are not set
174  if(fRecoParticles[i0].pfType > 3) pChi2 = 0;
175  }
176  //Fill and compute the PuppiWeight
177  int lNAlgos = fPuppiAlgo[pPupId].numAlgos();
178  for(int i1 = 0; i1 < lNAlgos; i1++) pVals.push_back(fVals[lNParticles*i1+i0]);
179  pWeight = fPuppiAlgo[pPupId].compute(pVals,pChi2);
180  //Apply the CHS weights
181  if(fRecoParticles[i0].id == 1 && fApplyCHS ) pWeight = 1;
182  if(fRecoParticles[i0].id == 2 && fApplyCHS ) pWeight = 0;
183  //Basic Weight Checks
184  if( ! edm::isFinite(pWeight)) {
185  pWeight = 0.0;
186  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;
187  }
188  //Basic Cuts
189  if(pWeight < fPuppiWeightCut) pWeight = 0; //==> Elminate the low Weight stuff
190  if(pWeight*fPFParticles[i0].pt() < fPuppiAlgo[pPupId].neutralPt(fNPV) && fRecoParticles[i0].id == 0 ) pWeight = 0; //threshold cut on the neutral Pt
191  fWeights .push_back(pWeight);
192  //Now get rid of the thrown out weights for the particle collection
193  if(std::abs(pWeight) < std::numeric_limits<double>::denorm_min() ) continue;
194  //Produce
195  PseudoJet curjet( pWeight*fPFParticles[i0].px(), pWeight*fPFParticles[i0].py(), pWeight*fPFParticles[i0].pz(), pWeight*fPFParticles[i0].e());
196  curjet.set_user_index(i0);
197  fPupParticles.push_back(curjet);
198  }
199  return fWeights;
200 }
201 
202 
#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)
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
Definition: DDAxes.h:10