CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Types | Private Member Functions | Private Attributes
MuonErrorMatrix Class Reference

#include <MuonErrorMatrix.h>

Public Types

enum  action { use, constructor }
 enum type to define if the class is used as a tool or to be created More...
 

Public Member Functions

void adjust (FreeTrajectoryState &state)
 adjust the error matrix on the state More...
 
void adjust (TrajectoryStateOnSurface &state)
 adjust the error matrix on the state More...
 
void close ()
 close the root file attached to the class More...
 
void complicatedTerm (const AlgebraicSymMatrix55 &input, AlgebraicSymMatrix55 &output)
 convert sigma/rho -> sigma2/COV More...
 
int findBin (TAxis *axis, double value)
 method to get the bin index, taking care of under/overlow: first(1)/last(GetNbins())returned More...
 
CurvilinearTrajectoryError get (GlobalVector momentum, bool convolute=true)
 main method to be used. Retrieve a 5x5 symetrical matrix according to parametrization of error or scale factor More...
 
TProfile3D * get (int i, int j)
 actually get access to the TProfile3D used for the parametrization More...
 
CurvilinearTrajectoryError getFast (GlobalVector momentum)
 
unsigned int index (int i, int j)
 
 MuonErrorMatrix (const edm::ParameterSet &pset)
 constructor from a parameter set More...
 
void simpleTerm (const AlgebraicSymMatrix55 &input, AlgebraicSymMatrix55 &output)
 convert sigma2/COV -> sigma/rho More...
 
 ~MuonErrorMatrix ()
 destructor More...
 

Static Public Member Functions

static bool divide (CurvilinearTrajectoryError &num_error, const CurvilinearTrajectoryError &denom_error)
 divide term by term the two matrix More...
 
static void multiply (CurvilinearTrajectoryError &initial_error, const CurvilinearTrajectoryError &scale_error)
 multiply term by term the two matrix More...
 
static double Term (const AlgebraicSymMatrix55 &curv, int i, int j)
 provide the numerical value used. sigma or correlation factor More...
 

Static Public Attributes

static const TString vars [5] ={"#frac{q}{|p|}","#lambda","#varphi_{0}","X_{T}","Y_{T}"}
 names of the variables of the 5x5 error matrix More...
 

Private Types

enum  TermAction { error, scale, assign }
 decide whether to scale of to assigne terms More...
 

Private Member Functions

TProfile3D * Index (int i, int j)
 internal method to get access to the profiles More...
 
int Pindex (int i, int j)
 internal methods to get the index of a matrix term. More...
 
double Rms (GlobalVector &momentum, int i, int j)
 internal method that retreives the error on the value of the parametrization for term i,j More...
 
double Value (GlobalVector &momentum, int i, int j, bool convolute=true)
 internal method that retreives the value of the parametrization for term i,j More...
 

Private Attributes

std::string theCategory
 log category: "MuonErrorMatrix" More...
 
TDirectory * theD
 the attached root file, where the parametrization is saved More...
 
TProfile3D * theData [15]
 15 TProfile, each holding he parametrization of each term of the 5x5 More...
 
TProfile3D * theData_fast [5][5]
 
TermAction theTermAction [15]
 

Detailed Description

This class holds an error matrix paramertization. Either the error matrix value or scale factors. The implementation uses TProfile3D with pt, eta, phi axis. Error/Scale factor matrix is obtained using get(GlobalVector momentum)

Author
Jean-Roch Vlimant UCSB
Finn Rebassoo UCSB

Definition at line 29 of file MuonErrorMatrix.h.

Member Enumeration Documentation

enum type to define if the class is used as a tool or to be created

Enumerator
use 
constructor 

Definition at line 32 of file MuonErrorMatrix.h.

decide whether to scale of to assigne terms

Enumerator
error 
scale 
assign 

Definition at line 88 of file MuonErrorMatrix.h.

Constructor & Destructor Documentation

MuonErrorMatrix::MuonErrorMatrix ( const edm::ParameterSet pset)

constructor from a parameter set

Definition at line 17 of file MuonErrorMatrix.cc.

References a, assign, constructor, data, TrackerOfflineValidation_Dqm_cff::dirName, error, MuonErrorMatrixValues_cff::errorMatrixValuesPSet, edm::ParameterSet::exists(), cmsRelvalreport::exit, MillePedeFileConverter_cfg::fileName, edm::FileInPath::fullPath(), edm::ParameterSet::getParameter(), mps_fire::i, LogDebug, LogTrace, maxEta, MuonErrorMatrixAnalyzer_cfi::maxPhi, MuonErrorMatrixAnalyzer_cfi::maxPt, cutBasedElectronID_CSA14_50ns_V0_cff::minEta, MuonErrorMatrixAnalyzer_cfi::minPhi, lostTracks_cfi::minPt, MuonErrorMatrixAnalyzer_cfi::NEta, MuonErrorMatrixAnalyzer_cfi::NPhi, MuonErrorMatrixAnalyzer_cfi::NPt, packedPFCandidateRefMixer_cfi::pf, Pi, Pindex(), scale, harvestTrackValidationPlots::str, AlCaHLTBitMon_QueryRunRegistry::string, theCategory, theD, theData, theData_fast, theTermAction, use, MuonErrorMatrixValues_cff::values, and vars.

17  :theD(0){
18  theCategory="MuonErrorMatrix";
19  std::string action = iConfig.getParameter<std::string>("action");
20 
21  bool madeFromCff=iConfig.exists("errorMatrixValuesPSet");
23 
25  if (!madeFromCff){ fileName = iConfig.getParameter<std::string>("rootFileName");}
26  else { errorMatrixValuesPSet =iConfig.getParameter<edm::ParameterSet>("errorMatrixValuesPSet");}
28 
29  int NPt=5;
30  std::vector<double> xBins;
31  double * xBinsArray = 0;
32  double minPt=1;
33  double maxPt=200;
34  int NEta=5;
35  std::vector<double> yBins;
36  double * yBinsArray =0;
37  double minEta=0;
38  double maxEta=2.5;
39  int NPhi=1;
40  double minPhi=-TMath::Pi();
41  double maxPhi=TMath::Pi();
42 
43  if (action!="use"){
44  a = constructor;
45 
46  NPt = iConfig.getParameter<int>("NPt");
47  if (NPt!=0){
48  minPt = iConfig.getParameter<double>("minPt");
49  maxPt = iConfig.getParameter<double>("maxPt");}
50  else{
51  xBins = iConfig.getParameter<std::vector<double> >("PtBins");
52  if (xBins.size()==0){edm::LogError( theCategory)<<"Npt=0 and no entries in the vector. I will do aseg fault soon.";}
53  NPt = xBins.size()-1;
54  xBinsArray = &(xBins.front());
55  minPt = xBins.front();
56  maxPt = xBins.back();}
57 
58  NEta = iConfig.getParameter<int>("NEta");
59  if (NEta!=0){
60  minEta = iConfig.getParameter<double>("minEta");
61  maxEta = iConfig.getParameter<double>("maxEta");}
62  else{
63  yBins = iConfig.getParameter<std::vector<double> >("EtaBins");
64  if (yBins.size()==0){edm::LogError( theCategory)<<"NEta=0 and no entries in the vector. I will do aseg fault soon.";}
65  NEta = yBins.size()-1;
66  yBinsArray = &(yBins.front());
67  minPt = yBins.front();
68  maxPt = yBins.back();}
69 
70  NPhi = iConfig.getParameter<int>("NPhi");
71  std::stringstream get(iConfig.getParameter<std::string>("minPhi"));
72  if (get.str() =="-Pi")
73  { minPhi =-TMath::Pi();}
74  else if(get.str() =="Pi")
75  { minPhi =TMath::Pi();}
76  else { get>>minPhi;}
77  get.str(iConfig.getParameter<std::string>("maxPhi"));
78  if (get.str() =="-Pi")
79  { maxPhi =-TMath::Pi();}
80  else if(get.str() =="Pi")
81  { maxPhi =TMath::Pi();}
82  else { get>>maxPhi;}
83  }//action!=use
84  else if (madeFromCff){
85 
86  xBins = errorMatrixValuesPSet.getParameter<std::vector<double> >("xAxis");
87  NPt = xBins.size()-1;
88  xBinsArray = &(xBins.front());
89  minPt = xBins.front();
90  maxPt = xBins.back();
91 
92  yBins = errorMatrixValuesPSet.getParameter<std::vector<double> >("yAxis");
93  NEta = yBins.size()-1;
94  yBinsArray = &(yBins.front());
95  minPt = yBins.front();
96  maxPt = yBins.back();
97 
98  std::vector<double> zBins = errorMatrixValuesPSet.getParameter<std::vector<double> >("zAxis");
99  NPhi=1;
100  if (zBins.size()!=2){
101  edm::LogError(theCategory)<<"please specify a zAxis with 2 entries only. more bins not implemented yet.";}
102  minPhi=zBins.front();
103  maxPhi=zBins.back();
104 
105  }
106 
107  if (a==use){
108  if (!madeFromCff){
109  edm::LogInfo(theCategory)<<"using an error matrix object from: "<<fileName;
110  edm::FileInPath data(fileName);
111  std::string fullpath = data.fullPath();
112  gROOT->cd();
113  theD = new TFile(fullpath.c_str());
114  theD->SetWritable(false);
115  }else{
116  static std::atomic<unsigned int> neverTheSame{0};
117  std::stringstream dirName("MuonErrorMatrixDirectory");
118  dirName<<neverTheSame++;
119  edm::LogInfo(theCategory)<<"using an error matrix object from configuration file. putting memory histograms to: "<<dirName.str();
120  gROOT->cd();
121  theD = new TDirectory(dirName.str().c_str(),"transient directory to host MuonErrorMatrix TProfile3D");
122  theD->SetWritable(false);
123  }
124  }
125  else{
126  edm::LogInfo(theCategory)<<"creating an error matrix object: "<<fileName;
127  theD = new TFile(fileName.c_str(),"recreate");
128  }
129 
130  if (a==use && !madeFromCff ){gROOT->cd();}
131  else {theD->cd();}
132 
133  for (int i=0;i!=5;i++){for (int j=i;j!=5;j++){
134  TString pfname(Form("pf3_V%1d%1d",i+1,j+1));
135  TProfile3D * pf =0;
136  if (a==use && !madeFromCff ){
137  //read from the rootfile
138  edm::LogVerbatim(theCategory)<<"getting "<<pfname<<" from "<<fileName;
139  pf = (TProfile3D *)theD->Get(pfname);
140  theData[Pindex(i,j)]=pf;
141  theData_fast[i][j]=pf; theData_fast[j][i]=pf;
142  }
143  else{
144  // curvilinear coordinate system
145  //need to make some input parameter to be to change the number of bins
146 
147  TString pftitle;
148  if (i==j){pftitle="#sigma_{"+vars[i]+"}";}
149  else{pftitle="#rho("+vars[i]+","+vars[j]+")";}
150  edm::LogVerbatim(theCategory)<<"booking "<<pfname<<" into "<<fileName;
151  pf = new TProfile3D(pfname,pftitle,NPt,minPt,maxPt,NEta,minEta,maxEta,NPhi,minPhi,maxPhi,"S");
152  pf->SetXTitle("muon p_{T} [GeV]");
153  pf->SetYTitle("muon |#eta|");
154  pf->SetZTitle("muon #varphi");
155 
156  //set variable size binning
157  if (xBinsArray){
158  pf->GetXaxis()->Set(NPt,xBinsArray);}
159  if (yBinsArray){
160  pf->GetYaxis()->Set(NEta,yBinsArray);}
161 
162  if (madeFromCff){
163  edm::ParameterSet pSet = errorMatrixValuesPSet.getParameter<edm::ParameterSet>(pfname.Data());
164  //set the values from the configuration file itself
165  std::vector<double> values = pSet.getParameter<std::vector<double> >("values");
166  unsigned int iX=pf->GetNbinsX();
167  unsigned int iY=pf->GetNbinsY();
168  unsigned int iZ=pf->GetNbinsZ();
169  unsigned int continuous_i=0;
170  for(unsigned int ix=1;ix<=iX;++ix){
171  for(unsigned int iy=1;iy<=iY;++iy){
172  for(unsigned int iz=1;iz<=iZ;++iz){
173  LogTrace(theCategory)<<"filling profile:"
174  <<"\n pt (x)= "<<pf->GetXaxis()->GetBinCenter(ix)
175  <<"\n eta (y)= "<<pf->GetYaxis()->GetBinCenter(iy)
176  <<"\n phi (z)= "<<pf->GetZaxis()->GetBinCenter(iz)
177  <<"\n value= "<<values[continuous_i];
178  pf->Fill(pf->GetXaxis()->GetBinCenter(ix),
179  pf->GetYaxis()->GetBinCenter(iy),
180  pf->GetZaxis()->GetBinCenter(iz),
181  values[continuous_i++]);
182  }}}
183  //term action
184  std::string tAction = pSet.getParameter<std::string>("action");
185  if (tAction == "scale")
186  theTermAction[Pindex(i,j)] = scale;
187  else if (tAction == "assign")
189  else{
190  edm::LogError(theCategory)<<" wrong configuration: term action: "<<tAction<<" is not recognized.";
191  theTermAction[Pindex(i,j)] = error;
192  }
193 
194  }
195  }
196 
197  LogDebug(theCategory)<<" index "<<i<<":"<<j<<" -> "<<Pindex(i,j);
198  theData[Pindex(i,j)]=pf;
199  theData_fast[i][j]=pf; theData_fast[j][i]=pf;
200  if (!pf){
201  edm::LogError(theCategory)<<" profile "<<pfname<<" in file "<<fileName<<" is not valid. exiting.";
202  exit(1);
203  }
204  }}
205 
206  //verify it
207  for (int i=0;i!=15;i++){
208  if (theData[i]) {edm::LogVerbatim(theCategory)<<i<<" :"<<theData[i]->GetName()
209  <<" "<< theData[i]->GetTitle()<<std::endl;}}
210 
211 }
#define LogDebug(id)
const double Pi
T getParameter(std::string const &) const
static const TString vars[5]
names of the variables of the 5x5 error matrix
double maxEta
action
enum type to define if the class is used as a tool or to be created
TermAction theTermAction[15]
int Pindex(int i, int j)
internal methods to get the index of a matrix term.
#define LogTrace(id)
TProfile3D * theData_fast[5][5]
NEta
set NEta=0 and the vector of double for variable size binning
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
TDirectory * theD
the attached root file, where the parametrization is saved
double a
Definition: hdecay.h:121
TProfile3D * theData[15]
15 TProfile, each holding he parametrization of each term of the 5x5
std::string theCategory
log category: "MuonErrorMatrix"
minPhi
set NPt=0 and the vector of double for variable size binning
MuonErrorMatrix::~MuonErrorMatrix ( )

destructor

Definition at line 225 of file MuonErrorMatrix.cc.

References close().

225  {
226  close();
227 }
void close()
close the root file attached to the class

Member Function Documentation

void MuonErrorMatrix::adjust ( FreeTrajectoryState state)

adjust the error matrix on the state

Definition at line 446 of file MuonErrorMatrix.cc.

References assign, complicatedTerm(), FreeTrajectoryState::curvilinearError(), error, mps_fire::i, LogDebug, makeMuonMisalignmentScenario::matrix, FreeTrajectoryState::momentum(), FreeTrajectoryState::parameters(), Pindex(), scale, simpleTerm(), theCategory, and theTermAction.

Referenced by TSGForRoadSearch::IPfts(), and TSGForRoadSearch::notAtIPtsos().

446  {
447 
448  LogDebug(theCategory+"|Adjust")<<"state: \n"<<state;
449  AlgebraicSymMatrix55 simpleTerms;
450  simpleTerm(state.curvilinearError(), simpleTerms);
451  //the above contains sigma(i), rho(i,j)
452  LogDebug(theCategory+"|Adjust")<<"state sigma(i), rho(i,j): \n"<<simpleTerms;
453 
454  AlgebraicSymMatrix55 simpleValues=get(state.momentum(),false).matrix();
455  LogDebug(theCategory+"|Adjust")<<"config: (i,i), (i,j): \n"<<simpleValues;
456 
457  for(int i = 0;i!=5;i++){for(int j = i;j!=5;j++){
458  //check on each term for desired action
459  switch(theTermAction[Pindex(i,j)]){
460  case scale:
461  {
462  simpleTerms(i,j) *= simpleValues(i,j);
463  break;
464  }
465  case assign:
466  {
467  simpleTerms(i,j) = simpleValues(i,j);
468  break;
469  }
470  case error:
471  {
472  edm::LogError(theCategory+"|Adjust")<<" cannot properly adjust for term: "<<i <<","<<j;
473  }
474  }
475  }}
476  LogDebug(theCategory+"|Adjust")<<"updated state sigma(i), rho(i,j): \n"<<simpleTerms;
477 
478  AlgebraicSymMatrix55 finalTerms;
479  complicatedTerm(simpleTerms, finalTerms);
480  LogDebug(theCategory+"|Adjust")<<"updated state COV(i,j): \n"<<finalTerms;
481 
482  CurvilinearTrajectoryError oMat(finalTerms);
483  state = FreeTrajectoryState(state.parameters(), oMat);
484  LogDebug(theCategory+"|Adjust")<<"updated state:\n"<<state;
485 }
#define LogDebug(id)
const GlobalTrajectoryParameters & parameters() const
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
const CurvilinearTrajectoryError & curvilinearError() const
void complicatedTerm(const AlgebraicSymMatrix55 &input, AlgebraicSymMatrix55 &output)
convert sigma/rho -> sigma2/COV
TermAction theTermAction[15]
int Pindex(int i, int j)
internal methods to get the index of a matrix term.
GlobalVector momentum() const
void simpleTerm(const AlgebraicSymMatrix55 &input, AlgebraicSymMatrix55 &output)
convert sigma2/COV -> sigma/rho
std::string theCategory
log category: "MuonErrorMatrix"
void MuonErrorMatrix::adjust ( TrajectoryStateOnSurface state)

adjust the error matrix on the state

Definition at line 489 of file MuonErrorMatrix.cc.

References assign, complicatedTerm(), TrajectoryStateOnSurface::curvilinearError(), error, TrajectoryStateOnSurface::globalMomentum(), TrajectoryStateOnSurface::globalParameters(), mps_fire::i, LogDebug, makeMuonMisalignmentScenario::matrix, Pindex(), scale, simpleTerm(), TrajectoryStateOnSurface::surface(), TrajectoryStateOnSurface::surfaceSide(), theCategory, theTermAction, and TrajectoryStateOnSurface::weight().

489  {
490 
491  AlgebraicSymMatrix55 simpleTerms;
492  simpleTerm(state.curvilinearError(), simpleTerms);
493  LogDebug(theCategory+"|Adjust")<<"state sigma(i), rho(i,j): \n"<<simpleTerms;
494 
495  AlgebraicSymMatrix55 simpleValues= get(state.globalMomentum(),false).matrix();
496  LogDebug(theCategory+"|Adjust")<<"config: (i,i), (i,j):\n"<<simpleValues;
497 
498  for(int i = 0;i!=5;i++){for(int j = i;j!=5;j++){
499  //check on each term for desired action
500  switch(theTermAction[Pindex(i,j)]){
501  case scale:
502  {
503  simpleTerms(i,j) *= simpleValues(i,j);
504  break;
505  }
506  case assign:
507  {
508  simpleTerms(i,j) = simpleValues(i,j);
509  break;
510  }
511  case error:
512  {
513  edm::LogError(theCategory+"|Adjust")<<" cannot properly adjust for term: "<<i <<","<<j;
514  }
515  }
516  }}
517  LogDebug(theCategory+"|Adjust")<<"updated state sigma(i), rho(i,j): \n"<<simpleTerms;
518 
519  AlgebraicSymMatrix55 finalTerms;
520  complicatedTerm(simpleTerms, finalTerms);
521  LogDebug(theCategory+"|Adjust")<<"updated state COV(i,j): \n"<<finalTerms;
522 
523  CurvilinearTrajectoryError oMat(finalTerms);
524  state = TrajectoryStateOnSurface(state.weight(),
525  state.globalParameters(),
526  oMat,
527  state.surface(),
528  state.surfaceSide()
529  );
530  LogDebug(theCategory+"|Adjust")<<"updated state:\n"<<state;
531 }
#define LogDebug(id)
const CurvilinearTrajectoryError & curvilinearError() const
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
const SurfaceType & surface() const
void complicatedTerm(const AlgebraicSymMatrix55 &input, AlgebraicSymMatrix55 &output)
convert sigma/rho -> sigma2/COV
TermAction theTermAction[15]
SurfaceSide surfaceSide() const
Position relative to material, defined relative to momentum vector.
int Pindex(int i, int j)
internal methods to get the index of a matrix term.
const GlobalTrajectoryParameters & globalParameters() const
void simpleTerm(const AlgebraicSymMatrix55 &input, AlgebraicSymMatrix55 &output)
convert sigma2/COV -> sigma/rho
GlobalVector globalMomentum() const
std::string theCategory
log category: "MuonErrorMatrix"
void MuonErrorMatrix::close ( void  )

close the root file attached to the class

Definition at line 215 of file MuonErrorMatrix.cc.

References mps_fire::i, theD, and theData.

Referenced by lumiQTWidget.ApplicationWindow::fileQuit(), esMonitoring.AsyncLineReaderMixin::handle_close(), esMonitoring.FDJsonServer::handle_close(), Vispa.Gui.BoxContentDialog.BoxContentDialog::keyPressEvent(), Vispa.Gui.FindDialog.FindDialog::keyPressEvent(), and ~MuonErrorMatrix().

215  {
216  //close the file
217  if (theD){
218  theD->cd();
219  //write to it first if constructor
220  if (theD->IsWritable()){
221  for (int i=0;i!=15;i++){ if (theData[i]) { theData[i]->Write();}}}
222  theD->Close();
223  }}
TDirectory * theD
the attached root file, where the parametrization is saved
TProfile3D * theData[15]
15 TProfile, each holding he parametrization of each term of the 5x5
void MuonErrorMatrix::complicatedTerm ( const AlgebraicSymMatrix55 input,
AlgebraicSymMatrix55 output 
)

convert sigma/rho -> sigma2/COV

Definition at line 434 of file MuonErrorMatrix.cc.

References mps_fire::i, input, and convertSQLitetoXML_cfg::output.

Referenced by adjust().

434  {
435 
436  for(int i = 0;i!=5;i++){for(int j = i;j!=5;j++){
437  if (i==j)
438  output(i,j) = input(i,j)*input(i,j); //sigma squared
439  else
440  output(i,j) = input(i,j)*input(i,i)*input(j,j); //rho*sigma*sigma
441  }}
442 }
static std::string const input
Definition: EdmProvDump.cc:44
bool MuonErrorMatrix::divide ( CurvilinearTrajectoryError num_error,
const CurvilinearTrajectoryError denom_error 
)
static

divide term by term the two matrix

Definition at line 411 of file MuonErrorMatrix.cc.

References mps_fire::i, and CurvilinearTrajectoryError::matrix().

411  {
412  //divide term by term the matrix
413  const AlgebraicSymMatrix55 & denom_matrix=denom_error.matrix();
414  AlgebraicSymMatrix55 num_matrix = num_error.matrix();
415  // the true type of the matrix is such that [i][j] is the same memory object as [j][i]: looping i:0-5, j:0-5 double multiply the terms
416  // need to loop only on i:0-5, j:i-5
417  for(int i = 0;i!=5;i++){for(int j = i;j!=5;j++){
418  if (denom_matrix(i,j)==0) return false;
419  num_matrix(i,j)/=denom_matrix(i,j);
420  }}
421  num_error = CurvilinearTrajectoryError(num_matrix);
422  return true;
423 }
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
const AlgebraicSymMatrix55 & matrix() const
int MuonErrorMatrix::findBin ( TAxis *  axis,
double  value 
)

method to get the bin index, taking care of under/overlow: first(1)/last(GetNbins())returned

Definition at line 289 of file MuonErrorMatrix.cc.

References mps_fire::result.

Referenced by getFast(), and Value().

289  {
290  //find the proper bin, protecting against under/over flow
291  int result = axis->FindBin(value);
292  if (result <= 0) result=1; //protect against under flow
293  else if (result > axis->GetNbins() ) result = axis->GetNbins();
294  return result;}
Definition: value.py:1
CurvilinearTrajectoryError MuonErrorMatrix::get ( GlobalVector  momentum,
bool  convolute = true 
)

main method to be used. Retrieve a 5x5 symetrical matrix according to parametrization of error or scale factor

Definition at line 231 of file MuonErrorMatrix.cc.

References mps_fire::i, and Value().

Referenced by Options.Options::__getitem__(), betterConfigParser.BetterConfigParser::__updateDict(), TSGFromPropagation::adjust(), MuonErrorMatrixAdjuster::fix_cov_matrix(), betterConfigParser.BetterConfigParser::getCompares(), betterConfigParser.BetterConfigParser::getGeneral(), betterConfigParser.BetterConfigParser::getResultingSection(), and TSGFromPropagation::init().

231  {
233  //retrieves a 55 matrix containing (i,i)^2 and (i,j)*(i,i)*(j,j)
234  for (int i=0;i!=5;i++){for (int j=i;j!=5;j++){
235  V(i,j) = Value(momentum,i,j,convolute);}}
236  return CurvilinearTrajectoryError(V);}
double Value(GlobalVector &momentum, int i, int j, bool convolute=true)
internal method that retreives the value of the parametrization for term i,j
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
TProfile3D* MuonErrorMatrix::get ( int  i,
int  j 
)
inline

actually get access to the TProfile3D used for the parametrization

Definition at line 53 of file MuonErrorMatrix.h.

References Index().

Referenced by Options.Options::__getitem__(), betterConfigParser.BetterConfigParser::__updateDict(), betterConfigParser.BetterConfigParser::getCompares(), betterConfigParser.BetterConfigParser::getGeneral(), and betterConfigParser.BetterConfigParser::getResultingSection().

53 {return Index(i,j);}
TProfile3D * Index(int i, int j)
internal method to get access to the profiles
CurvilinearTrajectoryError MuonErrorMatrix::getFast ( GlobalVector  momentum)

Definition at line 238 of file MuonErrorMatrix.cc.

References stringResolutionProvider_cfi::eta, PV3DBase< T, PVType, FrameType >::eta(), findBin(), mps_fire::i, PV3DBase< T, PVType, FrameType >::perp(), phi, PV3DBase< T, PVType, FrameType >::phi(), theData_fast, and MuonErrorMatrixValues_cff::values.

238  {
239  //will be faster but make assumptions that could be broken at some point
240  // same bining for all TProfile
242 
243  double pT = momentum.perp();
244  double eta = fabs(momentum.eta());
245  double phi = momentum.phi();
246 
247  //assume all the same axis in X,Y,Z
248  int iBin_x= findBin(theData_fast[0][0]->GetXaxis(),pT);
249  int iBin_y= findBin(theData_fast[0][0]->GetYaxis(),eta);
250  int iBin_z= findBin(theData_fast[0][0]->GetZaxis(),phi);
251 
252  //retreive values
253  double values[5][5]; //sigma_i and rho_ij
254  for (int i=0;i!=5;i++){for (int j=i;j!=5;j++){
255  values[i][j]=theData_fast[i][j]->GetBinContent(iBin_x, iBin_y, iBin_z);
256  }}
257 
258  for (int i=0;i!=5;i++){for (int j=i;j!=5;j++){
259  if (i==j){
260  //sigma_i * sigma_i
261  V(i,j) = values[i][j];
262  V(i,j)*=V(i,j);
263  }
264  else{
265  //sigma_i * sigma_j * rho_ij
266  V(i,j) = values[i][i] * values[j][j] * values[i][j];
267  }
268  }}
269 
270  return CurvilinearTrajectoryError(V);}
T perp() const
Definition: PV3DBase.h:72
int findBin(TAxis *axis, double value)
method to get the bin index, taking care of under/overlow: first(1)/last(GetNbins())returned ...
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
TProfile3D * theData_fast[5][5]
T eta() const
Definition: PV3DBase.h:76
unsigned int MuonErrorMatrix::index ( int  i,
int  j 
)
inline

Definition at line 54 of file MuonErrorMatrix.h.

References Pindex().

Referenced by BeautifulSoup.PageElement::insert().

54 {return Pindex(i,j);}
int Pindex(int i, int j)
internal methods to get the index of a matrix term.
TProfile3D* MuonErrorMatrix::Index ( int  i,
int  j 
)
inlineprivate

internal method to get access to the profiles

Definition at line 96 of file MuonErrorMatrix.h.

References Pindex(), Rms(), and Value().

Referenced by get(), Rms(), and Value().

96  {
97  return theData[Pindex(i,j)];}
int Pindex(int i, int j)
internal methods to get the index of a matrix term.
TProfile3D * theData[15]
15 TProfile, each holding he parametrization of each term of the 5x5
void MuonErrorMatrix::multiply ( CurvilinearTrajectoryError initial_error,
const CurvilinearTrajectoryError scale_error 
)
static

multiply term by term the two matrix

Definition at line 400 of file MuonErrorMatrix.cc.

References mps_fire::i, and CurvilinearTrajectoryError::matrix().

Referenced by TSGFromPropagation::adjust(), and FastTSGFromPropagation::adjust().

400  {
401  //scale term by term the matrix
402  const AlgebraicSymMatrix55 & scale_matrix=scale_error.matrix();
403  AlgebraicSymMatrix55 revised_matrix = initial_error.matrix();
404  // the true type of the matrix is such that [i][j] is the same memory object as [j][i]: looping i:0-5, j:0-5 double multiply the terms
405  // need to loop only on i:0-5, j:i-5
406  for(int i = 0;i!=5;i++){for(int j = i;j!=5;j++){
407  revised_matrix(i,j)*=scale_matrix(i,j);
408  }}
409  initial_error = CurvilinearTrajectoryError(revised_matrix);
410 }
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
const AlgebraicSymMatrix55 & matrix() const
int MuonErrorMatrix::Pindex ( int  i,
int  j 
)
inlineprivate

internal methods to get the index of a matrix term.

Definition at line 92 of file MuonErrorMatrix.h.

References funct::abs(), mps_fire::i, and PFRecoTauDiscriminationByIsolation_cfi::offset.

Referenced by adjust(), index(), Index(), and MuonErrorMatrix().

92  {
93  static const int offset[5]={0,5,5+4,5+4+3,5+4+3+2};
94  return offset[i]+abs(j-i);}
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double MuonErrorMatrix::Rms ( GlobalVector momentum,
int  i,
int  j 
)
private

internal method that retreives the error on the value of the parametrization for term i,j

Definition at line 354 of file MuonErrorMatrix.cc.

References stringResolutionProvider_cfi::eta, PV3DBase< T, PVType, FrameType >::eta(), Index(), LogDebug, PV3DBase< T, PVType, FrameType >::perp(), phi, PV3DBase< T, PVType, FrameType >::phi(), mps_fire::result, and theCategory.

Referenced by Index().

354  {
355  double result=0;
356  TProfile3D * ij = Index(i,j);
357  if (!ij){edm::LogError(theCategory)<<"cannot get the profile ("<<i<<":"<<i<<")"; return result;}
358  double pT = momentum.perp();
359  double eta = fabs(momentum.eta());
360  double phi = momentum.phi();
361 
362  int iBin_x= ij->GetXaxis()->FindBin(pT);
363  int iBin_y= ij->GetYaxis()->FindBin(eta);
364  int iBin_z= ij->GetZaxis()->FindBin(phi);
365  result=ij->GetBinError(iBin_x,iBin_y,iBin_z);
366 
367  LogDebug(theCategory)<<"for: (pT,eta,phi)=("<<pT<<", "<<eta<<", "<<phi<<") error["<<i<<","<<j<<"] is: "<<result;
368  return result;
369 }
#define LogDebug(id)
T perp() const
Definition: PV3DBase.h:72
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T eta() const
Definition: PV3DBase.h:76
std::string theCategory
log category: "MuonErrorMatrix"
TProfile3D * Index(int i, int j)
internal method to get access to the profiles
void MuonErrorMatrix::simpleTerm ( const AlgebraicSymMatrix55 input,
AlgebraicSymMatrix55 output 
)

convert sigma2/COV -> sigma/rho

Definition at line 425 of file MuonErrorMatrix.cc.

References mps_fire::i, input, convertSQLitetoXML_cfg::output, and mathSSE::sqrt().

Referenced by adjust().

425  {
426  for(int i = 0;i!=5;i++){for(int j = i;j!=5;j++){
427  if (i==j)
428  output(i,j) = sqrt(input(i,j)); //sigma
429  else
430  output(i,j) = input(i,j) / sqrt( input(i,i)* input(j,j)); //rho
431  }}
432 }
static std::string const input
Definition: EdmProvDump.cc:44
T sqrt(T t)
Definition: SSEVec.h:18
double MuonErrorMatrix::Term ( const AlgebraicSymMatrix55 curv,
int  i,
int  j 
)
static

provide the numerical value used. sigma or correlation factor

Definition at line 371 of file MuonErrorMatrix.cc.

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

371  {
372  //return sigma or correlation factor
373  double result=0;
374  if (i==j){
375  result = curv(i,j);
376  if (result<0){
377  //check validity of this guy
378  edm::LogError("MuonErrorMatrix")<<"invalid term in the error matrix.\n sii: "<< result;
379  return 0;}
380  return sqrt(result);}
381  else{
382  double si = curv(i,i);
383  double sj = curv(j,j);
384 
385  if (si<=0 || sj<=0){
386  edm::LogError("MuonErrorMatrix")<<"invalid term in the error matrix.\n si: "
387  <<si<<" sj: "<<sj<<". result will be corrupted\n"<<curv;
388  return 0;}
389 
390  si=sqrt(si);
391  sj=sqrt(sj);
392  //check validity
393 
394  return result = curv(i,j)/(si*sj);}
395  //by default
396  return 0;
397 }
T sqrt(T t)
Definition: SSEVec.h:18
double MuonErrorMatrix::Value ( GlobalVector momentum,
int  i,
int  j,
bool  convolute = true 
)
private

internal method that retreives the value of the parametrization for term i,j

Definition at line 297 of file MuonErrorMatrix.cc.

References corr, stringResolutionProvider_cfi::eta, PV3DBase< T, PVType, FrameType >::eta(), findBin(), cuy::ii, Index(), findQualityFiles::jj, LogDebug, PV3DBase< T, PVType, FrameType >::perp(), phi, PV3DBase< T, PVType, FrameType >::phi(), mps_fire::result, and theCategory.

Referenced by get(), and Index().

297  {
298  double result=0;
299  TProfile3D * ij = Index(i,j);
300  if (!ij) {edm::LogError(theCategory)<<"cannot get the profile ("<<i<<":"<<j<<")"; return result;}
301 
302  double pT = momentum.perp();
303  double eta = fabs(momentum.eta());
304  double phi = momentum.phi();
305 
306  int iBin_x= findBin(ij->GetXaxis(),pT);
307  int iBin_y= findBin(ij->GetYaxis(),eta);
308  int iBin_z= findBin(ij->GetZaxis(),phi);
309 
310  if (convolute){
311  if (i!=j){
312  //return the covariance = correlation*sigma_1 *sigma_2;
313  TProfile3D * ii = Index(i,i);
314  TProfile3D * jj = Index(j,j);
315  if (!ii){edm::LogError(theCategory)<<"cannot get the profile ("<<i<<":"<<i<<")"; return result;}
316  if (!jj){edm::LogError(theCategory)<<"cannot get the profile ("<<j<<":"<<j<<")"; return result;}
317 
318 
319  int iBin_i_x = findBin(ii->GetXaxis(),pT);
320  int iBin_i_y = findBin(ii->GetYaxis(),eta);
321  int iBin_i_z = findBin(ii->GetZaxis(),phi);
322  int iBin_j_x = findBin(jj->GetXaxis(),pT);
323  int iBin_j_y = findBin(jj->GetYaxis(),eta);
324  int iBin_j_z = findBin(jj->GetZaxis(),phi);
325 
326  double corr = ij->GetBinContent(iBin_x,iBin_y,iBin_z);
327  double sigma_1 = (ii->GetBinContent(iBin_i_x,iBin_i_y,iBin_i_z));
328  double sigma_2 = (jj->GetBinContent(iBin_j_x,iBin_j_y,iBin_j_z));
329 
330 
331  result=corr*sigma_1*sigma_2;
332  LogDebug(theCategory)<<"for: (pT,eta,phi)=("<<pT<<", "<<eta<<", "<<phi<<") nterms are"
333  <<"\nrho["<<i<<","<<j<<"]: "<<corr<<" ["<< iBin_x<<", "<<iBin_y<<", "<<iBin_z<<"]"
334  <<"\nsigma["<<i<<","<<i<<"]: "<< sigma_1
335  <<"\nsigma["<<j<<","<<j<<"]: "<< sigma_2
336  <<"Covariance["<<i<<","<<j<<"] is: "<<result;
337  return result;
338  }
339  else{
340  //return the variance = sigma_1 **2
341  // result=ij->GetBinContent(iBin);
342  result=ij->GetBinContent(iBin_x,iBin_y,iBin_z);
343  result*=result;
344  LogDebug(theCategory)<<"for: (pT,eta,phi)=("<<pT<<", "<<eta<<", "<<phi<<") sigma^2["<<i<<","<<j<<"] is: "<<result;
345  return result;
346  }
347  }else{
348  //do not convolute
349  result=ij->GetBinContent(iBin_x,iBin_y,iBin_z);
350  return result;
351  }
352 }
#define LogDebug(id)
T perp() const
Definition: PV3DBase.h:72
int findBin(TAxis *axis, double value)
method to get the bin index, taking care of under/overlow: first(1)/last(GetNbins())returned ...
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
JetCorrectorParameters corr
Definition: classes.h:5
ii
Definition: cuy.py:588
T eta() const
Definition: PV3DBase.h:76
std::string theCategory
log category: "MuonErrorMatrix"
TProfile3D * Index(int i, int j)
internal method to get access to the profiles

Member Data Documentation

std::string MuonErrorMatrix::theCategory
private

log category: "MuonErrorMatrix"

Definition at line 79 of file MuonErrorMatrix.h.

Referenced by adjust(), MuonErrorMatrix(), Rms(), and Value().

TDirectory* MuonErrorMatrix::theD
private

the attached root file, where the parametrization is saved

Definition at line 82 of file MuonErrorMatrix.h.

Referenced by close(), and MuonErrorMatrix().

TProfile3D* MuonErrorMatrix::theData[15]
private

15 TProfile, each holding he parametrization of each term of the 5x5

Definition at line 84 of file MuonErrorMatrix.h.

Referenced by close(), and MuonErrorMatrix().

TProfile3D* MuonErrorMatrix::theData_fast[5][5]
private

Definition at line 85 of file MuonErrorMatrix.h.

Referenced by getFast(), and MuonErrorMatrix().

TermAction MuonErrorMatrix::theTermAction[15]
private

Definition at line 89 of file MuonErrorMatrix.h.

Referenced by adjust(), and MuonErrorMatrix().

const TString MuonErrorMatrix::vars ={"#frac{q}{|p|}","#lambda","#varphi_{0}","X_{T}","Y_{T}"}
static

names of the variables of the 5x5 error matrix

Definition at line 57 of file MuonErrorMatrix.h.

Referenced by MuonErrorMatrix().