CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PuppiAlgo.cc
Go to the documentation of this file.
4 #include "fastjet/internal/base.hh"
5 #include "Math/QuantFuncMathCore.h"
6 #include "Math/SpecFuncMathCore.h"
7 #include "Math/ProbFunc.h"
8 #include "TMath.h"
9 
10 
12  fEtaMin = iConfig.getParameter<double>("etaMin");
13  fEtaMax = iConfig.getParameter<double>("etaMax");
14  fPtMin = iConfig.getParameter<double>("ptMin");
15  fNeutralPtMin = iConfig.getParameter<double>("MinNeutralPt"); // Weighted Neutral Pt Cut
16  fNeutralPtSlope = iConfig.getParameter<double>("MinNeutralPtSlope"); // Slope vs #pv
17  fRMSEtaSF = iConfig.getParameter<double>("RMSEtaSF");
18  fMedEtaSF = iConfig.getParameter<double>("MedEtaSF");
19  fEtaMaxExtrap = iConfig.getParameter<double>("EtaMaxExtrap");
20 
21  std::vector<edm::ParameterSet> lAlgos = iConfig.getParameter<std::vector<edm::ParameterSet> >("puppiAlgos");
22  fNAlgos = lAlgos.size();
23  //Uber Configurable Puppi
24  for(unsigned int i0 = 0; i0 < lAlgos.size(); i0++) {
25  int pAlgoId = lAlgos[i0].getParameter<int > ("algoId");
26  bool pCharged = lAlgos[i0].getParameter<bool> ("useCharged");
27  bool pWeight0 = lAlgos[i0].getParameter<bool> ("applyLowPUCorr");
28  int pComb = lAlgos[i0].getParameter<int> ("combOpt"); // 0=> add in chi2/1=>Multiply p-values
29  double pConeSize = lAlgos[i0].getParameter<double>("cone"); // Min Pt when computing pt and rms
30  double pRMSPtMin = lAlgos[i0].getParameter<double>("rmsPtMin"); // Min Pt when computing pt and rms
31  double pRMSSF = lAlgos[i0].getParameter<double>("rmsScaleFactor"); // Additional Tuning parameter for Jokers
32  fAlgoId .push_back(pAlgoId);
33  fCharged .push_back(pCharged);
34  fAdjust .push_back(pWeight0);
35  fCombId .push_back(pComb);
36  fConeSize .push_back(pConeSize);
37  fRMSPtMin .push_back(pRMSPtMin);
38  fRMSScaleFactor.push_back(pRMSSF);
39  double pRMS = 0;
40  double pMed = 0;
41  double pMean = 0;
42  int pNCount = 0;
43  fRMS .push_back(pRMS);
44  fMedian.push_back(pMed);
45  fMean .push_back(pMean);
46  fNCount.push_back(pNCount);
47  }
48 }
50  fPups .clear();
51  fPupsPV.clear();
52 }
54  fPups .clear();
55  fPupsPV.clear();
56  for(unsigned int i0 = 0; i0 < fNAlgos; i0++) {
57  fMedian[i0] = 0;
58  fRMS [i0] = 0;
59  fMean [i0] = 0;
60  fNCount[i0] = 0;
61  }
62 }
63 void PuppiAlgo::add(const fastjet::PseudoJet &iParticle,const double &iVal,const unsigned int iAlgo) {
64  if(iParticle.pt() < fRMSPtMin[iAlgo]) return;
65  // Change from SRR : Previously used fastjet::PseudoJet::user_index to decide the particle type.
66  // In CMSSW we use the user_index to specify the index in the input collection, so I invented
67  // a new mechanism using the fastjet UserInfo functionality. Of course, it's still just an integer
68  // but that interface could be changed (or augmented) if desired / needed.
69  int puppi_register = std::numeric_limits<int>::lowest();
70  if ( iParticle.has_user_info() ) {
71  PuppiContainer::PuppiUserInfo const * pInfo = dynamic_cast<PuppiContainer::PuppiUserInfo const *>( iParticle.user_info_ptr() );
72  if ( pInfo != 0 ) {
73  puppi_register = pInfo->puppi_register();
74  }
75  }
76  if ( puppi_register == std::numeric_limits<int>::lowest() ) {
77  throw cms::Exception("PuppiRegisterNotSet") << "The puppi register is not set. This must be set before use.\n";
78  }
79 
81  // if(fCharged[iAlgo] && std::abs(puppi_register) < 1) return;
82  // if(fCharged[iAlgo] && (std::abs(puppi_register) >=1 && std::abs(puppi_register) <=2)) fPupsPV.push_back(iVal);
83  //if(fCharged[iAlgo] && std::abs(puppi_register) < 3) return;
85  // fPups.push_back(iVal); //original
86  // fNCount[iAlgo]++;
87 
88  // added by Nhan -- for all eta regions, compute mean/RMS from the central charged PU
89  //std::cout << "std::abs(puppi_register) = " << std::abs(puppi_register) << std::endl;
90  if ((std::abs(iParticle.eta()) < fEtaMaxExtrap) && (std::abs(puppi_register) >= 3)){
91  fPups.push_back(iVal);
92  // fPupsPV.push_back(iVal);
93  fNCount[iAlgo]++;
94  }
95  // for the low PU case, correction. for checking that the PU-only median will be below the PV particles
96  if(std::abs(iParticle.eta()) < fEtaMaxExtrap && (std::abs(puppi_register) >=1 && std::abs(puppi_register) <=2)) fPupsPV.push_back(iVal);
97 
98 }
99 
102 //NHAN'S VERSION
103 void PuppiAlgo::computeMedRMS(const unsigned int &iAlgo,const double &iPVFrac) {
104 
105  //std::cout << "fNCount[iAlgo] = " << fNCount[iAlgo] << std::endl;
106  if(iAlgo >= fNAlgos ) return;
107  if(fNCount[iAlgo] == 0) return;
108 
109  // sort alphas
110  int lNBefore = 0;
111  for(unsigned int i0 = 0; i0 < iAlgo; i0++) lNBefore += fNCount[i0];
112  std::sort(fPups.begin()+lNBefore,fPups.begin()+lNBefore+fNCount[iAlgo]);
113 
114  // in case you have alphas == 0
115  int lNum0 = 0;
116  for(int i0 = lNBefore; i0 < lNBefore+fNCount[iAlgo]; i0++) {
117  if(fPups[i0] == 0) lNum0 = i0-lNBefore;
118  }
119 
120  // comput median, removed lCorr for now
121  int lNHalfway = lNBefore + lNum0 + int( double( fNCount[iAlgo]-lNum0 )*0.50);
122  fMedian[iAlgo] = fPups[lNHalfway];
123  double lMed = fMedian[iAlgo]; //Just to make the readability easier
124 
125  int lNRMS = 0;
126  for(int i0 = lNBefore; i0 < lNBefore+fNCount[iAlgo]; i0++) {
127  fMean[iAlgo] += fPups[i0];
128  if(fPups[i0] == 0) continue;
129  if(!fCharged[iAlgo] && fAdjust[iAlgo] && fPups[i0] > lMed) continue;
130  //if(fAdjust[iAlgo] && fPups[i0] > lMed) continue;
131  lNRMS++;
132  fRMS [iAlgo] += (fPups[i0]-lMed)*(fPups[i0]-lMed);
133  }
134  fMean[iAlgo]/=fNCount[iAlgo];
135  if(lNRMS > 0) fRMS [iAlgo]/=lNRMS;
136  if(fRMS[iAlgo] == 0) fRMS[iAlgo] = 1e-5;
137  // here is the raw RMS
138  fRMS [iAlgo] = sqrt(fRMS[iAlgo]);
139 
140  // some ways to do corrections to fRMS and fMedian
141  fRMS [iAlgo] *= fRMSScaleFactor[iAlgo];
142 
143  fRMS[iAlgo] *= fRMSEtaSF;
144  fMedian[iAlgo] *= fMedEtaSF;
145 
146  if(!fAdjust[iAlgo]) return;
147  //Adjust the p-value to correspond to the median
148  std::sort(fPupsPV.begin(),fPupsPV.end());
149  int lNPV = 0;
150  for(unsigned int i0 = 0; i0 < fPupsPV.size(); i0++) if(fPupsPV[i0] <= lMed ) lNPV++;
151  double lAdjust = double(lNPV)/double(fPupsPV.size()+fNCount[iAlgo]);
152  if(lAdjust > 0) fMedian[iAlgo] -= sqrt(ROOT::Math::chisquared_quantile(lAdjust,1.)*fRMS[iAlgo]);
153 
154 }
156 
157 //This code is probably a bit confusing
158 double PuppiAlgo::compute(std::vector<double> const &iVals,double iChi2) const {
159  if(fAlgoId[0] == -1) return 1;
160  double lVal = 0.;
161  double lPVal = 1.;
162  int lNDOF = 0;
163  for(unsigned int i0 = 0; i0 < fNAlgos; i0++) {
164  if(fNCount[i0] == 0) return 1.; //in the NoPU case return 1.
165  if(fCombId[i0] == 1 && i0 > 0) { //Compute the previous p-value so that p-values can be multiplieed
166  double pPVal = ROOT::Math::chisquared_cdf(lVal,lNDOF);
167  lPVal *= pPVal;
168  lNDOF = 0;
169  lVal = 0;
170  }
171  double pVal = iVals[i0];
172  //Special Check for any algo with log(0)
173  if(fAlgoId[i0] == 0 && iVals[i0] == 0) pVal = fMedian[i0];
174  if(fAlgoId[i0] == 3 && iVals[i0] == 0) pVal = fMedian[i0];
175  if(fAlgoId[i0] == 5 && iVals[i0] == 0) pVal = fMedian[i0];
176  lVal += (pVal-fMedian[i0])*(fabs(pVal-fMedian[i0]))/fRMS[i0]/fRMS[i0];
177  lNDOF++;
178  if(i0 == 0 && iChi2 != 0) lNDOF++; //Add external Chi2 to first element
179  if(i0 == 0 && iChi2 != 0) lVal+=iChi2; //Add external Chi2 to first element
180  }
181  //Top it off with the last calc
182  lPVal *= ROOT::Math::chisquared_cdf(lVal,lNDOF);
183  return lPVal;
184 }
T getParameter(std::string const &) const
void add(const fastjet::PseudoJet &iParticle, const double &iVal, const unsigned int iAlgo)
Definition: PuppiAlgo.cc:63
double fMedEtaSF
Definition: PuppiAlgo.h:42
double fNeutralPtMin
Definition: PuppiAlgo.h:38
float fPtMin
Definition: PuppiAlgo.h:37
std::vector< bool > fAdjust
Definition: PuppiAlgo.h:52
PuppiAlgo(edm::ParameterSet &iConfig)
Definition: PuppiAlgo.cc:11
double fRMSEtaSF
Definition: PuppiAlgo.h:41
float fEtaMin
Definition: PuppiAlgo.h:36
std::vector< double > fConeSize
Definition: PuppiAlgo.h:54
void reset()
Definition: PuppiAlgo.cc:53
T sqrt(T t)
Definition: SSEVec.h:48
~PuppiAlgo()
Definition: PuppiAlgo.cc:49
std::vector< double > fRMS
Definition: PuppiAlgo.h:45
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< float > fPupsPV
Definition: PuppiAlgo.h:49
double compute(std::vector< double > const &iVals, double iChi2) const
Definition: PuppiAlgo.cc:158
std::vector< float > fPups
Definition: PuppiAlgo.h:48
std::vector< double > fMean
Definition: PuppiAlgo.h:57
std::vector< int > fAlgoId
Definition: PuppiAlgo.h:50
std::vector< double > fRMSPtMin
Definition: PuppiAlgo.h:55
std::vector< int > fCombId
Definition: PuppiAlgo.h:53
std::vector< double > fRMSScaleFactor
Definition: PuppiAlgo.h:56
void computeMedRMS(const unsigned int &iAlgo, const double &iPVFrac)
Definition: PuppiAlgo.cc:103
double fEtaMaxExtrap
Definition: PuppiAlgo.h:43
double fNeutralPtSlope
Definition: PuppiAlgo.h:39
std::vector< bool > fCharged
Definition: PuppiAlgo.h:51
unsigned int fNAlgos
Definition: PuppiAlgo.h:34
float fEtaMax
Definition: PuppiAlgo.h:35
std::vector< int > fNCount
Definition: PuppiAlgo.h:58
std::vector< double > fMedian
Definition: PuppiAlgo.h:46