CMS 3D CMS Logo

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