CMS 3D CMS Logo

PFMEtSignInterfaceBase.cc
Go to the documentation of this file.
2 
5 
7 
8 #include <TVectorD.h>
9 
10 const double defaultPFMEtResolutionX = 10.;
11 const double defaultPFMEtResolutionY = 10.;
12 
13 const double epsilon = 1.e-9;
14 
16  : pfMEtResolution_(nullptr), inputFile_(nullptr), lut_(nullptr) {
18 
19  if (cfg.exists("addJERcorr")) {
20  edm::ParameterSet cfgJERcorr = cfg.getParameter<edm::ParameterSet>("addJERcorr");
21  edm::FileInPath inputFileName = cfgJERcorr.getParameter<edm::FileInPath>("inputFileName");
22  std::string lutName = cfgJERcorr.getParameter<std::string>("lutName");
23  if (inputFileName.location() != edm::FileInPath::Local)
24  throw cms::Exception("PFMEtSignInterfaceBase") << " Failed to find File = " << inputFileName << " !!\n";
25 
26  inputFile_ = new TFile(inputFileName.fullPath().data());
27  lut_ = dynamic_cast<TH2*>(inputFile_->Get(lutName.data()));
28  if (!lut_)
29  throw cms::Exception("PFMEtSignInterfaceBase") << " Failed to load LUT = " << lutName.data()
30  << " from file = " << inputFileName.fullPath().data() << " !!\n";
31  }
32 
33  verbosity_ = cfg.exists("verbosity") ? cfg.getParameter<int>("verbosity") : 0;
34 }
35 
37  delete pfMEtResolution_;
38  delete inputFile_;
39  delete lut_;
40 }
41 
42 reco::METCovMatrix PFMEtSignInterfaceBase::operator()(const std::vector<metsig::SigInputObj>& pfMEtSignObjects) const {
43  // if ( this->verbosity_ ) {
44  // std::cout << "<PFMEtSignInterfaceBase::operator()>:" << std::endl;
45  // std::cout << " pfMEtSignObjects: entries = " << pfMEtSignObjects.size() << std::endl;
46  // double dpt2Sum = 0.;
47  // for ( std::vector<metsig::SigInputObj>::const_iterator pfMEtSignObject = pfMEtSignObjects.begin();
48  // pfMEtSignObject != pfMEtSignObjects.end(); ++pfMEtSignObject ) {
49  // std::cout << pfMEtSignObject->get_type() << ": pt = " << pfMEtSignObject->get_energy() << ","
50  // << " phi = " << pfMEtSignObject->get_phi() << " --> dpt = " << pfMEtSignObject->get_sigma_e() << std::endl;
51  // dpt2Sum += pfMEtSignObject->get_sigma_e();
52  // }
53  // std::cout << "--> sqrt(sum(dpt^2)) = " << sqrt(dpt2Sum) << std::endl;
54  // }
55 
56  reco::METCovMatrix pfMEtCov;
57  if (pfMEtSignObjects.size() >= 2) {
58  metsig::significanceAlgo pfMEtSignAlgorithm;
59  pfMEtSignAlgorithm.addObjects(pfMEtSignObjects);
60  pfMEtCov = pfMEtSignAlgorithm.getSignifMatrix();
61 
62  double det = 0;
63  pfMEtCov.Det(det);
64 
65  // if ( this->verbosity_ && std::abs(det) > epsilon ) {
66  // //keep TMatrixD as it is much easier to find
67  // //eigenvectors and values than with SMatrix;
68  // //not used anyway, except for debugging
69  // TMatrixD tmpMatrix(2,2);
70  // tmpMatrix(0,0) = pfMEtCov(0,0);
71  // tmpMatrix(0,1) = pfMEtCov(0,1);
72  // tmpMatrix(1,0) = pfMEtCov(1,0);
73  // tmpMatrix(1,1) = pfMEtCov(1,1);
74 
75  // TVectorD eigenValues(2);
76  // TMatrixD eigenVectors = tmpMatrix.EigenVectors(eigenValues);
77  // // CV: eigenvectors are stored in columns
78  // // and are sorted such that the one corresponding to the highest eigenvalue is in the **first** column
79  // for ( unsigned iEigenVector = 0; iEigenVector < 2; ++iEigenVector ) {
80  // std::cout << "eigenVector #" << iEigenVector << " (eigenValue = " << eigenValues(iEigenVector) << "):"
81  // << " x = " << eigenVectors(0, iEigenVector) << ", y = " << eigenVectors(1, iEigenVector) << std::endl;
82  // }
83  // }
84 
85  //--- substitute (PF)MEt resolution matrix by default values
86  // in case resolution matrix cannot be inverted
87  if (std::abs(det) < epsilon) {
88  edm::LogWarning("PFMEtSignInterfaceBase::operator()")
89  << "Inversion of PFMEt covariance matrix failed, det = " << det
90  << " --> replacing covariance matrix by resolution defaults !!";
92  pfMEtCov(0, 1) = 0.;
93  pfMEtCov(1, 0) = 0.;
95  }
96  } else {
97  pfMEtCov(0, 0) = 0.;
98  pfMEtCov(0, 1) = 0.;
99  pfMEtCov(1, 0) = 0.;
100  pfMEtCov(1, 1) = 0.;
101  }
102 
103  return pfMEtCov;
104 }
significanceAlgo.h
MessageLogger.h
PFMEtSignInterfaceBase::~PFMEtSignInterfaceBase
~PFMEtSignInterfaceBase()
Definition: PFMEtSignInterfaceBase.cc:36
PFMEtSignInterfaceBase.h
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
PFMEtSignInterfaceBase::pfMEtResolution_
metsig::SignAlgoResolutions * pfMEtResolution_
Definition: PFMEtSignInterfaceBase.h:182
edm::FileInPath::Local
Definition: FileInPath.h:63
PFMEtSignInterfaceBase::lut_
TH2 * lut_
Definition: PFMEtSignInterfaceBase.h:188
edm::FileInPath
Definition: FileInPath.h:61
metsig::significanceAlgo::addObjects
const void addObjects(const std::vector< metsig::SigInputObj > &EventVec)
Definition: significanceAlgo.cc:109
InefficientDoubleROC.inputFileName
inputFileName
Definition: InefficientDoubleROC.py:437
defaultPFMEtResolutionX
const double defaultPFMEtResolutionX
Definition: PFMEtSignInterfaceBase.cc:10
PFMEtSignInterfaceBase::operator()
reco::METCovMatrix operator()(const std::vector< metsig::SigInputObj > &) const
Definition: PFMEtSignInterfaceBase.cc:42
PFMEtSignInterfaceBase::inputFile_
TFile * inputFile_
Definition: PFMEtSignInterfaceBase.h:187
edm::ParameterSet
Definition: ParameterSet.h:47
epsilon
const double epsilon
Definition: PFMEtSignInterfaceBase.cc:13
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
looper.cfg
cfg
Definition: looper.py:296
PFMEtSignInterfaceBase::verbosity_
int verbosity_
Definition: PFMEtSignInterfaceBase.h:190
metsig::significanceAlgo
Definition: significanceAlgo.h:86
metsig::SignAlgoResolutions
Definition: SignAlgoResolutions.h:61
PFMEtSignInterfaceBase::PFMEtSignInterfaceBase
PFMEtSignInterfaceBase(const edm::ParameterSet &)
Definition: PFMEtSignInterfaceBase.cc:15
Exception
Definition: hltDiff.cc:245
metsig::significanceAlgo::getSignifMatrix
reco::METCovMatrix getSignifMatrix() const
Definition: significanceAlgo.h:99
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
Exception.h
cms::Exception
Definition: Exception.h:70
defaultPFMEtResolutionY
const double defaultPFMEtResolutionY
Definition: PFMEtSignInterfaceBase.cc:11
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
reco::METCovMatrix
ROOT::Math::SMatrix< double, 2 > METCovMatrix
Definition: MET.h:39