CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Member Functions | Protected Member Functions | Protected Attributes
PuppiContainer Class Reference

#include <PuppiContainer.h>

Classes

class  PuppiUserInfo
 

Public Member Functions

void initialize (const std::vector< RecoObj > &iRecoObjects)
 
std::vector
< fastjet::PseudoJet > const & 
pfParticles () const
 
const std::vector< double > & puppiAlphas ()
 
const std::vector< double > & puppiAlphasMed ()
 
const std::vector< double > & puppiAlphasRMS ()
 
 PuppiContainer (const edm::ParameterSet &iConfig)
 
int puppiNAlgos ()
 
std::vector
< fastjet::PseudoJet > const & 
puppiParticles () const
 
const std::vector< double > & puppiRawAlphas ()
 
std::vector< double > const & puppiWeights ()
 
std::vector
< fastjet::PseudoJet > const & 
pvParticles () const
 
void setNPV (int iNPV)
 
 ~PuppiContainer ()
 

Protected Member Functions

double getChi2FromdZ (double iDZ)
 
int getPuppiId (float iPt, float iEta)
 
void getRawAlphas (int iOpt, std::vector< fastjet::PseudoJet > const &iConstits, std::vector< fastjet::PseudoJet > const &iParticles, std::vector< fastjet::PseudoJet > const &iChargeParticles)
 
void getRMSAvg (int iOpt, std::vector< fastjet::PseudoJet > const &iConstits, std::vector< fastjet::PseudoJet > const &iParticles, std::vector< fastjet::PseudoJet > const &iChargeParticles)
 
double goodVar (fastjet::PseudoJet const &iPart, std::vector< fastjet::PseudoJet > const &iParts, int iOpt, const double iRCone)
 
double var_within_R (int iId, const std::vector< fastjet::PseudoJet > &particles, const fastjet::PseudoJet &centre, const double R)
 

Protected Attributes

std::vector< double > fAlphaMed
 
std::vector< double > fAlphaRMS
 
bool fApplyCHS
 
std::vector< fastjet::PseudoJet > fChargedPV
 
bool fInvert
 
int fNAlgos
 
double fNeutralMinPt
 
double fNeutralSlope
 
int fNPV
 
std::vector< fastjet::PseudoJet > fPFParticles
 
std::vector< fastjet::PseudoJet > fPupParticles
 
std::vector< PuppiAlgofPuppiAlgo
 
bool fPuppiDiagnostics
 
double fPuppiWeightCut
 
double fPVFrac
 
std::vector< double > fRawAlphas
 
std::vector< RecoObjfRecoParticles
 
bool fUseExp
 
std::vector< double > fVals
 
std::vector< double > fWeights
 

Detailed Description

Definition at line 13 of file PuppiContainer.h.

Constructor & Destructor Documentation

PuppiContainer::PuppiContainer ( const edm::ParameterSet iConfig)

Definition at line 15 of file PuppiContainer.cc.

References edm::ParameterSet::getParameter().

15  {
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  std::vector<edm::ParameterSet> lAlgos = iConfig.getParameter<std::vector<edm::ParameterSet> >("algos");
22  fNAlgos = lAlgos.size();
23  for(unsigned int i0 = 0; i0 < lAlgos.size(); i0++) {
24  PuppiAlgo pPuppiConfig(lAlgos[i0]);
25  fPuppiAlgo.push_back(pPuppiConfig);
26  }
27 }
T getParameter(std::string const &) const
double fPuppiWeightCut
std::vector< PuppiAlgo > fPuppiAlgo
PuppiContainer::~PuppiContainer ( )

Definition at line 75 of file PuppiContainer.cc.

75 {}

Member Function Documentation

double PuppiContainer::getChi2FromdZ ( double  iDZ)
protected

Definition at line 198 of file PuppiContainer.cc.

References funct::abs(), and alignCSCRings::e.

198  {
199  //We need to obtain prob of PU + (1-Prob of LV)
200  // Prob(LV) = Gaus(dZ,sigma) where sigma = 1.5mm (its really more like 1mm)
201  //double lProbLV = ROOT::Math::normal_cdf_c(std::abs(iDZ),0.2)*2.; //*2 is to do it double sided
202  //Take iDZ to be corrected by sigma already
203  double lProbLV = ROOT::Math::normal_cdf_c(std::abs(iDZ),1.)*2.; //*2 is to do it double sided
204  double lProbPU = 1-lProbLV;
205  if(lProbPU <= 0) lProbPU = 1e-16; //Quick Trick to through out infs
206  if(lProbPU >= 0) lProbPU = 1-1e-16; //Ditto
207  double lChi2PU = TMath::ChisquareQuantile(lProbPU,1);
208  lChi2PU*=lChi2PU;
209  return lChi2PU;
210 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int PuppiContainer::getPuppiId ( float  iPt,
float  iEta 
)
protected

Definition at line 181 of file PuppiContainer.cc.

References funct::abs(), HLT_FULL_cff::etaMax, HLT_FULL_cff::etaMin, and PtMinSelector_cfg::ptMin.

181  {
182  int lId = -1;
183  for(int i0 = 0; i0 < fNAlgos; i0++) {
184  int nEtaBinsPerAlgo = fPuppiAlgo[i0].etaBins();
185  for (int i1 = 0; i1 < nEtaBinsPerAlgo; i1++){
186  if ( (std::abs(iEta) > fPuppiAlgo[i0].etaMin(i1)) && (std::abs(iEta) < fPuppiAlgo[i0].etaMax(i1)) ){
187  fPuppiAlgo[i0].fixAlgoEtaBin( i1 );
188  if(iPt > fPuppiAlgo[i0].ptMin()){
189  lId = i0;
190  break;
191  }
192  }
193  }
194  }
195  //if(lId == -1) std::cerr << "Error : Full fiducial range is not defined " << std::endl;
196  return lId;
197 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< PuppiAlgo > fPuppiAlgo
void PuppiContainer::getRawAlphas ( int  iOpt,
std::vector< fastjet::PseudoJet > const &  iConstits,
std::vector< fastjet::PseudoJet > const &  iParticles,
std::vector< fastjet::PseudoJet > const &  iChargeParticles 
)
protected

Definition at line 162 of file PuppiContainer.cc.

References edm::isFinite(), and LogDebug.

162  {
163  for(int j0 = 0; j0 < fNAlgos; j0++){
164  for(unsigned int i0 = 0; i0 < iConstits.size(); i0++ ) {
165  double pVal = -1;
166  //Get the Puppi Sub Algo (given iteration)
167  int pAlgo = fPuppiAlgo[j0].algoId (iOpt);
168  bool pCharged = fPuppiAlgo[j0].isCharged(iOpt);
169  double pCone = fPuppiAlgo[j0].coneSize (iOpt);
170  //Compute the Puppi Metric
171  if(!pCharged) pVal = goodVar(iConstits[i0],iParticles ,pAlgo,pCone);
172  if( pCharged) pVal = goodVar(iConstits[i0],iChargedParticles,pAlgo,pCone);
173  fRawAlphas.push_back(pVal);
174  if( ! edm::isFinite(pVal)) {
175  LogDebug( "NotFound" ) << "====> Value is Nan " << pVal << " == " << iConstits[i0].pt() << " -- " << iConstits[i0].eta() << endl;
176  continue;
177  }
178  }
179  }
180 }
#define LogDebug(id)
bool isFinite(T x)
std::vector< double > fRawAlphas
double goodVar(fastjet::PseudoJet const &iPart, std::vector< fastjet::PseudoJet > const &iParts, int iOpt, const double iRCone)
std::vector< PuppiAlgo > fPuppiAlgo
void PuppiContainer::getRMSAvg ( int  iOpt,
std::vector< fastjet::PseudoJet > const &  iConstits,
std::vector< fastjet::PseudoJet > const &  iParticles,
std::vector< fastjet::PseudoJet > const &  iChargeParticles 
)
protected

Definition at line 122 of file PuppiContainer.cc.

References eta, edm::isFinite(), LogDebug, and EnergyCorrector::pt.

122  {
123  for(unsigned int i0 = 0; i0 < iConstits.size(); i0++ ) {
124  double pVal = -1;
125  //Calculate the Puppi Algo to use
126  int pPupId = getPuppiId(iConstits[i0].pt(),iConstits[i0].eta());
127  if(pPupId == -1 || fPuppiAlgo[pPupId].numAlgos() <= iOpt){
128  fVals.push_back(-1);
129  continue;
130  }
131  //Get the Puppi Sub Algo (given iteration)
132  int pAlgo = fPuppiAlgo[pPupId].algoId (iOpt);
133  bool pCharged = fPuppiAlgo[pPupId].isCharged(iOpt);
134  double pCone = fPuppiAlgo[pPupId].coneSize (iOpt);
135  //Compute the Puppi Metric
136  if(!pCharged) pVal = goodVar(iConstits[i0],iParticles ,pAlgo,pCone);
137  if( pCharged) pVal = goodVar(iConstits[i0],iChargedParticles,pAlgo,pCone);
138  fVals.push_back(pVal);
139  //if(std::isnan(pVal) || std::isinf(pVal)) cerr << "====> Value is Nan " << pVal << " == " << iConstits[i0].pt() << " -- " << iConstits[i0].eta() << endl;
140  if( ! edm::isFinite(pVal)) {
141  LogDebug( "NotFound" ) << "====> Value is Nan " << pVal << " == " << iConstits[i0].pt() << " -- " << iConstits[i0].eta() << endl;
142  continue;
143  }
144 
145  // // fPuppiAlgo[pPupId].add(iConstits[i0],pVal,iOpt);
146  //code added by Nhan, now instead for every algorithm give it all the particles
147  for(int i1 = 0; i1 < fNAlgos; i1++){
148  pAlgo = fPuppiAlgo[i1].algoId (iOpt);
149  pCharged = fPuppiAlgo[i1].isCharged(iOpt);
150  pCone = fPuppiAlgo[i1].coneSize (iOpt);
151  double curVal = -1;
152  if(!pCharged) curVal = goodVar(iConstits[i0],iParticles ,pAlgo,pCone);
153  if( pCharged) curVal = goodVar(iConstits[i0],iChargedParticles,pAlgo,pCone);
154  //std::cout << "i1 = " << i1 << ", curVal = " << curVal << ", eta = " << iConstits[i0].eta() << ", pupID = " << pPupId << std::endl;
155  fPuppiAlgo[i1].add(iConstits[i0],curVal,iOpt);
156  }
157 
158  }
159  for(int i0 = 0; i0 < fNAlgos; i0++) fPuppiAlgo[i0].computeMedRMS(iOpt,fPVFrac);
160 }
#define LogDebug(id)
std::vector< double > fVals
bool isFinite(T x)
int getPuppiId(float iPt, float iEta)
double goodVar(fastjet::PseudoJet const &iPart, std::vector< fastjet::PseudoJet > const &iParts, int iOpt, const double iRCone)
std::vector< PuppiAlgo > fPuppiAlgo
double PuppiContainer::goodVar ( fastjet::PseudoJet const &  iPart,
std::vector< fastjet::PseudoJet > const &  iParts,
int  iOpt,
const double  iRCone 
)
protected

Definition at line 77 of file PuppiContainer.cc.

77  {
78  double lPup = 0;
79  lPup = var_within_R(iOpt,iParts,iPart,iRCone);
80  return lPup;
81 }
double var_within_R(int iId, const std::vector< fastjet::PseudoJet > &particles, const fastjet::PseudoJet &centre, const double R)
void PuppiContainer::initialize ( const std::vector< RecoObj > &  iRecoObjects)

Definition at line 29 of file PuppiContainer.cc.

References funct::abs(), i, edm::isFinite(), and or.

29  {
30  //Clear everything
31  fRecoParticles.resize(0);
32  fPFParticles .resize(0);
33  fChargedPV .resize(0);
34  fPupParticles .resize(0);
35  fWeights .resize(0);
36  fVals.resize(0);
37  fRawAlphas.resize(0);
38  fAlphaMed .resize(0);
39  fAlphaRMS .resize(0);
40  //fChargedNoPV.resize(0);
41  //Link to the RecoObjects
42  fPVFrac = 0.;
43  fNPV = 1.;
44  fRecoParticles = iRecoObjects;
45  for (unsigned int i = 0; i < fRecoParticles.size(); i++){
46  fastjet::PseudoJet curPseudoJet;
47  auto fRecoParticle = fRecoParticles[i];
48  // float nom = sqrt((fRecoParticle.m)*(fRecoParticle.m) + (fRecoParticle.pt)*(fRecoParticle.pt)*(cosh(fRecoParticle.eta))*(cosh(fRecoParticle.eta))) + (fRecoParticle.pt)*sinh(fRecoParticle.eta);//hacked
49  // float denom = sqrt((fRecoParticle.m)*(fRecoParticle.m) + (fRecoParticle.pt)*(fRecoParticle.pt));//hacked
50  // float rapidity = log(nom/denom);//hacked
51  if (edm::isFinite(fRecoParticle.rapidity)){
52  curPseudoJet.reset_PtYPhiM(fRecoParticle.pt,fRecoParticle.rapidity,fRecoParticle.phi,fRecoParticle.m);//hacked
53  } else {
54  curPseudoJet.reset_PtYPhiM(0, 99., 0, 0);//skipping may have been a better choice
55  }
56  //curPseudoJet.reset_PtYPhiM(fRecoParticle.pt,fRecoParticle.eta,fRecoParticle.phi,fRecoParticle.m);
57  int puppi_register = 0;
58  if(fRecoParticle.id == 0 or fRecoParticle.charge == 0) puppi_register = 0; // zero is neutral hadron
59  if(fRecoParticle.id == 1 and fRecoParticle.charge != 0) puppi_register = fRecoParticle.charge; // from PV use the
60  if(fRecoParticle.id == 2 and fRecoParticle.charge != 0) puppi_register = fRecoParticle.charge+5; // from NPV use the charge as key +5 as key
61  curPseudoJet.set_user_info( new PuppiUserInfo( puppi_register ) );
62  // fill vector of pseudojets for internal references
63  fPFParticles.push_back(curPseudoJet);
64  //Take Charged particles associated to PV
65  if(std::abs(fRecoParticle.id) == 1) fChargedPV.push_back(curPseudoJet);
66  if(std::abs(fRecoParticle.id) >= 1 ) fPVFrac+=1.;
67  //if((fRecoParticle.id == 0) && (inParticles[i].id == 2)) _genParticles.push_back( curPseudoJet);
68  //if(fRecoParticle.id <= 2 && !(inParticles[i].pt < fNeutralMinE && fRecoParticle.id < 2)) _pfchsParticles.push_back(curPseudoJet);
69  //if(fRecoParticle.id == 3) _chargedNoPV.push_back(curPseudoJet);
70  // if(fNPV < fRecoParticle.vtxId) fNPV = fRecoParticle.vtxId;
71  }
72  if (fPVFrac != 0) fPVFrac = double(fChargedPV.size())/fPVFrac;
73  else fPVFrac = 0;
74 }
int i
Definition: DBlmapReader.cc:9
std::vector< fastjet::PseudoJet > fPupParticles
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::EventIDconst &, edm::Timestampconst & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
std::vector< fastjet::PseudoJet > fChargedPV
std::vector< double > fVals
std::vector< fastjet::PseudoJet > fPFParticles
bool isFinite(T x)
std::vector< double > fRawAlphas
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< double > fAlphaMed
std::vector< double > fWeights
std::vector< RecoObj > fRecoParticles
std::vector< double > fAlphaRMS
std::vector<fastjet::PseudoJet> const& PuppiContainer::pfParticles ( ) const
inline

Definition at line 42 of file PuppiContainer.h.

References fPFParticles.

42 { return fPFParticles; }
std::vector< fastjet::PseudoJet > fPFParticles
const std::vector<double>& PuppiContainer::puppiAlphas ( )
inline

Definition at line 46 of file PuppiContainer.h.

References fVals.

46 { return fVals; }
std::vector< double > fVals
const std::vector<double>& PuppiContainer::puppiAlphasMed ( )
inline

Definition at line 48 of file PuppiContainer.h.

References fAlphaMed.

48 {return fAlphaMed;}
std::vector< double > fAlphaMed
const std::vector<double>& PuppiContainer::puppiAlphasRMS ( )
inline

Definition at line 49 of file PuppiContainer.h.

References fAlphaRMS.

49 {return fAlphaRMS;}
std::vector< double > fAlphaRMS
int PuppiContainer::puppiNAlgos ( )
inline

Definition at line 51 of file PuppiContainer.h.

References fNAlgos.

51 { return fNAlgos; }
std::vector<fastjet::PseudoJet> const& PuppiContainer::puppiParticles ( ) const
inline

Definition at line 52 of file PuppiContainer.h.

References fPupParticles.

52 { return fPupParticles;}
std::vector< fastjet::PseudoJet > fPupParticles
const std::vector<double>& PuppiContainer::puppiRawAlphas ( )
inline

Definition at line 45 of file PuppiContainer.h.

References fRawAlphas.

45 { return fRawAlphas; }
std::vector< double > fRawAlphas
std::vector< double > const & PuppiContainer::puppiWeights ( )

Definition at line 211 of file PuppiContainer.cc.

References alignCSCRings::e, eta, edm::isFinite(), LogDebug, bookConverter::max, EnergyCorrector::pt, reset(), and plotscripts::rms().

211  {
212  fPupParticles .resize(0);
213  fWeights .resize(0);
214  fVals .resize(0);
215  for(int i0 = 0; i0 < fNAlgos; i0++) fPuppiAlgo[i0].reset();
216 
217  int lNMaxAlgo = 1;
218  for(int i0 = 0; i0 < fNAlgos; i0++) lNMaxAlgo = std::max(fPuppiAlgo[i0].numAlgos(),lNMaxAlgo);
219  //Run through all compute mean and RMS
220  int lNParticles = fRecoParticles.size();
221  for(int i0 = 0; i0 < lNMaxAlgo; i0++) {
223  }
225 
226  std::vector<double> pVals;
227  for(int i0 = 0; i0 < lNParticles; i0++) {
228  //Refresh
229  pVals.clear();
230  double pWeight = 1;
231  //Get the Puppi Id and if ill defined move on
232  int pPupId = getPuppiId(fRecoParticles[i0].pt,fRecoParticles[i0].eta);
233  if(pPupId == -1) {
234  fWeights .push_back(pWeight);
235  fAlphaMed.push_back(-10);
236  fAlphaRMS.push_back(-10);
237  continue;
238  }
239  // fill the p-values
240  double pChi2 = 0;
241  if(fUseExp){
242  //Compute an Experimental Puppi Weight with delta Z info (very simple example)
243  pChi2 = getChi2FromdZ(fRecoParticles[i0].dZ);
244  //Now make sure Neutrals are not set
245  if(fRecoParticles[i0].pfType > 3) pChi2 = 0;
246  }
247  //Fill and compute the PuppiWeight
248  int lNAlgos = fPuppiAlgo[pPupId].numAlgos();
249  for(int i1 = 0; i1 < lNAlgos; i1++) pVals.push_back(fVals[lNParticles*i1+i0]);
250 
251  pWeight = fPuppiAlgo[pPupId].compute(pVals,pChi2);
252  //Apply the CHS weights
253  if(fRecoParticles[i0].id == 1 && fApplyCHS ) pWeight = 1;
254  if(fRecoParticles[i0].id == 2 && fApplyCHS ) pWeight = 0;
255  //Basic Weight Checks
256  if( ! edm::isFinite(pWeight)) {
257  pWeight = 0.0;
258  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;
259  }
260  //Basic Cuts
261  if(pWeight < fPuppiWeightCut) pWeight = 0; //==> Elminate the low Weight stuff
262  if(pWeight*fPFParticles[i0].pt() < fPuppiAlgo[pPupId].neutralPt(fNPV) && fRecoParticles[i0].id == 0 ) pWeight = 0; //threshold cut on the neutral Pt
263  if(fInvert) pWeight = 1.-pWeight;
264  //std::cout << "fRecoParticles[i0].pt = " << fRecoParticles[i0].pt << ", fRecoParticles[i0].charge = " << fRecoParticles[i0].charge << ", fRecoParticles[i0].id = " << fRecoParticles[i0].id << ", weight = " << pWeight << std::endl;
265 
266  fWeights .push_back(pWeight);
267  fAlphaMed.push_back(fPuppiAlgo[pPupId].median());
268  fAlphaRMS.push_back(fPuppiAlgo[pPupId].rms());
269  //Now get rid of the thrown out weights for the particle collection
270 
271  // leave these lines in, in case want to move eventually to having no 1-to-1 correspondence between puppi and pf cands
272  // if( std::abs(pWeight) < std::numeric_limits<double>::denorm_min() ) continue; // this line seems not to work like it's supposed to...
273  // if(std::abs(pWeight) <= 0. ) continue;
274 
275  //Produce
276  PseudoJet curjet( pWeight*fPFParticles[i0].px(), pWeight*fPFParticles[i0].py(), pWeight*fPFParticles[i0].pz(), pWeight*fPFParticles[i0].e() );
277  curjet.set_user_index(i0);
278  fPupParticles.push_back(curjet);
279  }
280  return fWeights;
281 }
#define LogDebug(id)
std::vector< fastjet::PseudoJet > fPupParticles
double getChi2FromdZ(double iDZ)
std::vector< fastjet::PseudoJet > fChargedPV
std::vector< double > fVals
std::vector< fastjet::PseudoJet > fPFParticles
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)
void getRawAlphas(int iOpt, std::vector< fastjet::PseudoJet > const &iConstits, std::vector< fastjet::PseudoJet > const &iParticles, std::vector< fastjet::PseudoJet > const &iChargeParticles)
std::vector< double > fAlphaMed
double fPuppiWeightCut
int getPuppiId(float iPt, float iEta)
std::vector< double > fWeights
std::vector< RecoObj > fRecoParticles
std::vector< double > fAlphaRMS
std::vector< PuppiAlgo > fPuppiAlgo
void reset(double vett[256])
Definition: TPedValues.cc:11
std::vector<fastjet::PseudoJet> const& PuppiContainer::pvParticles ( ) const
inline

Definition at line 43 of file PuppiContainer.h.

References fChargedPV.

43 { return fChargedPV; }
std::vector< fastjet::PseudoJet > fChargedPV
void PuppiContainer::setNPV ( int  iNPV)
inline

Definition at line 40 of file PuppiContainer.h.

References fNPV.

40 { fNPV = iNPV; }
double PuppiContainer::var_within_R ( int  iId,
const std::vector< fastjet::PseudoJet > &  particles,
const fastjet::PseudoJet &  centre,
const double  R 
)
protected

Definition at line 82 of file PuppiContainer.cc.

References reco::deltaR2(), i, dqm-mbProfile::log, min(), EnergyCorrector::pt, dttmaxenums::R, and MetTreeProducer::var().

82  {
83  if(iId == -1) return 1;
84 
85  //this is a circle in rapidity-phi
86  //it would make more sense to have var definition consistent
87  //fastjet::Selector sel = fastjet::SelectorCircle(R);
88  //sel.set_reference(centre);
89  //the original code used Selector infrastructure: it is too heavy here
90  //logic of SelectorCircle is preserved below
91 
92  vector<double > near_dR2s; near_dR2s.reserve(std::min(50UL, particles.size()));
93  vector<double > near_pts; near_pts.reserve(std::min(50UL, particles.size()));
94  for (auto const& part : particles){
95  if ( part.squared_distance(centre) < R*R ){
96  near_dR2s.push_back(reco::deltaR2(part, centre));
97  near_pts.push_back(part.pt());
98  }
99  }
100  double var = 0;
101  //double lSumPt = 0;
102  //if(iId == 1) for(auto pt : near_pts) lSumPt += pt;
103  auto nParts = near_dR2s.size();
104  for(auto i = 0UL; i < nParts; ++i){
105  auto dr2 = near_dR2s[i];
106  auto pt = near_pts[i];
107  if(dr2 < 0.0001) continue;
108  if(iId == 0) var += (pt/dr2);
109  else if(iId == 1) var += pt;
110  else if(iId == 2) var += (1./dr2);
111  else if(iId == 3) var += (1./dr2);
112  else if(iId == 4) var += pt;
113  else if(iId == 5) var += (pt * pt/dr2);
114  }
115  if(iId == 1) var += centre.pt(); //Sum in a cone
116  else if(iId == 0 && var != 0) var = log(var);
117  else if(iId == 3 && var != 0) var = log(var);
118  else if(iId == 5 && var != 0) var = log(var);
119  return var;
120 }
int i
Definition: DBlmapReader.cc:9
T min(T a, T b)
Definition: MathUtil.h:58
part
Definition: HCALResponse.h:20
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36

Member Data Documentation

std::vector<double> PuppiContainer::fAlphaMed
protected

Definition at line 70 of file PuppiContainer.h.

Referenced by puppiAlphasMed().

std::vector<double> PuppiContainer::fAlphaRMS
protected

Definition at line 71 of file PuppiContainer.h.

Referenced by puppiAlphasRMS().

bool PuppiContainer::fApplyCHS
protected

Definition at line 73 of file PuppiContainer.h.

std::vector<fastjet::PseudoJet> PuppiContainer::fChargedPV
protected

Definition at line 65 of file PuppiContainer.h.

Referenced by pvParticles().

bool PuppiContainer::fInvert
protected

Definition at line 74 of file PuppiContainer.h.

int PuppiContainer::fNAlgos
protected

Definition at line 79 of file PuppiContainer.h.

Referenced by puppiNAlgos().

double PuppiContainer::fNeutralMinPt
protected

Definition at line 76 of file PuppiContainer.h.

double PuppiContainer::fNeutralSlope
protected

Definition at line 77 of file PuppiContainer.h.

int PuppiContainer::fNPV
protected

Definition at line 80 of file PuppiContainer.h.

Referenced by setNPV().

std::vector<fastjet::PseudoJet> PuppiContainer::fPFParticles
protected

Definition at line 64 of file PuppiContainer.h.

Referenced by pfParticles().

std::vector<fastjet::PseudoJet> PuppiContainer::fPupParticles
protected

Definition at line 66 of file PuppiContainer.h.

Referenced by puppiParticles().

std::vector<PuppiAlgo> PuppiContainer::fPuppiAlgo
protected

Definition at line 82 of file PuppiContainer.h.

bool PuppiContainer::fPuppiDiagnostics
protected

Definition at line 62 of file PuppiContainer.h.

double PuppiContainer::fPuppiWeightCut
protected

Definition at line 78 of file PuppiContainer.h.

double PuppiContainer::fPVFrac
protected

Definition at line 81 of file PuppiContainer.h.

std::vector<double> PuppiContainer::fRawAlphas
protected

Definition at line 69 of file PuppiContainer.h.

Referenced by puppiRawAlphas().

std::vector<RecoObj> PuppiContainer::fRecoParticles
protected

Definition at line 63 of file PuppiContainer.h.

bool PuppiContainer::fUseExp
protected

Definition at line 75 of file PuppiContainer.h.

std::vector<double> PuppiContainer::fVals
protected

Definition at line 68 of file PuppiContainer.h.

Referenced by puppiAlphas().

std::vector<double> PuppiContainer::fWeights
protected

Definition at line 67 of file PuppiContainer.h.