CMS 3D CMS Logo

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