CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
CovarianceParameterization Class Reference

#include <CovarianceParameterization.h>

Classes

struct  CompressionSchema
 

Public Member Functions

 CovarianceParameterization ()
 
bool isValid () const
 
void load (int version)
 
int loadedVersion () const
 
float meanValue (int i, int j, int sign, float pt, float eta, int nHits, int pixelHits, float cii=1., float cjj=1.) const
 
float pack (float value, int schema, int i, int j, float pt, float eta, int nHits, int pixelHits, float cii=1., float cjj=1.) const
 
float unpack (uint16_t packed, int schema, int i, int j, float pt, float eta, int nHits, int pixelHits, float cii=1., float cjj=1.) const
 

Static Public Member Functions

static int index (int i, int j)
 

Private Member Functions

void addTheHistogram (std::vector< TH3D *> *HistoVector, std::string StringToAddInTheName, int i, int j, TFile &fileToRead)
 
void readFile (TFile &)
 

Private Attributes

std::vector< TH3D * > cov_elements_noPixelHit
 
std::vector< TH3D * > cov_elements_pixelHit
 
TFile * fileToRead_
 
int loadedVersion_
 
std::unordered_map< uint16_t, CompressionSchemaschemas
 

Detailed Description

Definition at line 24 of file CovarianceParameterization.h.

Constructor & Destructor Documentation

◆ CovarianceParameterization()

CovarianceParameterization::CovarianceParameterization ( )
inline

Member Function Documentation

◆ addTheHistogram()

void CovarianceParameterization::addTheHistogram ( std::vector< TH3D *> *  HistoVector,
std::string  StringToAddInTheName,
int  i,
int  j,
TFile &  fileToRead 
)
private

Definition at line 141 of file CovarianceParameterization.cc.

References mps_fire::i, dqmiolumiharvest::j, and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by readFile().

142  {
143  std::string List_covName[5] = {"qoverp", "lambda", "phi", "dxy", "dsz"};
144 
145  std::string histoNameString = "covariance_" + List_covName[i] + "_" + List_covName[j] + StringToAddInTheName +
146  "parametrization"; // + "_entries";
147  TH3D *matrixElememtHistogramm = (TH3D *)fileToRead.Get(histoNameString.c_str());
148  HistoVector->push_back(matrixElememtHistogramm);
149 }

◆ index()

static int CovarianceParameterization::index ( int  i,
int  j 
)
inlinestatic

Definition at line 26 of file CovarianceParameterization.h.

References mps_fire::i, and dqmiolumiharvest::j.

Referenced by CovarianceParameterization::CompressionSchema::operator()().

26  {
27  if (i >= j)
28  return j + i * (i + 1) / 2;
29  else
30  return i + j * (j + 1) / 2;
31  }

◆ isValid()

bool CovarianceParameterization::isValid ( void  ) const
inline

◆ load()

void CovarianceParameterization::load ( int  version)

Definition at line 81 of file CovarianceParameterization.cc.

References haddnano::cl, fileToRead_, printsummarytable::folder, edm::FileInPath::fullPath(), mps_fire::i, createfilelist::int, dqmiolumiharvest::j, dqmdumpme::k, crabWrapper::key, loadedVersion_, GetRecoTauVFromDQM_MC_cff::next, AlCaHLTBitMon_ParallelJobs::p, readFile(), schemas, AlCaHLTBitMon_QueryRunRegistry::string, and BeamSplash_cfg::version.

Referenced by pat::PackedCandidate::covarianceParameterization().

81  {
82  edm::FileInPath fip(
83  fmt::sprintf("DataFormats/PatCandidates/data/CovarianceParameterization_version%d.root", version));
84  fileToRead_ = TFile::Open(fip.fullPath().c_str());
85  TFile &fileToRead = *fileToRead_;
86  //Read files from here fip.fullPath().c_str();
87  if (fileToRead.IsOpen()) {
88  readFile(fileToRead);
89 
90  TIter next(((TDirectoryFile *)fileToRead.Get("schemas"))->GetListOfKeys());
91  TKey *key;
92  while ((key = (TKey *)next())) {
93  TClass *cl = gROOT->GetClass(key->GetClassName());
94  if (!cl->InheritsFrom("TDirectoryFile"))
95  continue;
96  std::string schemaNumber = key->ReadObj()->GetName();
97  uint16_t schemaN = std::stoi(schemaNumber);
98  //for (int folderNumber = 0; folderNumber < 6 ; folderNumber++) {
99  CompressionSchema schema;
100  for (int i = 0; i < 5; i++) {
101  for (int j = i; j < 5; j++) { //FILLING ONLY THE SCHEMA OF SOME ELEMENTS
102  std::string folder = "schemas/" + schemaNumber + "/" + char(48 + i) + char(48 + j);
103  std::string methodString = folder + "/method";
104  std::string targetString = folder + "/target";
105  std::string bitString = folder + "/bit";
106  std::vector<float> vParams;
107  TVector *p = (TVector *)fileToRead.Get((folder + "/param").c_str());
108  vParams.reserve(p->GetNoElements());
109  for (int k = 0; k < p->GetNoElements(); k++) {
110  vParams.push_back((*p)[k]);
111  }
112 
113  schema(i, j) = CompressionElement(
114  (CompressionElement::Method)((TParameter<int> *)fileToRead.Get(methodString.c_str()))->GetVal(),
115  (CompressionElement::Target)((TParameter<int> *)fileToRead.Get(targetString.c_str()))->GetVal(),
116  (int)((TParameter<int> *)fileToRead.Get(bitString.c_str()))->GetVal(),
117  vParams);
118  }
119  }
120  schemas[schemaN] = schema;
121  }
122 
124  } else {
125  loadedVersion_ = -1;
126  }
127 }
std::unordered_map< uint16_t, CompressionSchema > schemas

◆ loadedVersion()

int CovarianceParameterization::loadedVersion ( ) const
inline

◆ meanValue()

float CovarianceParameterization::meanValue ( int  i,
int  j,
int  sign,
float  pt,
float  eta,
int  nHits,
int  pixelHits,
float  cii = 1.,
float  cjj = 1. 
) const

Definition at line 151 of file CovarianceParameterization.cc.

References funct::abs(), cov_elements_noPixelHit, cov_elements_pixelHit, PVValHelper::eta, muonRecoAnalyzer_cfi::etaBin, mps_fire::i, dqmiolumiharvest::j, nHits, InitialStepPreSplitting_cff::pixelHits, DiDispStaMuonMonitor_cfi::pt, muonRecoAnalyzer_cfi::ptBin, and Validation_hcalonly_cfi::sign.

Referenced by pack(), and unpack().

152  {
153  int hitNumberToUse = nHits;
154  if (hitNumberToUse < 2)
155  hitNumberToUse = 2;
156  if (hitNumberToUse > 32)
157  hitNumberToUse = 32;
158  int ptBin = cov_elements_pixelHit[0]->GetXaxis()->FindBin(pt);
159  int etaBin = cov_elements_pixelHit[0]->GetYaxis()->FindBin(std::abs(eta));
160  int hitBin = cov_elements_pixelHit[0]->GetZaxis()->FindBin(hitNumberToUse);
161  int min_idx = i;
162  int max_idx = j;
163 
164  if (i > j) {
165  min_idx = j;
166  max_idx = i;
167  }
168 
169  int indexOfTheHitogramInTheList = ((9 - min_idx) * min_idx) / 2 + max_idx;
170 
171  double meanValue = 0.;
172  if (pixelHits > 0) {
173  meanValue = sign * cov_elements_pixelHit[indexOfTheHitogramInTheList]->GetBinContent(ptBin, etaBin, hitBin);
174  } else {
175  meanValue = sign * cov_elements_noPixelHit[indexOfTheHitogramInTheList]->GetBinContent(ptBin, etaBin, hitBin);
176  }
177  return meanValue;
178 }
float meanValue(int i, int j, int sign, float pt, float eta, int nHits, int pixelHits, float cii=1., float cjj=1.) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< TH3D * > cov_elements_pixelHit
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
std::vector< TH3D * > cov_elements_noPixelHit

◆ pack()

float CovarianceParameterization::pack ( float  value,
int  schema,
int  i,
int  j,
float  pt,
float  eta,
int  nHits,
int  pixelHits,
float  cii = 1.,
float  cjj = 1. 
) const

Definition at line 180 of file CovarianceParameterization.cc.

References PVValHelper::eta, mps_fire::i, dqmiolumiharvest::j, meanValue(), nHits, InitialStepPreSplitting_cff::pixelHits, DiDispStaMuonMonitor_cfi::pt, schemas, edm::second(), and std::swap().

Referenced by pat::PackedCandidate::packCovarianceElement().

181  {
182  if (i > j)
183  std::swap(i, j);
184  float ref = meanValue(i, j, 1., pt, eta, nHits, pixelHits, cii, cjj);
185  if (ref == 0) {
186  schema = 0;
187  }
188  if (schema == 0 && i == j && (i == 2 || i == 0))
189  ref = 1. / (pt * pt);
190  /* //Used for debugging, to be later removed
191  uint16_t p=(*schemas.find(schema)).second(i,j).pack(value,ref);
192  float up=(*schemas.find(schema)).second(i,j).unpack(p,ref);
193  std::cout << "check " << i << " " << j << " " << value << " " << up << " " << p << " " << ref << " " << schema<< std::endl;*/
194  return (*schemas.find(schema)).second(i, j).pack(value, ref);
195 }
U second(std::pair< T, U > const &p)
float meanValue(int i, int j, int sign, float pt, float eta, int nHits, int pixelHits, float cii=1., float cjj=1.) const
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
Definition: value.py:1
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
std::unordered_map< uint16_t, CompressionSchema > schemas

◆ readFile()

void CovarianceParameterization::readFile ( TFile &  f)
private

Definition at line 129 of file CovarianceParameterization.cc.

References addTheHistogram(), cov_elements_noPixelHit, cov_elements_pixelHit, f, mps_fire::i, dqmiolumiharvest::j, and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by load().

129  {
130  for (int i = 0; i < 5; i++) {
131  for (int j = i; j < 5; j++) {
132  std::string String_first_positive = "_pixel_";
133  std::string String_second_positive = "_noPixel_";
134 
135  addTheHistogram(&cov_elements_pixelHit, String_first_positive, i, j, f);
136  addTheHistogram(&cov_elements_noPixelHit, String_second_positive, i, j, f);
137  }
138  }
139 }
double f[11][100]
std::vector< TH3D * > cov_elements_pixelHit
void addTheHistogram(std::vector< TH3D *> *HistoVector, std::string StringToAddInTheName, int i, int j, TFile &fileToRead)
std::vector< TH3D * > cov_elements_noPixelHit

◆ unpack()

float CovarianceParameterization::unpack ( uint16_t  packed,
int  schema,
int  i,
int  j,
float  pt,
float  eta,
int  nHits,
int  pixelHits,
float  cii = 1.,
float  cjj = 1. 
) const

Definition at line 196 of file CovarianceParameterization.cc.

References MillePedeFileConverter_cfg::e, PVValHelper::eta, mps_fire::i, dqmiolumiharvest::j, meanValue(), nHits, InitialStepPreSplitting_cff::pixelHits, DiDispStaMuonMonitor_cfi::pt, schemas, edm::second(), and std::swap().

Referenced by pat::PackedCandidate::unpackCovarianceElement().

198  {
199  if (i > j)
200  std::swap(i, j);
201  float ref = meanValue(i, j, 1., pt, eta, nHits, pixelHits, cii, cjj);
202  if (ref == 0) {
203  schema = 0;
204  }
205  if (schema == 0 && i == j && (i == 2 || i == 0))
206  ref = 1. / (pt * pt);
207  if (i == j && (*schemas.find(schema)).second(i, j).unpack(packed, ref) == 0)
208  return 1e-9;
209  else
210  return (*schemas.find(schema)).second(i, j).unpack(packed, ref);
211 }
U second(std::pair< T, U > const &p)
float meanValue(int i, int j, int sign, float pt, float eta, int nHits, int pixelHits, float cii=1., float cjj=1.) const
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
std::unordered_map< uint16_t, CompressionSchema > schemas

Member Data Documentation

◆ cov_elements_noPixelHit

std::vector<TH3D *> CovarianceParameterization::cov_elements_noPixelHit
private

Definition at line 73 of file CovarianceParameterization.h.

Referenced by meanValue(), and readFile().

◆ cov_elements_pixelHit

std::vector<TH3D *> CovarianceParameterization::cov_elements_pixelHit
private

Definition at line 72 of file CovarianceParameterization.h.

Referenced by meanValue(), and readFile().

◆ fileToRead_

TFile* CovarianceParameterization::fileToRead_
private

Definition at line 70 of file CovarianceParameterization.h.

Referenced by load().

◆ loadedVersion_

int CovarianceParameterization::loadedVersion_
private

Definition at line 69 of file CovarianceParameterization.h.

Referenced by isValid(), load(), and loadedVersion().

◆ schemas

std::unordered_map<uint16_t, CompressionSchema> CovarianceParameterization::schemas
private

Definition at line 71 of file CovarianceParameterization.h.

Referenced by load(), pack(), and unpack().