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