CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Types | Private Member Functions | Private Attributes
metsig::SignAlgoResolutions Class Reference

#include <SignAlgoResolutions.h>

Public Member Functions

void addResolutions (const edm::ParameterSet &iConfig)
 
double eval (const resolutionType &type, const resolutionFunc &func, const double &et, const double &phi, const double &eta, const double &p) const
 
double eval (const resolutionType &type, const resolutionFunc &func, const double &et, const double &phi, const double &eta) const
 
metsig::SigInputObj evalPF (const reco::PFCandidate *candidate) const
 
metsig::SigInputObj evalPFJet (const reco::PFJet *jet) const
 
bool isFilled () const
 
 SignAlgoResolutions ()
 
 SignAlgoResolutions (const edm::ParameterSet &iConfig)
 

Private Types

typedef std::pair
< metsig::resolutionType,
metsig::resolutionFunc
functionCombo
 
typedef std::vector< double > functionPars
 

Private Member Functions

void addfunction (const resolutionType type, const resolutionFunc func, std::vector< double > parameters)
 
double ElectronPtResolution (const reco::PFCandidate *c) const
 
double EtFunction (const functionPars &x, const functionPars &par) const
 
double getfunc (const resolutionType &type, const resolutionFunc &func, std::vector< double > &x) const
 
void initializeJetResolutions (const edm::ParameterSet &iConfig)
 
double PFunction (const functionPars &x, const functionPars &par) const
 
double PhiConstFunction (const functionPars &x, const functionPars &par) const
 
double PhiFunction (const functionPars &x, const functionPars &par) const
 

Private Attributes

std::map< functionCombo,
functionPars
functionmap_
 
std::vector< double > jdphi [10]
 
std::vector< double > jdpt [10]
 
PFEnergyResolutionpfresol_
 
JetResolutionphiResol_
 
JetResolutionptResol_
 
double ptResolThreshold_
 

Detailed Description

Definition at line 44 of file SignAlgoResolutions.h.

Member Typedef Documentation

Definition at line 62 of file SignAlgoResolutions.h.

typedef std::vector<double> metsig::SignAlgoResolutions::functionPars
private

Definition at line 63 of file SignAlgoResolutions.h.

Constructor & Destructor Documentation

metsig::SignAlgoResolutions::SignAlgoResolutions ( )
inline

Definition at line 47 of file SignAlgoResolutions.h.

47 :functionmap_(){;}
std::map< functionCombo, functionPars > functionmap_
metsig::SignAlgoResolutions::SignAlgoResolutions ( const edm::ParameterSet iConfig)

Definition at line 35 of file SignAlgoResolutions.cc.

References addResolutions().

35  :
36  functionmap_(),
37  ptResol_(0),
38  phiResol_(0)
39 {
40  addResolutions(iConfig);
41 }
std::map< functionCombo, functionPars > functionmap_
void addResolutions(const edm::ParameterSet &iConfig)

Member Function Documentation

void metsig::SignAlgoResolutions::addfunction ( const resolutionType  type,
const resolutionFunc  func,
std::vector< double >  parameters 
)
private

Definition at line 308 of file SignAlgoResolutions.cc.

References Parameters::parameters.

308  {
309 
310  functionCombo mypair(type,func);
311  functionmap_[mypair]=parameters;
312 
313 }
type
Definition: HCALResponse.h:22
dictionary parameters
Definition: Parameters.py:2
std::map< functionCombo, functionPars > functionmap_
std::pair< metsig::resolutionType, metsig::resolutionFunc > functionCombo
void metsig::SignAlgoResolutions::addResolutions ( const edm::ParameterSet iConfig)

Definition at line 144 of file SignAlgoResolutions.cc.

References metsig::caloEB, metsig::caloEE, metsig::caloHB, metsig::caloHE, metsig::caloHF, metsig::caloHO, metsig::ET, edm::ParameterSet::getParameter(), initializeJetResolutions(), metsig::PFtype1, metsig::PFtype2, metsig::PFtype3, metsig::PFtype4, metsig::PFtype5, metsig::PFtype6, metsig::PFtype7, and metsig::PHI.

Referenced by SignAlgoResolutions().

144  {
145  using namespace std;
146 
147  // Jet Resolutions - for now load from the files. Migrate to EventSetup asap.
149 
150  ptResolThreshold_ = iConfig.getParameter<double>("ptresolthreshold");
151 
152 
153  //get temporary low pT pfjet resolutions
154  for (int ieta=0; ieta<10; ieta++){
155  jdpt[ieta] = iConfig.getParameter<std::vector<double> >(Form("jdpt%d", ieta));
156  jdphi[ieta] = iConfig.getParameter<std::vector<double> >(Form("jdphi%d", ieta));
157  }
158 
159 
160  // for now: do this by hand - this can obviously also be done via ESSource etc.
161  functionPars etparameters(3,0);
162  functionPars phiparameters(1,0);
163  // set the parameters per function:
164  // ECAL, BARREL:
165  std::vector<double> ebet = iConfig.getParameter<std::vector<double> >("EB_EtResPar");
166  std::vector<double> ebphi = iConfig.getParameter<std::vector<double> >("EB_PhiResPar");
167 
168  etparameters[0]=ebet[0];
169  etparameters[1]=ebet[1];
170  etparameters[2]=ebet[2];
171  phiparameters[0]=ebphi[0];
172  addfunction(caloEB,ET,etparameters);
173  addfunction(caloEB,PHI,phiparameters);
174  // ECAL, ENDCAP:
175  std::vector<double> eeet = iConfig.getParameter<std::vector<double> >("EE_EtResPar");
176  std::vector<double> eephi = iConfig.getParameter<std::vector<double> >("EE_PhiResPar");
177 
178  etparameters[0]=eeet[0];
179  etparameters[1]=eeet[1];
180  etparameters[2]=eeet[2];
181  phiparameters[0]=eephi[0];
182  addfunction(caloEE,ET,etparameters);
183  addfunction(caloEE,PHI,phiparameters);
184  // HCAL, BARREL:
185  std::vector<double> hbet = iConfig.getParameter<std::vector<double> >("HB_EtResPar");
186  std::vector<double> hbphi = iConfig.getParameter<std::vector<double> >("HB_PhiResPar");
187 
188  etparameters[0]=hbet[0];
189  etparameters[1]=hbet[1];
190  etparameters[2]=hbet[2];
191  phiparameters[0]=hbphi[0];
192  addfunction(caloHB,ET,etparameters);
193  addfunction(caloHB,PHI,phiparameters);
194  // HCAL, ENDCAP:
195  std::vector<double> heet = iConfig.getParameter<std::vector<double> >("HE_EtResPar");
196  std::vector<double> hephi = iConfig.getParameter<std::vector<double> >("HE_PhiResPar");
197 
198  etparameters[0]=heet[0];
199  etparameters[1]=heet[1];
200  etparameters[2]=heet[2];
201  phiparameters[0]=hephi[0];
202  addfunction(caloHE,ET,etparameters);
203  addfunction(caloHE,PHI,phiparameters);
204  // HCAL, Outer
205  std::vector<double> hoet = iConfig.getParameter<std::vector<double> >("HO_EtResPar");
206  std::vector<double> hophi = iConfig.getParameter<std::vector<double> >("HO_PhiResPar");
207 
208 
209  etparameters[0]=hoet[0];
210  etparameters[1]=hoet[1];
211  etparameters[2]=hoet[2];
212  phiparameters[0]=hophi[0];
213  addfunction(caloHO,ET,etparameters);
214  addfunction(caloHO,PHI,phiparameters);
215  // HCAL, Forward
216  std::vector<double> hfet = iConfig.getParameter<std::vector<double> >("HF_EtResPar");
217  std::vector<double> hfphi = iConfig.getParameter<std::vector<double> >("HF_PhiResPar");
218 
219  etparameters[0]=hfet[0];
220  etparameters[1]=hfet[1];
221  etparameters[2]=hfet[2];
222  phiparameters[0]=hfphi[0];
223  addfunction(caloHF,ET,etparameters);
224  addfunction(caloHF,PHI,phiparameters);
225 
226 
227  // PF objects:
228  // type 1:
229  std::vector<double> pf1et = iConfig.getParameter<std::vector<double> >("PF_EtResType1");
230  std::vector<double> pf1phi = iConfig.getParameter<std::vector<double> >("PF_PhiResType1");
231  etparameters[0]=pf1et[0];
232  etparameters[1]=pf1et[1];
233  etparameters[2]=pf1et[2];
234  phiparameters[0]=pf1phi[0];
235  addfunction(PFtype1,ET,etparameters);
236  addfunction(PFtype1,PHI,phiparameters);
237 
238  // PF objects:
239  // type 2:
240  std::vector<double> pf2et = iConfig.getParameter<std::vector<double> >("PF_EtResType2");
241  std::vector<double> pf2phi = iConfig.getParameter<std::vector<double> >("PF_PhiResType2");
242  etparameters[0]=pf2et[0];
243  etparameters[1]=pf2et[1];
244  etparameters[2]=pf2et[2];
245  phiparameters[0]=pf2phi[0];
246  addfunction(PFtype2,ET,etparameters);
247  addfunction(PFtype2,PHI,phiparameters);
248 
249  // PF objects:
250  // type 3:
251  std::vector<double> pf3et = iConfig.getParameter<std::vector<double> >("PF_EtResType3");
252  std::vector<double> pf3phi = iConfig.getParameter<std::vector<double> >("PF_PhiResType3");
253  etparameters[0]=pf3et[0];
254  etparameters[1]=pf3et[1];
255  etparameters[2]=pf3et[2];
256  phiparameters[0]=pf3phi[0];
257  addfunction(PFtype3,ET,etparameters);
258  addfunction(PFtype3,PHI,phiparameters);
259 
260  // PF objects:
261  // type 4:
262  std::vector<double> pf4et = iConfig.getParameter<std::vector<double> >("PF_EtResType4");
263  std::vector<double> pf4phi = iConfig.getParameter<std::vector<double> >("PF_PhiResType4");
264  etparameters[0]=pf4et[0];
265  etparameters[1]=pf4et[1];
266  etparameters[2]=pf4et[2];
267  //phiparameters[0]=pf4phi[0];
268  addfunction(PFtype4,ET,etparameters);
269  addfunction(PFtype4,PHI,pf4phi); //use the same functional form for photon phi error as for pT, pass whole vector
270 
271  // PF objects:
272  // type 5:
273  std::vector<double> pf5et = iConfig.getParameter<std::vector<double> >("PF_EtResType5");
274  std::vector<double> pf5phi = iConfig.getParameter<std::vector<double> >("PF_PhiResType5");
275  etparameters[0]=pf5et[0];
276  etparameters[1]=pf5et[1];
277  etparameters[2]=pf5et[2];
278  phiparameters[0]=pf5phi[0];
279  addfunction(PFtype5,ET,etparameters);
280  addfunction(PFtype5,PHI,pf5phi);
281 
282  // PF objects:
283  // type 6:
284  std::vector<double> pf6et = iConfig.getParameter<std::vector<double> >("PF_EtResType6");
285  std::vector<double> pf6phi = iConfig.getParameter<std::vector<double> >("PF_PhiResType6");
286  etparameters[0]=pf6et[0];
287  etparameters[1]=pf6et[1];
288  etparameters[2]=pf6et[2];
289  phiparameters[0]=pf6phi[0];
290  addfunction(PFtype6,ET,etparameters);
291  addfunction(PFtype6,PHI,phiparameters);
292 
293 
294  // PF objects:
295  // type 7:
296  std::vector<double> pf7et = iConfig.getParameter<std::vector<double> >("PF_EtResType7");
297  std::vector<double> pf7phi = iConfig.getParameter<std::vector<double> >("PF_PhiResType7");
298  etparameters[0]=pf7et[0];
299  etparameters[1]=pf7et[1];
300  etparameters[2]=pf7et[2];
301  phiparameters[0]=pf7phi[0];
302  addfunction(PFtype7,ET,etparameters);
303  addfunction(PFtype7,PHI,phiparameters);
304 
305  return;
306 }
T getParameter(std::string const &) const
std::vector< double > jdpt[10]
std::vector< double > functionPars
void initializeJetResolutions(const edm::ParameterSet &iConfig)
void addfunction(const resolutionType type, const resolutionFunc func, std::vector< double > parameters)
std::vector< double > jdphi[10]
double metsig::SignAlgoResolutions::ElectronPtResolution ( const reco::PFCandidate c) const
private

Definition at line 416 of file SignAlgoResolutions.cc.

References relval_parameters_module::energy, reco::LeafCandidate::energy(), reco::LeafCandidate::eta(), and eta().

416  {
417 
418  double eta = c->eta();
419  double energy = c->energy();
420  double dEnergy = pfresol_->getEnergyResolutionEm(energy, eta);
421 
422  return dEnergy/cosh(eta);
423 }
T eta() const
virtual double eta() const
momentum pseudorapidity
virtual double energy() const
energy
double getEnergyResolutionEm(double CorrectedEnergy, double eta) const
PFEnergyResolution * pfresol_
double metsig::SignAlgoResolutions::EtFunction ( const functionPars x,
const functionPars par 
) const
private

Definition at line 341 of file SignAlgoResolutions.cc.

References query::result, and mathSSE::sqrt().

342 {
343  if(par.size()<3)
344  return 0.;
345  if(x.size()<1)
346  return 0.;
347  double et=x[0];
348  if(et<=0.)
349  return 0.;
350  double result = et*sqrt((par[2]*par[2])+(par[1]*par[1]/et)+(par[0]*par[0]/(et*et)));
351  return result;
352 }
T sqrt(T t)
Definition: SSEVec.h:46
tuple result
Definition: query.py:137
x
Definition: VDTMath.h:216
double metsig::SignAlgoResolutions::eval ( const resolutionType type,
const resolutionFunc func,
const double &  et,
const double &  phi,
const double &  eta,
const double &  p 
) const

Definition at line 53 of file SignAlgoResolutions.cc.

References eta(), AlCaHLTBitMon_ParallelJobs::p, phi, and vdt::x.

53  {
54 
55  functionPars x(4);
56  x[0]=et;
57  x[1]=phi;
58  x[2]=eta;
59  x[3]=p;
60  // std::cout << "getting function of type " << type << " " << func << " " << x[0] << " " << x[1] << " " << x[2] << " " << x[3] << std::endl;
61  return getfunc(type,func,x);
62 
63 }
type
Definition: HCALResponse.h:22
double getfunc(const resolutionType &type, const resolutionFunc &func, std::vector< double > &x) const
T eta() const
std::vector< double > functionPars
x
Definition: VDTMath.h:216
Definition: DDAxes.h:10
double metsig::SignAlgoResolutions::eval ( const resolutionType type,
const resolutionFunc func,
const double &  et,
const double &  phi,
const double &  eta 
) const

Definition at line 45 of file SignAlgoResolutions.cc.

References create_public_lumi_plots::exp, AlCaHLTBitMon_ParallelJobs::p, funct::sin(), and theta().

45  {
46  // derive p from et and eta;
47  double theta = 2*atan(exp(-eta));
48  double p = et / sin(theta); // rough assumption: take e and p equivalent...
49  return eval(type,func,et,phi,eta,p);
50 
51 }
type
Definition: HCALResponse.h:22
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
T eta() const
double eval(const resolutionType &type, const resolutionFunc &func, const double &et, const double &phi, const double &eta, const double &p) const
Definition: DDAxes.h:10
metsig::SigInputObj metsig::SignAlgoResolutions::evalPF ( const reco::PFCandidate candidate) const

Definition at line 65 of file SignAlgoResolutions.cc.

References reco::LeafCandidate::energy(), metsig::ET, reco::LeafCandidate::eta(), eta(), edm::Ref< C, T, F >::isNull(), mergeVDriftHistosByStation::name, reco::PFCandidate::particleId(), metsig::PFtype1, metsig::PFtype2, metsig::PFtype3, metsig::PFtype4, metsig::PFtype5, metsig::PFtype6, metsig::PFtype7, phi, metsig::PHI, reco::LeafCandidate::phi(), funct::sin(), reco::LeafCandidate::theta(), and reco::PFCandidate::trackRef().

65  {
66  double eta = candidate->eta();
67  double phi = candidate->phi();
68  double et = candidate->energy()*sin(candidate->theta());
69  resolutionType thetype;
70  std::string name;
71  int type = candidate->particleId();
72  switch (type)
73  {
74  case 1:
75  thetype=PFtype1;name="PFChargedHadron"; break;
76  case 2:
77  thetype=PFtype2;name="PFChargedEM"; break;
78  case 3:
79  thetype=PFtype3;name="PFMuon"; break;
80  case 4:
81  thetype=PFtype4;name="PFNeutralEM"; break;
82  case 5:
83  thetype=PFtype5;name="PFNeutralHadron"; break;
84  case 6:
85  thetype=PFtype6;name="PFtype6"; break;
86  case 7:
87  thetype=PFtype7;name="PFtype7"; break;
88  default:
89  thetype=PFtype7;name="PFunknown"; break;
90  }
91 
92  double d_et=0, d_phi=0; //d_phi here is the error on phi component of the et
93  reco::TrackRef trackRef = candidate->trackRef();
94  if(!trackRef.isNull()){
95  d_phi = et*trackRef->phiError();
96  d_et = (type==2) ? ElectronPtResolution(candidate) : trackRef->ptError();
97  //if(type==2) std::cout << eval(thetype,ET,et,phi,eta) << " " << trackRef->ptError() << " "<< ElectronPtResolution(candidate) << std::endl;
98  }
99  else{
100  d_et = eval(thetype,ET,et,phi,eta);
101  d_phi = eval(thetype,PHI,et,phi,eta);
102  }
103 
104  metsig::SigInputObj resultingobj(name,et,phi,d_et,d_phi);
105  return resultingobj;
106 }
type
Definition: HCALResponse.h:22
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T eta() const
double eval(const resolutionType &type, const resolutionFunc &func, const double &et, const double &phi, const double &eta, const double &p) const
virtual double eta() const
momentum pseudorapidity
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:339
double ElectronPtResolution(const reco::PFCandidate *c) const
virtual double energy() const
energy
bool isNull() const
Checks for null.
Definition: Ref.h:247
virtual double theta() const
momentum polar angle
virtual ParticleType particleId() const
Definition: PFCandidate.h:324
virtual double phi() const
momentum azimuthal angle
Definition: DDAxes.h:10
metsig::SigInputObj metsig::SignAlgoResolutions::evalPFJet ( const reco::PFJet jet) const

Definition at line 110 of file SignAlgoResolutions.cc.

References reco::LeafCandidate::eta(), reco::LeafCandidate::phi(), and reco::LeafCandidate::pt().

110  {
111 
112  double jpt = jet->pt();
113  double jphi = jet->phi();
114  double jeta = jet->eta();
115  double jdeltapt = 999.;
116  double jdeltapphi = 999.;
117 
118  if(jpt<ptResolThreshold_ && jpt<20.){ //use temporary fix for low pT jets
119  double feta = TMath::Abs(jeta);
120  int ieta = feta<5.? int(feta/0.5) : 9; //bin size = 0.5
121  int ipt = jpt>3. ? int(jpt-3./2) : 0; //bin size =2, starting from ptmin=3GeV
122  jdeltapt = jdpt[ieta][ipt];
123  jdeltapphi = jpt*jdphi[ieta][ipt];
124  }
125  else{
126  //use the resolution functions at |eta|=5 to avoid crash for jets with large eta.
127  if(jeta>5) jeta=5;
128  if(jeta<-5) jeta=-5;
129  TF1* fPtEta = ptResol_->parameterEta("sigma",jeta);
130  TF1* fPhiEta = phiResol_->parameterEta("sigma",jeta);
131  jdeltapt = jpt>ptResolThreshold_ ? jpt*fPtEta->Eval(jpt) : jpt*fPtEta->Eval(ptResolThreshold_);
132  jdeltapphi = jpt>ptResolThreshold_ ? jpt*fPhiEta->Eval(jpt) : jpt*fPhiEta->Eval(ptResolThreshold_);
133  delete fPtEta;
134  delete fPhiEta;
135  }
136 
137  std::string inputtype = "jet";
138  metsig::SigInputObj obj_jet(inputtype,jpt,jphi,jdeltapt,jdeltapphi);
139  //std::cout << "RESOLUTIONS JET: " << jpt << " " << jphi<< " " <<jdeltapt << " " << jdeltapphi << std::endl;
140  return obj_jet;
141 }
TF1 * parameterEta(const std::string &parameterName, float eta)
std::vector< double > jdpt[10]
virtual double eta() const
momentum pseudorapidity
virtual double pt() const
transverse momentum
std::vector< double > jdphi[10]
virtual double phi() const
momentum azimuthal angle
double metsig::SignAlgoResolutions::getfunc ( const resolutionType type,
const resolutionFunc func,
std::vector< double > &  x 
) const
private

Definition at line 315 of file SignAlgoResolutions.cc.

References metsig::CONSTPHI, metsig::ET, metsig::PHI, query::result, edm::second(), metsig::TRACKP, and makeHLTPrescaleTable::values.

315  {
316 
317  double result=0;
318  functionCombo mypair(type,func);
319 
320  if(functionmap_.count(mypair)==0){
321  return result;
322  }
323 
324  functionPars values = (functionmap_.find(mypair))->second;
325  switch ( func ){
326  case metsig::ET :
327  return EtFunction(x,values);
328  case metsig::PHI :
329  return PhiFunction(x,values);
330  case metsig::TRACKP :
331  return PFunction(x,values);
332  case metsig::CONSTPHI :
333  return PhiConstFunction(x,values);
334  }
335 
336  // std::cout << "returning function " << type << " " << func << " " << result << " " << x[0] << std::endl;
337 
338  return result;
339 }
type
Definition: HCALResponse.h:22
double PhiFunction(const functionPars &x, const functionPars &par) const
U second(std::pair< T, U > const &p)
tuple result
Definition: query.py:137
double EtFunction(const functionPars &x, const functionPars &par) const
std::vector< double > functionPars
double PhiConstFunction(const functionPars &x, const functionPars &par) const
std::map< functionCombo, functionPars > functionmap_
x
Definition: VDTMath.h:216
std::pair< metsig::resolutionType, metsig::resolutionFunc > functionCombo
double PFunction(const functionPars &x, const functionPars &par) const
void metsig::SignAlgoResolutions::initializeJetResolutions ( const edm::ParameterSet iConfig)
private

Definition at line 386 of file SignAlgoResolutions.cc.

References dtNoiseDBValidation_cfg::cerr, launcher::cmssw_base, cmsBenchmark::cmssw_release_base, edm::ParameterSet::getParameter(), and scaleCards::path.

Referenced by addResolutions().

386  {
387 
388  using namespace std;
389 
390  // only reinitialize the resolutsion if the pointers are zero
391  if ( ptResol_ == 0 ) {
392  string resolutionsAlgo = iConfig.getParameter<std::string>("resolutionsAlgo");
393  string resolutionsEra = iConfig.getParameter<std::string>("resolutionsEra");
394 
395  string cmssw_base(getenv("CMSSW_BASE"));
396  string cmssw_release_base(getenv("CMSSW_RELEASE_BASE"));
397  string path = cmssw_base + "/src/CondFormats/JetMETObjects/data";
398  struct stat st;
399  if (stat(path.c_str(),&st)!=0) {
400  path = cmssw_release_base + "/src/CondFormats/JetMETObjects/data";
401  }
402  if (stat(path.c_str(),&st)!=0) {
403  cerr<<"ERROR: tried to set path but failed, abort."<<endl;
404  }
405  string era(resolutionsEra);
406  string alg(resolutionsAlgo);
407  string ptFileName = path + "/" + era + "_PtResolution_" +alg+".txt";
408  string phiFileName = path + "/" + era + "_PhiResolution_"+alg+".txt";
409 
410  ptResol_ = new JetResolution(ptFileName,false);
411  phiResol_ = new JetResolution(phiFileName,false);
412  }
413 }
T getParameter(std::string const &) const
list path
Definition: scaleCards.py:51
list cmssw_release_base
Definition: cmsBenchmark.py:61
list cmssw_base
Definition: launcher.py:19
bool metsig::SignAlgoResolutions::isFilled ( ) const
inline

Definition at line 55 of file SignAlgoResolutions.h.

References functionmap_.

55 {return functionmap_.size()>0;}
std::map< functionCombo, functionPars > functionmap_
double metsig::SignAlgoResolutions::PFunction ( const functionPars x,
const functionPars par 
) const
private

Definition at line 374 of file SignAlgoResolutions.cc.

375 {
376  // not currently implemented
377  return 0;
378 }
double metsig::SignAlgoResolutions::PhiConstFunction ( const functionPars x,
const functionPars par 
) const
private

Definition at line 380 of file SignAlgoResolutions.cc.

381 {
382  return par[0];
383 }
double metsig::SignAlgoResolutions::PhiFunction ( const functionPars x,
const functionPars par 
) const
private

Definition at line 355 of file SignAlgoResolutions.cc.

References mathSSE::sqrt().

356 {
357  double et=x[0];
358  if(et<=0.){
359  return 0.;
360  }
361 
362  //if 1 parameter is C provided, returns C*pT, if three parameters N, S, C are provided, it returns the usual resolution value, as for sigmaPt
363  if(par.size()!=1 && par.size()!=3){//only 1 or 3 parameters supported for phi function
364  return 0.;
365  }
366  else if(par.size()==1){
367  return par[0]*et;
368  }
369  else{
370  return et*sqrt((par[2]*par[2])+(par[1]*par[1]/et)+(par[0]*par[0]/(et*et)));
371  }
372 
373 }
T sqrt(T t)
Definition: SSEVec.h:46
x
Definition: VDTMath.h:216

Member Data Documentation

std::map<functionCombo,functionPars> metsig::SignAlgoResolutions::functionmap_
private

Definition at line 64 of file SignAlgoResolutions.h.

Referenced by isFilled().

std::vector<double> metsig::SignAlgoResolutions::jdphi[10]
private

Definition at line 77 of file SignAlgoResolutions.h.

std::vector<double> metsig::SignAlgoResolutions::jdpt[10]
private

Definition at line 76 of file SignAlgoResolutions.h.

PFEnergyResolution* metsig::SignAlgoResolutions::pfresol_
private

Definition at line 81 of file SignAlgoResolutions.h.

JetResolution* metsig::SignAlgoResolutions::phiResol_
private

Definition at line 80 of file SignAlgoResolutions.h.

JetResolution* metsig::SignAlgoResolutions::ptResol_
private

Definition at line 79 of file SignAlgoResolutions.h.

double metsig::SignAlgoResolutions::ptResolThreshold_
private

Definition at line 72 of file SignAlgoResolutions.h.