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