CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MuonErrorMatrix.cc
Go to the documentation of this file.
4 
5 #include "TROOT.h"
6 #include "TString.h"
7 #include "TRandom2.h"
8 #include "TMath.h"
9 
10 #include <sstream>
11 
12 using namespace std;
13 
14 const TString MuonErrorMatrix::vars[5]={"#frac{q}{|p|}","#lambda","#varphi_{0}","X_{T}","Y_{T}"};
15 
17  theCategory="MuonErrorMatrix";
18  std::string action = iConfig.getParameter<std::string>("action");
19 
20  bool madeFromCff=iConfig.exists("errorMatrixValuesPSet");
21  edm::ParameterSet errorMatrixValuesPSet;
22 
24  if (!madeFromCff){ fileName = iConfig.getParameter<std::string>("rootFileName");}
25  else { errorMatrixValuesPSet =iConfig.getParameter<edm::ParameterSet>("errorMatrixValuesPSet");}
27 
28  int NPt=5;
29  std::vector<double> xBins;
30  double * xBinsArray = 0;
31  double minPt=1;
32  double maxPt=200;
33  int NEta=5;
34  std::vector<double> yBins;
35  double * yBinsArray =0;
36  double minEta=0;
37  double maxEta=2.5;
38  int NPhi=1;
39  double minPhi=-TMath::Pi();
40  double maxPhi=TMath::Pi();
41 
42  if (action!="use"){
43  a = constructor;
44 
45  NPt = iConfig.getParameter<int>("NPt");
46  if (NPt!=0){
47  minPt = iConfig.getParameter<double>("minPt");
48  maxPt = iConfig.getParameter<double>("maxPt");}
49  else{
50  xBins = iConfig.getParameter<std::vector<double> >("PtBins");
51  if (xBins.size()==0){edm::LogError( theCategory)<<"Npt=0 and no entries in the vector. I will do aseg fault soon.";}
52  NPt = xBins.size()-1;
53  xBinsArray = &(xBins.front());
54  minPt = xBins.front();
55  maxPt = xBins.back();}
56 
57  NEta = iConfig.getParameter<int>("NEta");
58  if (NEta!=0){
59  minEta = iConfig.getParameter<double>("minEta");
60  maxEta = iConfig.getParameter<double>("maxEta");}
61  else{
62  yBins = iConfig.getParameter<std::vector<double> >("EtaBins");
63  if (yBins.size()==0){edm::LogError( theCategory)<<"NEta=0 and no entries in the vector. I will do aseg fault soon.";}
64  NEta = yBins.size()-1;
65  yBinsArray = &(yBins.front());
66  minPt = yBins.front();
67  maxPt = yBins.back();}
68 
69  NPhi = iConfig.getParameter<int>("NPhi");
70  std::stringstream get(iConfig.getParameter<std::string>("minPhi"));
71  if (get.str() =="-Pi")
72  { minPhi =-TMath::Pi();}
73  else if(get.str() =="Pi")
74  { minPhi =TMath::Pi();}
75  else { get>>minPhi;}
76  get.str(iConfig.getParameter<std::string>("maxPhi"));
77  if (get.str() =="-Pi")
78  { maxPhi =-TMath::Pi();}
79  else if(get.str() =="Pi")
80  { maxPhi =TMath::Pi();}
81  else { get>>maxPhi;}
82  }//action!=use
83  else if (madeFromCff){
84 
85  xBins = errorMatrixValuesPSet.getParameter<std::vector<double> >("xAxis");
86  NPt = xBins.size()-1;
87  xBinsArray = &(xBins.front());
88  minPt = xBins.front();
89  maxPt = xBins.back();
90 
91  yBins = errorMatrixValuesPSet.getParameter<std::vector<double> >("yAxis");
92  NEta = yBins.size()-1;
93  yBinsArray = &(yBins.front());
94  minPt = yBins.front();
95  maxPt = yBins.back();
96 
97  std::vector<double> zBins = errorMatrixValuesPSet.getParameter<std::vector<double> >("zAxis");
98  NPhi=1;
99  if (zBins.size()!=2){
100  edm::LogError(theCategory)<<"please specify a zAxis with 2 entries only. more bins not implemented yet.";}
101  minPhi=zBins.front();
102  maxPhi=zBins.back();
103 
104  }
105 
106  if (a==use){
107  if (!madeFromCff){
108  edm::LogInfo(theCategory)<<"using an error matrix object from: "<<fileName;
109  edm::FileInPath data(fileName);
110  std::string fullpath = data.fullPath();
111  gROOT->cd();
112  theD = new TFile(fullpath.c_str());
113  theD->SetWritable(false);
114  }else{
115  static unsigned int neverTheSame=0;
116  std::stringstream dirName("MuonErrorMatrixDirectory");
117  dirName<<neverTheSame++;
118  edm::LogInfo(theCategory)<<"using an error matrix object from configuration file. putting memory histograms to: "<<dirName.str();
119  gROOT->cd();
120  theD = new TDirectory(dirName.str().c_str(),"transient directory to host MuonErrorMatrix TProfile3D");
121  theD->SetWritable(false);
122  }
123  }
124  else{
125  edm::LogInfo(theCategory)<<"creating an error matrix object: "<<fileName;
126  theD = new TFile(fileName.c_str(),"recreate");
127  }
128 
129  if (a==use && !madeFromCff ){gROOT->cd();}
130  else {theD->cd();}
131 
132  for (int i=0;i!=5;i++){for (int j=i;j!=5;j++){
133  TString pfname(Form("pf3_V%1d%1d",i+1,j+1));
134  TProfile3D * pf =0;
135  if (a==use && !madeFromCff ){
136  //read from the rootfile
137  edm::LogVerbatim(theCategory)<<"getting "<<pfname<<" from "<<fileName;
138  pf = (TProfile3D *)theD->Get(pfname);
139  theData[Pindex(i,j)]=pf;
140  theData_fast[i][j]=pf; theData_fast[j][i]=pf;
141  }
142  else{
143  // curvilinear coordinate system
144  //need to make some input parameter to be to change the number of bins
145 
146  TString pftitle;
147  if (i==j){pftitle="#sigma_{"+vars[i]+"}";}
148  else{pftitle="#rho("+vars[i]+","+vars[j]+")";}
149  edm::LogVerbatim(theCategory)<<"booking "<<pfname<<" into "<<fileName;
150  pf = new TProfile3D(pfname,pftitle,NPt,minPt,maxPt,NEta,minEta,maxEta,NPhi,minPhi,maxPhi,"S");
151  pf->SetXTitle("muon p_{T} [GeV]");
152  pf->SetYTitle("muon |#eta|");
153  pf->SetZTitle("muon #varphi");
154 
155  //set variable size binning
156  if (xBinsArray){
157  pf->GetXaxis()->Set(NPt,xBinsArray);}
158  if (yBinsArray){
159  pf->GetYaxis()->Set(NEta,yBinsArray);}
160 
161  if (madeFromCff){
162  edm::ParameterSet pSet = errorMatrixValuesPSet.getParameter<edm::ParameterSet>(pfname.Data());
163  //set the values from the configuration file itself
164  std::vector<double> values = pSet.getParameter<std::vector<double> >("values");
165  unsigned int iX=pf->GetNbinsX();
166  unsigned int iY=pf->GetNbinsY();
167  unsigned int iZ=pf->GetNbinsZ();
168  unsigned int continuous_i=0;
169  for(unsigned int ix=1;ix<=iX;++ix){
170  for(unsigned int iy=1;iy<=iY;++iy){
171  for(unsigned int iz=1;iz<=iZ;++iz){
172  LogTrace(theCategory)<<"filling profile:"
173  <<"\n pt (x)= "<<pf->GetXaxis()->GetBinCenter(ix)
174  <<"\n eta (y)= "<<pf->GetYaxis()->GetBinCenter(iy)
175  <<"\n phi (z)= "<<pf->GetZaxis()->GetBinCenter(iz)
176  <<"\n value= "<<values[continuous_i];
177  pf->Fill(pf->GetXaxis()->GetBinCenter(ix),
178  pf->GetYaxis()->GetBinCenter(iy),
179  pf->GetZaxis()->GetBinCenter(iz),
180  values[continuous_i++]);
181  }}}
182  //term action
183  std::string tAction = pSet.getParameter<std::string>("action");
184  if (tAction == "scale")
186  else if (tAction == "assign")
188  else{
189  edm::LogError(theCategory)<<" wrong configuration: term action: "<<tAction<<" is not recognized.";
191  }
192 
193  }
194  }
195 
196  LogDebug(theCategory)<<" index "<<i<<":"<<j<<" -> "<<Pindex(i,j);
197  theData[Pindex(i,j)]=pf;
198  theData_fast[i][j]=pf; theData_fast[j][i]=pf;
199  if (!pf){
200  edm::LogError(theCategory)<<" profile "<<pfname<<" in file "<<fileName<<" is not valid. exiting.";
201  exit(1);
202  }
203  }}
204 
205  //verify it
206  for (int i=0;i!=15;i++){
207  if (theData[i]) {edm::LogVerbatim(theCategory)<<i<<" :"<<theData[i]->GetName()
208  <<" "<< theData[i]->GetTitle()<<std::endl;}}
209 
210 }
211 
212 //void MuonErrorMatrix::writeIntoCFF(){}
213 
215  //close the file
216  if (theD){
217  theD->cd();
218  //write to it first if constructor
219  if (theD->IsWritable()){
220  for (int i=0;i!=15;i++){ if (theData[i]) { theData[i]->Write();}}}
221  theD->Close();
222  }}
223 
225  close();
226 }
227 
228 
229 
232  //retrieves a 55 matrix containing (i,i)^2 and (i,j)*(i,i)*(j,j)
233  for (int i=0;i!=5;i++){for (int j=i;j!=5;j++){
234  V(i,j) = Value(momentum,i,j,convolute);}}
235  return CurvilinearTrajectoryError(V);}
236 
238  //will be faster but make assumptions that could be broken at some point
239  // same bining for all TProfile
241 
242  double pT = momentum.perp();
243  double eta = fabs(momentum.eta());
244  double phi = momentum.phi();
245 
246  //assume all the same axis in X,Y,Z
247  int iBin_x= findBin(theData_fast[0][0]->GetXaxis(),pT);
248  int iBin_y= findBin(theData_fast[0][0]->GetYaxis(),eta);
249  int iBin_z= findBin(theData_fast[0][0]->GetZaxis(),phi);
250 
251  //retreive values
252  double values[5][5]; //sigma_i and rho_ij
253  for (int i=0;i!=5;i++){for (int j=i;j!=5;j++){
254  values[i][j]=theData_fast[i][j]->GetBinContent(iBin_x, iBin_y, iBin_z);
255  }}
256 
257  for (int i=0;i!=5;i++){for (int j=i;j!=5;j++){
258  if (i==j){
259  //sigma_i * sigma_i
260  V(i,j) = values[i][j];
261  V(i,j)*=V(i,j);
262  }
263  else{
264  //sigma_i * sigma_j * rho_ij
265  V(i,j) = values[i][i] * values[j][j] * values[i][j];
266  }
267  }}
268 
269  return CurvilinearTrajectoryError(V);}
270 
271 
272 
273 /*CurvilinearTrajectoryError MuonErrorMatrix::get_random(GlobalVector momentum) {
274  static TRandom2 rand;
275  AlgebraicSymMatrix55 V;//result
276  //first proceed with diagonal elements
277  for (int i=0;i!=5;i++){
278  V(i,i)=rand.Gaus( Value(momentum,i,i), Rms(momentum,i,i));}
279 
280  //now proceed with the correlations
281  for (int i=0;i!=5;i++){for (int j=i+1;j<5;j++){
282  double corr = rand.Gaus( Value(momentum,i,j), Rms(momentum,i,j));
283  //assign the covariance from correlation and sigmas
284  V(i,j)= corr * sqrt( V[i][i] * V[j][j]);}}
285  return CurvilinearTrajectoryError(V); }
286 */
287 
288 int MuonErrorMatrix::findBin(TAxis * axis, double value){
289  //find the proper bin, protecting against under/over flow
290  int result = axis->FindBin(value);
291  if (result <= 0) result=1; //protect against under flow
292  else if (result > axis->GetNbins() ) result = axis->GetNbins();
293  return result;}
294 
295 
296 double MuonErrorMatrix::Value(GlobalVector & momentum, int i, int j,bool convolute) {
297  double result=0;
298  TProfile3D * ij = Index(i,j);
299  if (!ij) {edm::LogError(theCategory)<<"cannot get the profile ("<<i<<":"<<j<<")"; return result;}
300 
301  double pT = momentum.perp();
302  double eta = fabs(momentum.eta());
303  double phi = momentum.phi();
304 
305  int iBin_x= findBin(ij->GetXaxis(),pT);
306  int iBin_y= findBin(ij->GetYaxis(),eta);
307  int iBin_z= findBin(ij->GetZaxis(),phi);
308 
309  if (convolute){
310  if (i!=j){
311  //return the covariance = correlation*sigma_1 *sigma_2;
312  TProfile3D * ii = Index(i,i);
313  TProfile3D * jj = Index(j,j);
314  if (!ii){edm::LogError(theCategory)<<"cannot get the profile ("<<i<<":"<<i<<")"; return result;}
315  if (!jj){edm::LogError(theCategory)<<"cannot get the profile ("<<j<<":"<<j<<")"; return result;}
316 
317 
318  int iBin_i_x = findBin(ii->GetXaxis(),pT);
319  int iBin_i_y = findBin(ii->GetYaxis(),eta);
320  int iBin_i_z = findBin(ii->GetZaxis(),phi);
321  int iBin_j_x = findBin(jj->GetXaxis(),pT);
322  int iBin_j_y = findBin(jj->GetYaxis(),eta);
323  int iBin_j_z = findBin(jj->GetZaxis(),phi);
324 
325  double corr = ij->GetBinContent(iBin_x,iBin_y,iBin_z);
326  double sigma_1 = (ii->GetBinContent(iBin_i_x,iBin_i_y,iBin_i_z));
327  double sigma_2 = (jj->GetBinContent(iBin_j_x,iBin_j_y,iBin_j_z));
328 
329 
330  result=corr*sigma_1*sigma_2;
331  LogDebug(theCategory)<<"for: (pT,eta,phi)=("<<pT<<", "<<eta<<", "<<phi<<") nterms are"
332  <<"\nrho["<<i<<","<<j<<"]: "<<corr<<" ["<< iBin_x<<", "<<iBin_y<<", "<<iBin_z<<"]"
333  <<"\nsigma["<<i<<","<<i<<"]: "<< sigma_1
334  <<"\nsigma["<<j<<","<<j<<"]: "<< sigma_2
335  <<"Covariance["<<i<<","<<j<<"] is: "<<result;
336  return result;
337  }
338  else{
339  //return the variance = sigma_1 **2
340  // result=ij->GetBinContent(iBin);
341  result=ij->GetBinContent(iBin_x,iBin_y,iBin_z);
342  result*=result;
343  LogDebug(theCategory)<<"for: (pT,eta,phi)=("<<pT<<", "<<eta<<", "<<phi<<") sigma^2["<<i<<","<<j<<"] is: "<<result;
344  return result;
345  }
346  }else{
347  //do not convolute
348  result=ij->GetBinContent(iBin_x,iBin_y,iBin_z);
349  return result;
350  }
351 }
352 
353 double MuonErrorMatrix::Rms(GlobalVector & momentum, int i, int j) {
354  double result=0;
355  TProfile3D * ij = Index(i,j);
356  if (!ij){edm::LogError(theCategory)<<"cannot get the profile ("<<i<<":"<<i<<")"; return result;}
357  double pT = momentum.perp();
358  double eta = fabs(momentum.eta());
359  double phi = momentum.phi();
360 
361  int iBin_x= ij->GetXaxis()->FindBin(pT);
362  int iBin_y= ij->GetYaxis()->FindBin(eta);
363  int iBin_z= ij->GetZaxis()->FindBin(phi);
364  result=ij->GetBinError(iBin_x,iBin_y,iBin_z);
365 
366  LogDebug(theCategory)<<"for: (pT,eta,phi)=("<<pT<<", "<<eta<<", "<<phi<<") error["<<i<<","<<j<<"] is: "<<result;
367  return result;
368 }
369 
370 double MuonErrorMatrix::Term(const AlgebraicSymMatrix55 & curv, int i, int j){
371  //return sigma or correlation factor
372  double result=0;
373  if (i==j){
374  result = curv(i,j);
375  if (result<0){
376  //check validity of this guy
377  edm::LogError("MuonErrorMatrix")<<"invalid term in the error matrix.\n sii: "<< result;
378  return 0;}
379  return sqrt(result);}
380  else{
381  double si = curv(i,i);
382  double sj = curv(j,j);
383 
384  if (si<=0 || sj<=0){
385  edm::LogError("MuonErrorMatrix")<<"invalid term in the error matrix.\n si: "
386  <<si<<" sj: "<<sj<<". result will be corrupted\n"<<curv;
387  return 0;}
388 
389  si=sqrt(si);
390  sj=sqrt(sj);
391  //check validity
392 
393  return result = curv(i,j)/(si*sj);}
394  //by default
395  return 0;
396 }
397 
398 
400  //scale term by term the matrix
401  const AlgebraicSymMatrix55 & scale_matrix=scale_error.matrix();
402  AlgebraicSymMatrix55 revised_matrix = initial_error.matrix();
403  // 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
404  // need to loop only on i:0-5, j:i-5
405  for(int i = 0;i!=5;i++){for(int j = i;j!=5;j++){
406  revised_matrix(i,j)*=scale_matrix(i,j);
407  }}
408  initial_error = CurvilinearTrajectoryError(revised_matrix);
409 }
411  //divide term by term the matrix
412  const AlgebraicSymMatrix55 & denom_matrix=denom_error.matrix();
413  AlgebraicSymMatrix55 num_matrix = num_error.matrix();
414  // 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
415  // need to loop only on i:0-5, j:i-5
416  for(int i = 0;i!=5;i++){for(int j = i;j!=5;j++){
417  if (denom_matrix(i,j)==0) return false;
418  num_matrix(i,j)/=denom_matrix(i,j);
419  }}
420  num_error = CurvilinearTrajectoryError(num_matrix);
421  return true;
422 }
423 
425  for(int i = 0;i!=5;i++){for(int j = i;j!=5;j++){
426  if (i==j)
427  output(i,j) = sqrt(input(i,j)); //sigma
428  else
429  output(i,j) = input(i,j) / sqrt( input(i,i)* input(j,j)); //rho
430  }}
431 }
432 
434 
435  for(int i = 0;i!=5;i++){for(int j = i;j!=5;j++){
436  if (i==j)
437  output(i,j) = input(i,j)*input(i,j); //sigma squared
438  else
439  output(i,j) = input(i,j)*input(i,i)*input(j,j); //rho*sigma*sigma
440  }}
441 }
442 
443 
444 
446 
447  LogDebug(theCategory+"|Adjust")<<"state: \n"<<state;
448  AlgebraicSymMatrix55 simpleTerms;
449  simpleTerm(state.curvilinearError(), simpleTerms);
450  //the above contains sigma(i), rho(i,j)
451  LogDebug(theCategory+"|Adjust")<<"state sigma(i), rho(i,j): \n"<<simpleTerms;
452 
453  AlgebraicSymMatrix55 simpleValues=get(state.momentum(),false).matrix();
454  LogDebug(theCategory+"|Adjust")<<"config: (i,i), (i,j): \n"<<simpleValues;
455 
456  for(int i = 0;i!=5;i++){for(int j = i;j!=5;j++){
457  //check on each term for desired action
458  switch(theTermAction[Pindex(i,j)]){
459  case scale:
460  {
461  simpleTerms(i,j) *= simpleValues(i,j);
462  break;
463  }
464  case assign:
465  {
466  simpleTerms(i,j) = simpleValues(i,j);
467  break;
468  }
469  case error:
470  {
471  edm::LogError(theCategory+"|Adjust")<<" cannot properly adjust for term: "<<i <<","<<j;
472  }
473  }
474  }}
475  LogDebug(theCategory+"|Adjust")<<"updated state sigma(i), rho(i,j): \n"<<simpleTerms;
476 
477  AlgebraicSymMatrix55 finalTerms;
478  complicatedTerm(simpleTerms, finalTerms);
479  LogDebug(theCategory+"|Adjust")<<"updated state COV(i,j): \n"<<finalTerms;
480 
481  CurvilinearTrajectoryError oMat(finalTerms);
482  state = FreeTrajectoryState(state.parameters(), oMat);
483  LogDebug(theCategory+"|Adjust")<<"updated state:\n"<<state;
484 }
485 
486 
487 
489 
490  AlgebraicSymMatrix55 simpleTerms;
491  simpleTerm(state.curvilinearError(), simpleTerms);
492  LogDebug(theCategory+"|Adjust")<<"state sigma(i), rho(i,j): \n"<<simpleTerms;
493 
494  AlgebraicSymMatrix55 simpleValues= get(state.globalMomentum(),false).matrix();
495  LogDebug(theCategory+"|Adjust")<<"config: (i,i), (i,j):\n"<<simpleValues;
496 
497  for(int i = 0;i!=5;i++){for(int j = i;j!=5;j++){
498  //check on each term for desired action
499  switch(theTermAction[Pindex(i,j)]){
500  case scale:
501  {
502  simpleTerms(i,j) *= simpleValues(i,j);
503  break;
504  }
505  case assign:
506  {
507  simpleTerms(i,j) = simpleValues(i,j);
508  break;
509  }
510  case error:
511  {
512  edm::LogError(theCategory+"|Adjust")<<" cannot properly adjust for term: "<<i <<","<<j;
513  }
514  }
515  }}
516  LogDebug(theCategory+"|Adjust")<<"updated state sigma(i), rho(i,j): \n"<<simpleTerms;
517 
518  AlgebraicSymMatrix55 finalTerms;
519  complicatedTerm(simpleTerms, finalTerms);
520  LogDebug(theCategory+"|Adjust")<<"updated state COV(i,j): \n"<<finalTerms;
521 
522  CurvilinearTrajectoryError oMat(finalTerms);
523  state = TrajectoryStateOnSurface(state.weight(),
524  state.globalParameters(),
525  oMat,
526  state.surface(),
527  state.surfaceSide()
528  );
529  LogDebug(theCategory+"|Adjust")<<"updated state:\n"<<state;
530 }
#define LogDebug(id)
const double Pi
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
void adjust(FreeTrajectoryState &state)
adjust the error matrix on the state
static bool divide(CurvilinearTrajectoryError &num_error, const CurvilinearTrajectoryError &denom_error)
divide term by term the two matrix
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 ...
const GlobalTrajectoryParameters & parameters() const
static const TString vars[5]
names of the variables of the 5x5 error matrix
const CurvilinearTrajectoryError & curvilinearError() const
double Value(GlobalVector &momentum, int i, int j, bool convolute=true)
internal method that retreives the value of the parametrization for term i,j
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
bool exists(std::string const &parameterName) const
checks if a parameter exists
CurvilinearTrajectoryError get(GlobalVector momentum, bool convolute=true)
main method to be used. Retrieve a 5x5 symetrical matrix according to parametrization of error or sca...
double maxEta
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
~MuonErrorMatrix()
destructor
T eta() const
static double Term(const AlgebraicSymMatrix55 &curv, int i, int j)
provide the numerical value used. sigma or correlation factor
int ii
Definition: cuy.py:588
const CurvilinearTrajectoryError & curvilinearError() const
action
enum type to define if the class is used as a tool or to be created
static std::string const input
Definition: EdmProvDump.cc:44
const SurfaceType & surface() const
void complicatedTerm(const AlgebraicSymMatrix55 &input, AlgebraicSymMatrix55 &output)
convert sigma/rho -&gt; sigma2/COV
TermAction theTermAction[15]
T sqrt(T t)
Definition: SSEVec.h:48
SurfaceSide surfaceSide() const
Position relative to material, defined relative to momentum vector.
tuple result
Definition: query.py:137
MuonErrorMatrix(const edm::ParameterSet &pset)
constructor from a parameter set
int j
Definition: DBlmapReader.cc:9
static void multiply(CurvilinearTrajectoryError &initial_error, const CurvilinearTrajectoryError &scale_error)
multiply term by term the two matrix
int Pindex(int i, int j)
internal methods to get the index of a matrix term.
GlobalVector momentum() const
#define LogTrace(id)
JetCorrectorParameters corr
Definition: classes.h:5
TProfile3D * theData_fast[5][5]
const GlobalTrajectoryParameters & globalParameters() const
void simpleTerm(const AlgebraicSymMatrix55 &input, AlgebraicSymMatrix55 &output)
convert sigma2/COV -&gt; sigma/rho
double Rms(GlobalVector &momentum, int i, int j)
internal method that retreives the error on the value of the parametrization for term i...
T eta() const
Definition: PV3DBase.h:76
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
const AlgebraicSymMatrix55 & matrix() const
GlobalVector globalMomentum() const
TDirectory * theD
the attached root file, where the parametrization is saved
double a
Definition: hdecay.h:121
CurvilinearTrajectoryError getFast(GlobalVector momentum)
std::string fullPath() const
Definition: FileInPath.cc:165
TProfile3D * theData[15]
15 TProfile, each holding he parametrization of each term of the 5x5
std::string theCategory
log category: &quot;MuonErrorMatrix&quot;
void close()
close the root file attached to the class
TProfile3D * Index(int i, int j)
internal method to get access to the profiles
Definition: DDAxes.h:10