00001 #include "RecoMuon/TrackingTools/interface/MuonErrorMatrix.h"
00002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00003 #include "FWCore/ParameterSet/interface/FileInPath.h"
00004
00005 #include "TROOT.h"
00006 #include "TString.h"
00007 #include "TRandom2.h"
00008 #include "TMath.h"
00009
00010 #include <sstream>
00011
00012 using namespace std;
00013
00014 const TString MuonErrorMatrix::vars[5]={"#frac{q}{|p|}","#lambda","#varphi_{0}","X_{T}","Y_{T}"};
00015
00016 MuonErrorMatrix::MuonErrorMatrix(const edm::ParameterSet & iConfig):theD(0){
00017 theCategory="MuonErrorMatrix";
00018 std::string action = iConfig.getParameter<std::string>("action");
00019
00020 bool madeFromCff=iConfig.exists("errorMatrixValuesPSet");
00021 edm::ParameterSet errorMatrixValuesPSet;
00022
00023 std::string fileName;
00024 if (!madeFromCff){ fileName = iConfig.getParameter<std::string>("rootFileName");}
00025 else { errorMatrixValuesPSet =iConfig.getParameter<edm::ParameterSet>("errorMatrixValuesPSet");}
00026 MuonErrorMatrix::action a= use;
00027
00028 int NPt=5;
00029 std::vector<double> xBins;
00030 double * xBinsArray = 0;
00031 double minPt=1;
00032 double maxPt=200;
00033 int NEta=5;
00034 std::vector<double> yBins;
00035 double * yBinsArray =0;
00036 double minEta=0;
00037 double maxEta=2.5;
00038 int NPhi=1;
00039 double minPhi=-TMath::Pi();
00040 double maxPhi=TMath::Pi();
00041
00042 if (action!="use"){
00043 a = constructor;
00044
00045 NPt = iConfig.getParameter<int>("NPt");
00046 if (NPt!=0){
00047 minPt = iConfig.getParameter<double>("minPt");
00048 maxPt = iConfig.getParameter<double>("maxPt");}
00049 else{
00050 xBins = iConfig.getParameter<std::vector<double> >("PtBins");
00051 if (xBins.size()==0){edm::LogError( theCategory)<<"Npt=0 and no entries in the vector. I will do aseg fault soon.";}
00052 NPt = xBins.size()-1;
00053 xBinsArray = &(xBins.front());
00054 minPt = xBins.front();
00055 maxPt = xBins.back();}
00056
00057 NEta = iConfig.getParameter<int>("NEta");
00058 if (NEta!=0){
00059 minEta = iConfig.getParameter<double>("minEta");
00060 maxEta = iConfig.getParameter<double>("maxEta");}
00061 else{
00062 yBins = iConfig.getParameter<std::vector<double> >("EtaBins");
00063 if (yBins.size()==0){edm::LogError( theCategory)<<"NEta=0 and no entries in the vector. I will do aseg fault soon.";}
00064 NEta = yBins.size()-1;
00065 yBinsArray = &(yBins.front());
00066 minPt = yBins.front();
00067 maxPt = yBins.back();}
00068
00069 NPhi = iConfig.getParameter<int>("NPhi");
00070 std::stringstream get(iConfig.getParameter<std::string>("minPhi"));
00071 if (get.str() =="-Pi")
00072 { minPhi =-TMath::Pi();}
00073 else if(get.str() =="Pi")
00074 { minPhi =TMath::Pi();}
00075 else { get>>minPhi;}
00076 get.str(iConfig.getParameter<std::string>("maxPhi"));
00077 if (get.str() =="-Pi")
00078 { maxPhi =-TMath::Pi();}
00079 else if(get.str() =="Pi")
00080 { maxPhi =TMath::Pi();}
00081 else { get>>maxPhi;}
00082 }
00083 else if (madeFromCff){
00084
00085 xBins = errorMatrixValuesPSet.getParameter<std::vector<double> >("xAxis");
00086 NPt = xBins.size()-1;
00087 xBinsArray = &(xBins.front());
00088 minPt = xBins.front();
00089 maxPt = xBins.back();
00090
00091 yBins = errorMatrixValuesPSet.getParameter<std::vector<double> >("yAxis");
00092 NEta = yBins.size()-1;
00093 yBinsArray = &(yBins.front());
00094 minPt = yBins.front();
00095 maxPt = yBins.back();
00096
00097 std::vector<double> zBins = errorMatrixValuesPSet.getParameter<std::vector<double> >("zAxis");
00098 NPhi=1;
00099 if (zBins.size()!=2){
00100 edm::LogError(theCategory)<<"please specify a zAxis with 2 entries only. more bins not implemented yet.";}
00101 minPhi=zBins.front();
00102 maxPhi=zBins.back();
00103
00104 }
00105
00106 if (a==use){
00107 if (!madeFromCff){
00108 edm::LogInfo(theCategory)<<"using an error matrix object from: "<<fileName;
00109 edm::FileInPath data(fileName);
00110 std::string fullpath = data.fullPath();
00111 gROOT->cd();
00112 theD = new TFile(fullpath.c_str());
00113 theD->SetWritable(false);
00114 }else{
00115 static unsigned int neverTheSame=0;
00116 std::stringstream dirName("MuonErrorMatrixDirectory");
00117 dirName<<neverTheSame++;
00118 edm::LogInfo(theCategory)<<"using an error matrix object from configuration file. putting memory histograms to: "<<dirName.str();
00119 gROOT->cd();
00120 theD = new TDirectory(dirName.str().c_str(),"transient directory to host MuonErrorMatrix TProfile3D");
00121 theD->SetWritable(false);
00122 }
00123 }
00124 else{
00125 edm::LogInfo(theCategory)<<"creating an error matrix object: "<<fileName;
00126 theD = new TFile(fileName.c_str(),"recreate");
00127 }
00128
00129 if (a==use && !madeFromCff ){gROOT->cd();}
00130 else {theD->cd();}
00131
00132 for (int i=0;i!=5;i++){for (int j=i;j!=5;j++){
00133 TString pfname(Form("pf3_V%1d%1d",i+1,j+1));
00134 TProfile3D * pf =0;
00135 if (a==use && !madeFromCff ){
00136
00137 edm::LogVerbatim(theCategory)<<"getting "<<pfname<<" from "<<fileName;
00138 pf = (TProfile3D *)theD->Get(pfname);
00139 theData[Pindex(i,j)]=pf;
00140 theData_fast[i][j]=pf; theData_fast[j][i]=pf;
00141 }
00142 else{
00143
00144
00145
00146 TString pftitle;
00147 if (i==j){pftitle="#sigma_{"+vars[i]+"}";}
00148 else{pftitle="#rho("+vars[i]+","+vars[j]+")";}
00149 edm::LogVerbatim(theCategory)<<"booking "<<pfname<<" into "<<fileName;
00150 pf = new TProfile3D(pfname,pftitle,NPt,minPt,maxPt,NEta,minEta,maxEta,NPhi,minPhi,maxPhi,"S");
00151 pf->SetXTitle("muon p_{T} [GeV]");
00152 pf->SetYTitle("muon |#eta|");
00153 pf->SetZTitle("muon #varphi");
00154
00155
00156 if (xBinsArray){
00157 pf->GetXaxis()->Set(NPt,xBinsArray);}
00158 if (yBinsArray){
00159 pf->GetYaxis()->Set(NEta,yBinsArray);}
00160
00161 if (madeFromCff){
00162 edm::ParameterSet pSet = errorMatrixValuesPSet.getParameter<edm::ParameterSet>(pfname.Data());
00163
00164 std::vector<double> values = pSet.getParameter<std::vector<double> >("values");
00165 unsigned int iX=pf->GetNbinsX();
00166 unsigned int iY=pf->GetNbinsY();
00167 unsigned int iZ=pf->GetNbinsZ();
00168 unsigned int continuous_i=0;
00169 for(unsigned int ix=1;ix<=iX;++ix){
00170 for(unsigned int iy=1;iy<=iY;++iy){
00171 for(unsigned int iz=1;iz<=iZ;++iz){
00172 LogTrace(theCategory)<<"filling profile:"
00173 <<"\n pt (x)= "<<pf->GetXaxis()->GetBinCenter(ix)
00174 <<"\n eta (y)= "<<pf->GetYaxis()->GetBinCenter(iy)
00175 <<"\n phi (z)= "<<pf->GetZaxis()->GetBinCenter(iz)
00176 <<"\n value= "<<values[continuous_i];
00177 pf->Fill(pf->GetXaxis()->GetBinCenter(ix),
00178 pf->GetYaxis()->GetBinCenter(iy),
00179 pf->GetZaxis()->GetBinCenter(iz),
00180 values[continuous_i++]);
00181 }}}
00182
00183 std::string tAction = pSet.getParameter<std::string>("action");
00184 if (tAction == "scale")
00185 theTermAction[Pindex(i,j)] = scale;
00186 else if (tAction == "assign")
00187 theTermAction[Pindex(i,j)] = assign;
00188 else{
00189 edm::LogError(theCategory)<<" wrong configuration: term action: "<<tAction<<" is not recognized.";
00190 theTermAction[Pindex(i,j)] = error;
00191 }
00192
00193 }
00194 }
00195
00196 LogDebug(theCategory)<<" index "<<i<<":"<<j<<" -> "<<Pindex(i,j);
00197 theData[Pindex(i,j)]=pf;
00198 theData_fast[i][j]=pf; theData_fast[j][i]=pf;
00199 if (!pf){
00200 edm::LogError(theCategory)<<" profile "<<pfname<<" in file "<<fileName<<" is not valid. exiting.";
00201 exit(1);
00202 }
00203 }}
00204
00205
00206 for (int i=0;i!=15;i++){
00207 if (theData[i]) {edm::LogVerbatim(theCategory)<<i<<" :"<<theData[i]->GetName()
00208 <<" "<< theData[i]->GetTitle()<<std::endl;}}
00209
00210 }
00211
00212
00213
00214 void MuonErrorMatrix::close(){
00215
00216 if (theD){
00217 theD->cd();
00218
00219 if (theD->IsWritable()){
00220 for (int i=0;i!=15;i++){ if (theData[i]) { theData[i]->Write();}}}
00221 theD->Close();
00222 }}
00223
00224 MuonErrorMatrix::~MuonErrorMatrix() {
00225 close();
00226 }
00227
00228
00229
00230 CurvilinearTrajectoryError MuonErrorMatrix::get(GlobalVector momentum,bool convolute) {
00231 AlgebraicSymMatrix55 V;
00232
00233 for (int i=0;i!=5;i++){for (int j=i;j!=5;j++){
00234 V(i,j) = Value(momentum,i,j,convolute);}}
00235 return CurvilinearTrajectoryError(V);}
00236
00237 CurvilinearTrajectoryError MuonErrorMatrix::getFast(GlobalVector momentum) {
00238
00239
00240 AlgebraicSymMatrix55 V;
00241
00242 double pT = momentum.perp();
00243 double eta = fabs(momentum.eta());
00244 double phi = momentum.phi();
00245
00246
00247 int iBin_x= findBin(theData_fast[0][0]->GetXaxis(),pT);
00248 int iBin_y= findBin(theData_fast[0][0]->GetYaxis(),eta);
00249 int iBin_z= findBin(theData_fast[0][0]->GetZaxis(),phi);
00250
00251
00252 double values[5][5];
00253 for (int i=0;i!=5;i++){for (int j=i;j!=5;j++){
00254 values[i][j]=theData_fast[i][j]->GetBinContent(iBin_x, iBin_y, iBin_z);
00255 }}
00256
00257 for (int i=0;i!=5;i++){for (int j=i;j!=5;j++){
00258 if (i==j){
00259
00260 V(i,j) = values[i][j];
00261 V(i,j)*=V(i,j);
00262 }
00263 else{
00264
00265 V(i,j) = values[i][i] * values[j][j] * values[i][j];
00266 }
00267 }}
00268
00269 return CurvilinearTrajectoryError(V);}
00270
00271
00272
00273
00274
00275
00276
00277
00278
00279
00280
00281
00282
00283
00284
00285
00286
00287
00288 int MuonErrorMatrix::findBin(TAxis * axis, double value){
00289
00290 int result = axis->FindBin(value);
00291 if (result <= 0) result=1;
00292 else if (result > axis->GetNbins() ) result = axis->GetNbins();
00293 return result;}
00294
00295
00296 double MuonErrorMatrix::Value(GlobalVector & momentum, int i, int j,bool convolute) {
00297 double result=0;
00298 TProfile3D * ij = Index(i,j);
00299 if (!ij) {edm::LogError(theCategory)<<"cannot get the profile ("<<i<<":"<<j<<")"; return result;}
00300
00301 double pT = momentum.perp();
00302 double eta = fabs(momentum.eta());
00303 double phi = momentum.phi();
00304
00305 int iBin_x= findBin(ij->GetXaxis(),pT);
00306 int iBin_y= findBin(ij->GetYaxis(),eta);
00307 int iBin_z= findBin(ij->GetZaxis(),phi);
00308
00309 if (convolute){
00310 if (i!=j){
00311
00312 TProfile3D * ii = Index(i,i);
00313 TProfile3D * jj = Index(j,j);
00314 if (!ii){edm::LogError(theCategory)<<"cannot get the profile ("<<i<<":"<<i<<")"; return result;}
00315 if (!jj){edm::LogError(theCategory)<<"cannot get the profile ("<<j<<":"<<j<<")"; return result;}
00316
00317
00318 int iBin_i_x = findBin(ii->GetXaxis(),pT);
00319 int iBin_i_y = findBin(ii->GetYaxis(),eta);
00320 int iBin_i_z = findBin(ii->GetZaxis(),phi);
00321 int iBin_j_x = findBin(jj->GetXaxis(),pT);
00322 int iBin_j_y = findBin(jj->GetYaxis(),eta);
00323 int iBin_j_z = findBin(jj->GetZaxis(),phi);
00324
00325 double corr = ij->GetBinContent(iBin_x,iBin_y,iBin_z);
00326 double sigma_1 = (ii->GetBinContent(iBin_i_x,iBin_i_y,iBin_i_z));
00327 double sigma_2 = (jj->GetBinContent(iBin_j_x,iBin_j_y,iBin_j_z));
00328
00329
00330 result=corr*sigma_1*sigma_2;
00331 LogDebug(theCategory)<<"for: (pT,eta,phi)=("<<pT<<", "<<eta<<", "<<phi<<") nterms are"
00332 <<"\nrho["<<i<<","<<j<<"]: "<<corr<<" ["<< iBin_x<<", "<<iBin_y<<", "<<iBin_z<<"]"
00333 <<"\nsigma["<<i<<","<<i<<"]: "<< sigma_1
00334 <<"\nsigma["<<j<<","<<j<<"]: "<< sigma_2
00335 <<"Covariance["<<i<<","<<j<<"] is: "<<result;
00336 return result;
00337 }
00338 else{
00339
00340
00341 result=ij->GetBinContent(iBin_x,iBin_y,iBin_z);
00342 result*=result;
00343 LogDebug(theCategory)<<"for: (pT,eta,phi)=("<<pT<<", "<<eta<<", "<<phi<<") sigma^2["<<i<<","<<j<<"] is: "<<result;
00344 return result;
00345 }
00346 }else{
00347
00348 result=ij->GetBinContent(iBin_x,iBin_y,iBin_z);
00349 return result;
00350 }
00351 }
00352
00353 double MuonErrorMatrix::Rms(GlobalVector & momentum, int i, int j) {
00354 double result=0;
00355 TProfile3D * ij = Index(i,j);
00356 if (!ij){edm::LogError(theCategory)<<"cannot get the profile ("<<i<<":"<<i<<")"; return result;}
00357 double pT = momentum.perp();
00358 double eta = fabs(momentum.eta());
00359 double phi = momentum.phi();
00360
00361 int iBin_x= ij->GetXaxis()->FindBin(pT);
00362 int iBin_y= ij->GetYaxis()->FindBin(eta);
00363 int iBin_z= ij->GetZaxis()->FindBin(phi);
00364 result=ij->GetBinError(iBin_x,iBin_y,iBin_z);
00365
00366 LogDebug(theCategory)<<"for: (pT,eta,phi)=("<<pT<<", "<<eta<<", "<<phi<<") error["<<i<<","<<j<<"] is: "<<result;
00367 return result;
00368 }
00369
00370 double MuonErrorMatrix::Term(const AlgebraicSymMatrix55 & curv, int i, int j){
00371
00372 double result=0;
00373 if (i==j){
00374 result = curv(i,j);
00375 if (result<0){
00376
00377 edm::LogError("MuonErrorMatrix")<<"invalid term in the error matrix.\n sii: "<< result;
00378 return 0;}
00379 return sqrt(result);}
00380 else{
00381 double si = curv(i,i);
00382 double sj = curv(j,j);
00383
00384 if (si<=0 || sj<=0){
00385 edm::LogError("MuonErrorMatrix")<<"invalid term in the error matrix.\n si: "
00386 <<si<<" sj: "<<sj<<". result will be corrupted\n"<<curv;
00387 return 0;}
00388
00389 si=sqrt(si);
00390 sj=sqrt(sj);
00391
00392
00393 return result = curv(i,j)/(si*sj);}
00394
00395 return 0;
00396 }
00397
00398
00399 void MuonErrorMatrix::multiply(CurvilinearTrajectoryError & initial_error, const CurvilinearTrajectoryError & scale_error){
00400
00401 const AlgebraicSymMatrix55 & scale_matrix=scale_error.matrix();
00402 AlgebraicSymMatrix55 revised_matrix = initial_error.matrix();
00403
00404
00405 for(int i = 0;i!=5;i++){for(int j = i;j!=5;j++){
00406 revised_matrix(i,j)*=scale_matrix(i,j);
00407 }}
00408 initial_error = CurvilinearTrajectoryError(revised_matrix);
00409 }
00410 bool MuonErrorMatrix::divide(CurvilinearTrajectoryError & num_error, const CurvilinearTrajectoryError & denom_error){
00411
00412 const AlgebraicSymMatrix55 & denom_matrix=denom_error.matrix();
00413 AlgebraicSymMatrix55 num_matrix = num_error.matrix();
00414
00415
00416 for(int i = 0;i!=5;i++){for(int j = i;j!=5;j++){
00417 if (denom_matrix(i,j)==0) return false;
00418 num_matrix(i,j)/=denom_matrix(i,j);
00419 }}
00420 num_error = CurvilinearTrajectoryError(num_matrix);
00421 return true;
00422 }
00423
00424 void MuonErrorMatrix::simpleTerm(const AlgebraicSymMatrix55 & input, AlgebraicSymMatrix55 & output){
00425 for(int i = 0;i!=5;i++){for(int j = i;j!=5;j++){
00426 if (i==j)
00427 output(i,j) = sqrt(input(i,j));
00428 else
00429 output(i,j) = input(i,j) / sqrt( input(i,i)* input(j,j));
00430 }}
00431 }
00432
00433 void MuonErrorMatrix::complicatedTerm(const AlgebraicSymMatrix55 & input, AlgebraicSymMatrix55 & output){
00434
00435 for(int i = 0;i!=5;i++){for(int j = i;j!=5;j++){
00436 if (i==j)
00437 output(i,j) = input(i,j)*input(i,j);
00438 else
00439 output(i,j) = input(i,j)*input(i,i)*input(j,j);
00440 }}
00441 }
00442
00443
00444
00445 void MuonErrorMatrix::adjust(FreeTrajectoryState & state){
00446
00447 LogDebug(theCategory+"|Adjust")<<"state: \n"<<state;
00448 AlgebraicSymMatrix55 simpleTerms;
00449 simpleTerm(state.curvilinearError(), simpleTerms);
00450
00451 LogDebug(theCategory+"|Adjust")<<"state sigma(i), rho(i,j): \n"<<simpleTerms;
00452
00453 AlgebraicSymMatrix55 simpleValues=get(state.momentum(),false).matrix();
00454 LogDebug(theCategory+"|Adjust")<<"config: (i,i), (i,j): \n"<<simpleValues;
00455
00456 for(int i = 0;i!=5;i++){for(int j = i;j!=5;j++){
00457
00458 switch(theTermAction[Pindex(i,j)]){
00459 case scale:
00460 {
00461 simpleTerms(i,j) *= simpleValues(i,j);
00462 break;
00463 }
00464 case assign:
00465 {
00466 simpleTerms(i,j) = simpleValues(i,j);
00467 break;
00468 }
00469 case error:
00470 {
00471 edm::LogError(theCategory+"|Adjust")<<" cannot properly adjust for term: "<<i <<","<<j;
00472 }
00473 }
00474 }}
00475 LogDebug(theCategory+"|Adjust")<<"updated state sigma(i), rho(i,j): \n"<<simpleTerms;
00476
00477 AlgebraicSymMatrix55 finalTerms;
00478 complicatedTerm(simpleTerms, finalTerms);
00479 LogDebug(theCategory+"|Adjust")<<"updated state COV(i,j): \n"<<finalTerms;
00480
00481 CurvilinearTrajectoryError oMat(finalTerms);
00482 state = FreeTrajectoryState(state.parameters(), oMat);
00483 LogDebug(theCategory+"|Adjust")<<"updated state:\n"<<state;
00484 }
00485
00486
00487
00488 void MuonErrorMatrix::adjust(TrajectoryStateOnSurface & state){
00489
00490 AlgebraicSymMatrix55 simpleTerms;
00491 simpleTerm(state.curvilinearError(), simpleTerms);
00492 LogDebug(theCategory+"|Adjust")<<"state sigma(i), rho(i,j): \n"<<simpleTerms;
00493
00494 AlgebraicSymMatrix55 simpleValues= get(state.globalMomentum(),false).matrix();
00495 LogDebug(theCategory+"|Adjust")<<"config: (i,i), (i,j):\n"<<simpleValues;
00496
00497 for(int i = 0;i!=5;i++){for(int j = i;j!=5;j++){
00498
00499 switch(theTermAction[Pindex(i,j)]){
00500 case scale:
00501 {
00502 simpleTerms(i,j) *= simpleValues(i,j);
00503 break;
00504 }
00505 case assign:
00506 {
00507 simpleTerms(i,j) = simpleValues(i,j);
00508 break;
00509 }
00510 case error:
00511 {
00512 edm::LogError(theCategory+"|Adjust")<<" cannot properly adjust for term: "<<i <<","<<j;
00513 }
00514 }
00515 }}
00516 LogDebug(theCategory+"|Adjust")<<"updated state sigma(i), rho(i,j): \n"<<simpleTerms;
00517
00518 AlgebraicSymMatrix55 finalTerms;
00519 complicatedTerm(simpleTerms, finalTerms);
00520 LogDebug(theCategory+"|Adjust")<<"updated state COV(i,j): \n"<<finalTerms;
00521
00522 CurvilinearTrajectoryError oMat(finalTerms);
00523 state = TrajectoryStateOnSurface(state.globalParameters(),
00524 oMat,
00525 state.surface(),
00526 state.surfaceSide(),
00527 state.weight());
00528 LogDebug(theCategory+"|Adjust")<<"updated state:\n"<<state;
00529 }