CMS 3D CMS Logo

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