Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-03 08:20:17

0001 /*
0002  * This file is part of KFParticle package
0003  * Copyright (C) 2007-2019 FIAS Frankfurt Institute for Advanced Studies
0004  *               2007-2019 Goethe University of Frankfurt
0005  *               2007-2019 Ivan Kisel <I.Kisel@compeng.uni-frankfurt.de>
0006  *               2007-2019 Maksym Zyzak
0007  *
0008  * KFParticle is free software: you can redistribute it and/or modify
0009  * it under the terms of the GNU General Public License as published by
0010  * the Free Software Foundation, either version 3 of the License, or
0011  * (at your option) any later version.
0012  *
0013  * KFParticle is distributed in the hope that it will be useful,
0014  * but WITHOUT ANY WARRANTY; without even the implied warranty of
0015  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
0016  * GNU General Public License for more details.
0017  *
0018  * You should have received a copy of the GNU General Public License
0019  * along with this program.  If not, see <https://www.gnu.org/licenses/>.
0020  */
0021 
0022 #ifdef DO_TPCCATRACKER_EFF_PERFORMANCE
0023 
0024 #include "KFParticlePerformanceBase.h"
0025 
0026 #include "TDirectory.h"
0027 #include "TH1.h"
0028 #include "TH2.h"
0029 #include "TH3.h"
0030 #include "TProfile.h"
0031 #include "TProfile2D.h"
0032 
0033 
0034 KFParticlePerformanceBase::KFParticlePerformanceBase():
0035   fParteff(), fPVeff(), fPVeffMCReconstructable(), outfileName(), histodir(0), fNEvents(0), fStoreMCHistograms(1), 
0036   fStorePrimSecHistograms(1), fStoreZRHistograms(1),fHistoDir(0)
0037 {
0038   /** The default constructor. Initialises all pointers to nullptr.  **/
0039   for(int iParticle=0; iParticle<KFPartEfficiencies::nParticles; iParticle++)
0040   {
0041     for(int iFitQA=0; iFitQA<nFitQA; iFitQA++)
0042     {
0043       hFitDaughtersQA         [iParticle][iFitQA] = nullptr;
0044       hFitQA                  [iParticle][iFitQA] = nullptr;
0045       hFitQANoConstraint      [iParticle][iFitQA] = nullptr;
0046       hFitQAMassConstraint    [iParticle][iFitQA] = nullptr;
0047       hFitQATopoConstraint    [iParticle][iFitQA] = nullptr;
0048       hFitQATopoMassConstraint[iParticle][iFitQA] = nullptr;
0049     }
0050   }
0051 
0052   for(int iParticle=0; iParticle<KFPartEfficiencies::nParticles; iParticle++)
0053     for(int iQA=0; iQA<nDSToParticleQA; iQA++)
0054       hDSToParticleQA[iParticle][iQA] = nullptr;
0055   
0056   for(int iParameterSet=0; iParameterSet<nParametersSet; iParameterSet++)
0057   {
0058     for(int iParticle=0; iParticle<KFPartEfficiencies::nParticles; iParticle++)
0059     {
0060       for(int iHisto=0; iHisto<nHistoPartParam; iHisto++)
0061       {
0062         hPartParam               [iParameterSet][iParticle][iHisto] = nullptr;
0063         hPartParamPrimary        [iParameterSet][iParticle][iHisto] = nullptr;
0064         hPartParamPrimaryMass    [iParameterSet][iParticle][iHisto] = nullptr;
0065         hPartParamPrimaryTopo    [iParameterSet][iParticle][iHisto] = nullptr;
0066         hPartParamPrimaryTopoMass[iParameterSet][iParticle][iHisto] = nullptr;
0067         hPartParamSecondary      [iParameterSet][iParticle][iHisto] = nullptr;
0068         hPartParamSecondaryMass  [iParameterSet][iParticle][iHisto] = nullptr;
0069       }
0070     }
0071   }
0072 
0073   for(int iParameterSet=0; iParameterSet<nParametersSet; iParameterSet++)
0074   {
0075     for(int iParticle=0; iParticle<KFPartEfficiencies::nParticles; iParticle++)
0076     {
0077       for(int iHisto=0; iHisto<nHistoPartParam2D; iHisto++)
0078       {
0079         hPartParam2D               [iParameterSet][iParticle][iHisto] = nullptr;
0080         hPartParam2DPrimary        [iParameterSet][iParticle][iHisto] = nullptr;
0081         hPartParam2DPrimaryMass    [iParameterSet][iParticle][iHisto] = nullptr;
0082         hPartParam2DPrimaryTopo    [iParameterSet][iParticle][iHisto] = nullptr;
0083         hPartParam2DPrimaryTopoMass[iParameterSet][iParticle][iHisto] = nullptr;
0084         hPartParam2DSecondary      [iParameterSet][iParticle][iHisto] = nullptr;
0085         hPartParam2DSecondaryMass  [iParameterSet][iParticle][iHisto] = nullptr;
0086       }
0087     }
0088   }
0089 
0090   for(int iParticle=0; iParticle<KFPartEfficiencies::nParticles; iParticle++)
0091     for(int iHisto=0; iHisto<nHistoPartParam3D; iHisto++)
0092       hPartParam3D[0][iParticle][iHisto] = nullptr;
0093 
0094   for(int iParticle=0; iParticle<KFPartEfficiencies::nParticles; iParticle++)
0095     for(int iEffSet=0; iEffSet<3; iEffSet++)
0096       for(int iEff=0; iEff<nPartEfficiency; iEff++)
0097         hPartEfficiency[iParticle][iEffSet][iEff] = nullptr;
0098         
0099   for(int iParticle=0; iParticle<KFPartEfficiencies::nParticles; iParticle++)
0100     for(int iEffSet=0; iEffSet<3; iEffSet++)
0101       for(int iEff=0; iEff<nPartEfficiency2D; iEff++)
0102         hPartEfficiency2D[iParticle][iEffSet][iEff] = nullptr;
0103   
0104   for(int iEffSet=0; iEffSet<2; iEffSet++)
0105     for(int iHistoPV=0; iHistoPV<nHistosPV; iHistoPV++)
0106       hPVFitQa[iEffSet][iHistoPV] = nullptr;
0107 
0108   for(int iEffSet1=0; iEffSet1<2; iEffSet1++)
0109     for(int iEffSet2=0; iEffSet2<2; iEffSet2++)
0110       for(int iHistoPV=0; iHistoPV<nHistosPV-1; iHistoPV++)
0111         hPVFitQa2D[iEffSet1][iEffSet2][iHistoPV] = nullptr;
0112 
0113   for(int iParam=0; iParam<nHistosPVParam; iParam++)
0114   {
0115     hPVParam      [iParam] = nullptr;
0116     hPVParamGhost [iParam] = nullptr;
0117     hPVParamSignal[iParam] = nullptr;
0118     hPVParamPileup[iParam] = nullptr;
0119     hPVParamBG    [iParam] = nullptr;
0120   }
0121   
0122   for(int iParam=0; iParam<nHistosPVParam2D; iParam++)
0123     hPVParam2D[iParam] = nullptr;
0124   
0125   for(int iFitQA=0; iFitQA<nFitPVTracksQA; iFitQA++)
0126     hFitPVTracksQA[iFitQA] = nullptr;
0127   
0128   for(int iTP=0; iTP<nHistosTP; iTP++)
0129     hTrackParameters[iTP] = nullptr;
0130 
0131   for(int iEffSet=0; iEffSet<4; iEffSet++)
0132     for(int iEff=0; iEff<nPVefficiency; iEff++)
0133       hPVefficiency[iEffSet][iEff] = nullptr;
0134 }
0135 
0136 void KFParticlePerformanceBase::CreateHistos(std::string histoDir, TDirectory* outFile, std::map<int,bool> decays)
0137 {
0138   /** Creates all histograms. If "outFile" is provided - creates a new ROOT directory and stores all 
0139    ** histograms there. Otherwise histograms are stored in TDirectory::CurrentDirectory().
0140    ** \param[in] histoDir - name of the ROOT directory with histograms
0141    ** \param[in] outFile - pointer to the external ROOT directory or file, where histograms should be stored
0142    ** \param[in] decays - a list of decays, for which histograms are created, if empty - histograms are
0143    ** created for all decay channels from the KF Particle Finder reconstruction scheme
0144    **/
0145   TDirectory *curdir = gDirectory;
0146   if (outFile) {
0147     outFile->cd();
0148     fHistoDir = outFile;
0149     if (histoDir != "") {
0150       fHistoDir = outFile->mkdir( TString(histoDir) );
0151       fHistoDir->cd();
0152     }
0153   } else {
0154     fHistoDir = TDirectory::CurrentDirectory();
0155   }
0156   {
0157     gDirectory->mkdir("KFParticlesFinder");
0158     gDirectory->cd("KFParticlesFinder");
0159     histodir = gDirectory;
0160     gDirectory->mkdir("Particles");
0161     gDirectory->cd("Particles");
0162     for(int iPart=0; iPart<fParteff.nParticles; ++iPart)
0163     {
0164       if(!(decays.empty()) && (iPart < fParteff.fFirstStableParticleIndex || iPart > fParteff.fLastStableParticleIndex))
0165         if(decays.find(fParteff.partPDG[iPart]) == decays.end()) continue;
0166         
0167       gDirectory->mkdir(fParteff.partName[iPart].data());
0168       gDirectory->cd(fParteff.partName[iPart].data());
0169       {
0170         if(fStoreMCHistograms)
0171         {
0172           TString res = "res";
0173           TString pull = "pull";
0174 
0175           gDirectory->mkdir("DaughtersQA");
0176           gDirectory->cd("DaughtersQA");
0177           {
0178             TString parName[nFitQA/2] = {"X","Y","Z","Px","Py","Pz","E","M"};
0179             int nBins = 100;
0180             float xMax[nFitQA/2] = {0.15,0.15,0.03,0.01,0.01,0.06,0.06,0.01};
0181   //             float xMax[nFitQA/2] = {2.,2.,5.,0.3,0.3,0.3,0.03,0.03};
0182 
0183             for( int iH=0; iH<nFitQA/2; iH++ ){
0184               hFitDaughtersQA[iPart][iH]   = new TH1F((res+parName[iH]).Data(),
0185                                                       (GetDirectoryPath()+res+parName[iH]).Data(), 
0186                                                       nBins, -xMax[iH],xMax[iH]);
0187               hFitDaughtersQA[iPart][iH+8] = new TH1F((pull+parName[iH]).Data(),
0188                                                       (GetDirectoryPath()+pull+parName[iH]).Data(), 
0189                                                       nBins, -6,6);
0190             }
0191           }
0192           gDirectory->cd(".."); //particle directory
0193 
0194           gDirectory->mkdir("DSToParticleQA");
0195           gDirectory->cd("DSToParticleQA");
0196           {
0197             TString parName[3] = {"X","Y","Z"};
0198             int nBins = 100;
0199             float xMax[3] = {0.5, 0.5, 2.};
0200 
0201             for( int iH=0; iH<3; iH++ ){
0202               hDSToParticleQA[iPart][iH]   = new TH1F((res+parName[iH]).Data(),
0203                                                       (GetDirectoryPath()+res+parName[iH]).Data(), 
0204                                                       nBins, -xMax[iH],xMax[iH]);
0205               hDSToParticleQA[iPart][iH+3] = new TH1F((pull+parName[iH]).Data(),
0206                                                       (GetDirectoryPath()+pull+parName[iH]).Data(), 
0207                                                       nBins, -6,6);
0208             }
0209             
0210             hDSToParticleQA[iPart][6] = new TH1F("r", (GetDirectoryPath()+TString("r")).Data(), 1000, 0.0, 20.0);
0211           }
0212           gDirectory->cd(".."); //particle directory
0213           
0214           CreateFitHistograms(hFitQA[iPart], iPart);
0215           CreateEfficiencyHistograms(hPartEfficiency[iPart],hPartEfficiency2D[iPart]);
0216         }
0217         gDirectory->mkdir("Parameters");
0218         gDirectory->cd("Parameters");
0219         {
0220           const bool drawZR = IsCollectZRHistogram(iPart);
0221           CreateParameterHistograms(hPartParam[0], hPartParam2D[0], hPartParam3D[0], iPart, drawZR);
0222 
0223           if(IsCollect3DHistogram(iPart))
0224           {
0225             gDirectory->mkdir("SignalReco");
0226             gDirectory->cd("SignalReco");
0227             {
0228               CreateParameterHistograms(hPartParam[4], hPartParam2D[4], 0, iPart, drawZR);
0229             }
0230             gDirectory->cd(".."); // Parameters
0231             gDirectory->mkdir("BGReco");
0232             gDirectory->cd("BGReco");
0233             {
0234               CreateParameterHistograms(hPartParam[5], hPartParam2D[5], 0, iPart, drawZR);
0235             }
0236             gDirectory->cd(".."); // Parameters
0237           }
0238           
0239           if(fStoreMCHistograms)
0240           {
0241             gDirectory->mkdir("Signal");
0242             gDirectory->cd("Signal");
0243             {
0244               CreateParameterHistograms(hPartParam[1], hPartParam2D[1], 0, iPart, drawZR);
0245             }
0246             gDirectory->cd(".."); // particle directory / Parameters
0247             gDirectory->mkdir("Background");
0248             gDirectory->cd("Background");
0249             {
0250               CreateParameterHistograms(hPartParam[2], hPartParam2D[2], 0, iPart, drawZR);
0251             }
0252             gDirectory->cd(".."); // particle directory
0253             gDirectory->mkdir("Ghost");
0254             gDirectory->cd("Ghost");
0255             {
0256               CreateParameterHistograms(hPartParam[3], hPartParam2D[3], 0, iPart, drawZR);
0257             }
0258             gDirectory->cd(".."); // Parameters
0259             gDirectory->mkdir("MCSignal");
0260             gDirectory->cd("MCSignal");
0261             {
0262               CreateParameterHistograms(hPartParam[6], hPartParam2D[6], 0, iPart, drawZR);
0263             }
0264             gDirectory->cd(".."); // Parameters
0265             
0266             
0267             bool plotPrimaryHistograms = abs(fParteff.partPDG[iPart]) == 310 ||
0268                                          abs(fParteff.partPDG[iPart]) == 3122 ||
0269                                          abs(fParteff.partPDG[iPart]) == 22 ||
0270                                          abs(fParteff.partPDG[iPart]) == 111 ||
0271                                          abs(fParteff.partPDG[iPart]) == 3312 ||
0272                                          abs(fParteff.partPDG[iPart]) == 3334;  
0273                                           
0274             bool plotSecondaryHistograms = abs(fParteff.partPDG[iPart]) == 310 ||
0275                                            abs(fParteff.partPDG[iPart]) == 3122 ||
0276                                            abs(fParteff.partPDG[iPart]) == 22 ||
0277                                            abs(fParteff.partPDG[iPart]) == 111;
0278                                             
0279             if(fStorePrimSecHistograms && plotPrimaryHistograms)
0280             {
0281               gDirectory->mkdir("Primary");
0282               gDirectory->cd("Primary");
0283               {
0284                 CreateParameterSubfolder("NoConstraint (1C-Fit)", hPartParamPrimary, hPartParam2DPrimary, hFitQANoConstraint, iPart, true);
0285                 CreateParameterSubfolder("MassConstraint (2C-Fit)", hPartParamPrimaryMass, hPartParam2DPrimaryMass, hFitQAMassConstraint, iPart, true);
0286                 CreateParameterSubfolder("PVConstraint (3C-Fit)", hPartParamPrimaryTopo, hPartParam2DPrimaryTopo, hFitQATopoConstraint, iPart, true);
0287                 CreateParameterSubfolder("PVMassConstraint (4C-Fit)", hPartParamPrimaryTopoMass, hPartParam2DPrimaryTopoMass, hFitQATopoMassConstraint, iPart, true);
0288               }
0289               gDirectory->cd(".."); // particle directory / Parameters
0290             }
0291             
0292             if(fStorePrimSecHistograms && plotSecondaryHistograms)
0293             {
0294               gDirectory->mkdir("Secondary");
0295               gDirectory->cd("Secondary");
0296               {
0297                 CreateParameterSubfolder("NoConstraint (1C-Fit)", hPartParamSecondary, hPartParam2DSecondary, 0, iPart, true);
0298                 CreateParameterSubfolder("MassConstraint (2C-Fit)", hPartParamSecondaryMass, hPartParam2DSecondaryMass, 0, iPart, true);
0299               }
0300               gDirectory->cd(".."); // particle directory / Parameters
0301             }
0302           }
0303         }
0304         gDirectory->cd(".."); //particle directory
0305       }
0306       gDirectory->cd(".."); //Particles
0307     }
0308     gDirectory->cd(".."); //main
0309     gDirectory->mkdir("PrimaryVertexQA");
0310     gDirectory->cd("PrimaryVertexQA");
0311     {
0312       struct {
0313         TString name;
0314         TString title;
0315         Int_t n;
0316         Double_t l,r;
0317       } Table[nHistosPV]=
0318       {
0319         {"PVResX",  "x_{rec}-x_{mc}, cm", 100, -0.1f, 0.1f},
0320         {"PVResY",  "y_{rec}-y_{mc}, cm", 100, -0.1f, 0.1f},
0321         {"PVResZ",  "z_{rec}-z_{mc}, cm", 100, -1.f, 1.f},
0322         {"PVPullX", "Pull X",             100, -6.f, 6.f},
0323         {"PVPullY", "Pull Y",             100, -6.f, 6.f},
0324         {"PVPullZ", "Pull Z",             100, -6.f, 6.f},
0325         {"Lost",    "Lost tracks",        102, -0.01f, 1.01f}        
0326       };
0327       
0328       TString parName[nHistosPVParam] = {"x","y","z","r","Ntracks","Chi2","NDF","Chi2NDF","prob", "PVpurity", 
0329                                          "ghostTr", "triggerTr", "pileupTr", "bgTr", "dzSamePV"};
0330       TString parAxisName[nHistosPVParam] = {"x [cm]","y [cm]","z [cm]","r [cm]","N tracks","Chi2","NDF","Chi2NDF","prob","purity",
0331                                              "ghost tracks [%]", "trigger tracks [%]", "pileup tracks [%]", "bg tracks [%]", "dz [cm]"};
0332       int nBins[nHistosPVParam] = {1000,1000,1000,1000,1001,10000,1001,10000,100,102,102,102,102,102,1000};
0333       float xMin[nHistosPVParam] = {-1., -1., -10.,  0,   -0.5,    0.,   -0.5,    0., 0., -0.01, -0.01, -0.01, -0.01, -0.01, 0.};
0334       float xMax[nHistosPVParam] = { 1.,  1.,  10., 10, 1000.5, 1000., 1000.5, 1000., 1.,  1.01,  1.01,  1.01,  1.01,  1.01, 100.};
0335       
0336       TString parName2D[nHistosPVParam2D] = {"xy"};
0337       TString parXAxisName2D[nHistosPVParam2D] = {"x [cm]"};
0338       TString parYAxisName2D[nHistosPVParam2D] = {"y [cm]"};
0339       int nBinsX2D[nHistosPVParam2D] = {1000};
0340       float xMin2D[nHistosPVParam2D] = {-1.};
0341       float xMax2D[nHistosPVParam2D] = { 1.};
0342       int nBinsY2D[nHistosPVParam2D] = {1000};
0343       float yMin2D[nHistosPVParam2D] = {-1.};
0344       float yMax2D[nHistosPVParam2D] = { 1.};
0345       
0346       for(int iH=0; iH<nHistosPVParam; iH++)
0347       {
0348         hPVParam[iH]       = new TH1F(parName[iH].Data(),(GetDirectoryPath()+parName[iH]).Data(),
0349                                         nBins[iH],xMin[iH],xMax[iH]);
0350         hPVParam[iH]->GetXaxis()->SetTitle(parAxisName[iH].Data());
0351       }
0352 
0353       for(int iH=0; iH<nHistosPVParam2D; iH++)
0354       {
0355         hPVParam2D[iH]       = new TH2F(parName2D[iH].Data(),(GetDirectoryPath()+parName2D[iH]).Data(),
0356                                         nBinsX2D[iH],xMin2D[iH],xMax2D[iH],
0357                                         nBinsY2D[iH],yMin2D[iH],yMax2D[iH]);
0358         hPVParam2D[iH]->GetXaxis()->SetTitle(parXAxisName2D[iH].Data());
0359         hPVParam2D[iH]->GetYaxis()->SetTitle(parYAxisName2D[iH].Data());
0360         hPVParam2D[iH]->GetYaxis()->SetTitleOffset(1.0);
0361       }
0362       
0363       gDirectory->mkdir("Efficiency");
0364       gDirectory->cd("Efficiency");
0365       {
0366         TString effName[nPVefficiency] = {"effVsNMCPVTracks","effVsNMCPV","effVsNMCTracks","effVsNPVTracks","effVsNPV","effVsNTracks"};
0367         int nBinsEff[nPVefficiency]  = { 100 , 100 ,  100 ,  100 , 100 , 1000 };
0368         float xMinEff[nPVefficiency] = {   0.,   0.,    0.,    0.,   0.,    0.};
0369         float xMaxEff[nPVefficiency] = { 100., 100., 1000.,  100., 100., 1000.};
0370 
0371         gDirectory->mkdir("Signal");
0372         gDirectory->cd("Signal");
0373         {
0374           for( int iH=0; iH<nPVefficiency; iH++ ){
0375             hPVefficiency[0][iH]   = new TProfile( effName[iH].Data(), (GetDirectoryPath()+effName[iH]).Data(), nBinsEff[iH], xMinEff[iH], xMaxEff[iH]);
0376           }
0377         }
0378         gDirectory->cd(".."); //L1
0379 
0380         gDirectory->mkdir("Pileup");
0381         gDirectory->cd("Pileup");
0382         {
0383           for( int iH=0; iH<nPVefficiency; iH++ ){
0384             hPVefficiency[1][iH]   = new TProfile( effName[iH].Data(), (GetDirectoryPath()+effName[iH].Data()), nBinsEff[iH], xMinEff[iH], xMaxEff[iH]);
0385           }
0386         }
0387         gDirectory->cd(".."); //L1
0388         
0389         gDirectory->mkdir("Signal_MCReconstructable");
0390         gDirectory->cd("Signal_MCReconstructable");
0391         {
0392           for( int iH=0; iH<nPVefficiency; iH++ ){
0393             hPVefficiency[2][iH]   = new TProfile( effName[iH].Data(), (GetDirectoryPath()+effName[iH].Data()), nBinsEff[iH], xMinEff[iH], xMaxEff[iH]);
0394           }
0395         }
0396         gDirectory->cd(".."); //L1
0397         
0398         gDirectory->mkdir("Pileup_MCReconstructable");
0399         gDirectory->cd("Pileup_MCReconstructable");
0400         {
0401           for( int iH=0; iH<nPVefficiency; iH++ ){
0402             hPVefficiency[3][iH]   = new TProfile( effName[iH].Data(), (GetDirectoryPath()+effName[iH].Data()), nBinsEff[iH], xMinEff[iH], xMaxEff[iH]);
0403           }
0404         }
0405         gDirectory->cd(".."); //L1
0406       }      
0407       gDirectory->cd(".."); //L1
0408       
0409       gDirectory->mkdir("PVTracksQA");
0410       gDirectory->cd("PVTracksQA");
0411       {
0412         TString resTrPV = "resTrPV";
0413         TString pullTrPV = "pullTrPV";
0414         TString parNameTrPV[nFitPVTracksQA/2] = {"X","Y","Z","Px","Py","Pz"};
0415         int nBinsTrPV = 100;
0416         float xMaxTrPV[nFitPVTracksQA/2] = {0.5,0.5,0.5,0.05,0.05,0.05};
0417 
0418         for( int iH=0; iH<nFitPVTracksQA/2; iH++ ){
0419           hFitPVTracksQA[iH]   = new TH1F((resTrPV+parNameTrPV[iH]).Data(),
0420                                                   (GetDirectoryPath()+resTrPV+parNameTrPV[iH]).Data(), 
0421                                                   nBinsTrPV, -xMaxTrPV[iH],xMaxTrPV[iH]);
0422           hFitPVTracksQA[iH+nFitPVTracksQA/2] = new TH1F((pullTrPV+parNameTrPV[iH]).Data(),
0423                                                   (GetDirectoryPath()+pullTrPV+parNameTrPV[iH]).Data(), 
0424                                                   nBinsTrPV, -6,6);
0425         }
0426       }      
0427       gDirectory->cd(".."); //L1
0428       
0429       gDirectory->mkdir("Signal");
0430       gDirectory->cd("Signal");
0431       {
0432         gDirectory->mkdir("FitQA");
0433         gDirectory->cd("FitQA");
0434         {
0435           gDirectory->mkdir("FitQAvcNMCPVTracks");
0436           gDirectory->cd("FitQAvcNMCPVTracks");
0437           {
0438             for(int iHPV=0; iHPV<nHistosPV-1; ++iHPV){
0439               hPVFitQa2D[0][0][iHPV] = new TH2F(Table[iHPV].name.Data(),(GetDirectoryPath()+Table[iHPV].title).Data(),
0440                                                 500, 0., 5000.,
0441                                                 Table[iHPV].n, Table[iHPV].l, Table[iHPV].r);
0442             }
0443           }
0444           gDirectory->cd(".."); //FitQA
0445 
0446           gDirectory->mkdir("FitQAvcNPVTracks");
0447           gDirectory->cd("FitQAvcNPVTracks");
0448           {
0449             for(int iHPV=0; iHPV<nHistosPV-1; ++iHPV){
0450               hPVFitQa2D[0][1][iHPV] = new TH2F(Table[iHPV].name.Data(),(GetDirectoryPath()+Table[iHPV].title).Data(),
0451                                                 500, 0., 5000.,
0452                                                 Table[iHPV].n, Table[iHPV].l, Table[iHPV].r);
0453             }
0454           }
0455           gDirectory->cd(".."); //FitQA
0456           
0457           for(int iHPV=0; iHPV<nHistosPV; ++iHPV){
0458             hPVFitQa[0][iHPV] = new TH1F(Table[iHPV].name.Data(),(GetDirectoryPath()+Table[iHPV].title).Data(),
0459                                          Table[iHPV].n, Table[iHPV].l, Table[iHPV].r);
0460           }
0461         }
0462         gDirectory->cd(".."); //Signal
0463 
0464         for(int iH=0; iH<nHistosPVParam; iH++)
0465         {
0466           hPVParamSignal[iH] = new TH1F((parName[iH]).Data(),(GetDirectoryPath()+parName[iH]).Data(),
0467                                         nBins[iH],xMin[iH],xMax[iH]);
0468           hPVParamSignal[iH]->GetXaxis()->SetTitle(parAxisName[iH].Data());
0469         }
0470       }      
0471       gDirectory->cd(".."); //L1
0472 
0473       gDirectory->mkdir("Pileup");
0474       gDirectory->cd("Pileup");
0475       {
0476         gDirectory->mkdir("FitQA");
0477         gDirectory->cd("FitQA");
0478         {
0479           gDirectory->mkdir("FitQAvcNMCPVTracks");
0480           gDirectory->cd("FitQAvcNMCPVTracks");
0481           {
0482             for(int iHPV=0; iHPV<nHistosPV-1; ++iHPV){
0483               hPVFitQa2D[1][0][iHPV] = new TH2F(Table[iHPV].name.Data(),(GetDirectoryPath()+Table[iHPV].title).Data(),
0484                                                 500, 0., 5000.,
0485                                                 Table[iHPV].n, Table[iHPV].l, Table[iHPV].r);
0486             }
0487           }
0488           gDirectory->cd(".."); //FitQA
0489 
0490           gDirectory->mkdir("FitQAvcNPVTracks");
0491           gDirectory->cd("FitQAvcNPVTracks");
0492           {
0493             for(int iHPV=0; iHPV<nHistosPV-1; ++iHPV){
0494               hPVFitQa2D[1][1][iHPV] = new TH2F(Table[iHPV].name.Data(),(GetDirectoryPath()+Table[iHPV].title).Data(),
0495                                                 500, 0., 5000.,
0496                                                 Table[iHPV].n, Table[iHPV].l, Table[iHPV].r);
0497             }
0498           }
0499           gDirectory->cd(".."); //FitQA
0500           
0501           for(int iHPV=0; iHPV<nHistosPV; ++iHPV){
0502             hPVFitQa[1][iHPV] = new TH1F(Table[iHPV].name.Data(),(GetDirectoryPath()+Table[iHPV].title).Data(),
0503                                          Table[iHPV].n, Table[iHPV].l, Table[iHPV].r);
0504           }
0505         }
0506         gDirectory->cd(".."); //Signal
0507         
0508         for(int iH=0; iH<nHistosPVParam; iH++)
0509         {
0510           hPVParamPileup[iH] = new TH1F((parName[iH]).Data(),(GetDirectoryPath()+parName[iH]).Data(),
0511                                         nBins[iH],xMin[iH],xMax[iH]);
0512           hPVParamPileup[iH]->GetXaxis()->SetTitle(parAxisName[iH].Data());
0513         }
0514       }      
0515       gDirectory->cd(".."); //L1
0516       
0517       gDirectory->mkdir("Background");
0518       gDirectory->cd("Background");
0519       {
0520         for(int iH=0; iH<nHistosPVParam; iH++)
0521         {
0522           hPVParamBG[iH] = new TH1F((parName[iH]).Data(),(GetDirectoryPath()+parName[iH]).Data(),
0523                                         nBins[iH],xMin[iH],xMax[iH]);
0524           hPVParamBG[iH]->GetXaxis()->SetTitle(parAxisName[iH].Data());
0525         }
0526       }      
0527       gDirectory->cd(".."); //L1
0528       
0529       gDirectory->mkdir("Ghost");
0530       gDirectory->cd("Ghost");
0531       {
0532         for(int iH=0; iH<nHistosPVParam; iH++)
0533         {
0534           hPVParamGhost[iH] = new TH1F((parName[iH]).Data(),(GetDirectoryPath()+parName[iH]).Data(),
0535                                         nBins[iH],xMin[iH],xMax[iH]);
0536           hPVParamGhost[iH]->GetXaxis()->SetTitle(parAxisName[iH].Data());
0537         }
0538       }
0539       gDirectory->cd(".."); //L1
0540     }
0541     gDirectory->cd(".."); //L1
0542     gDirectory->mkdir("TrackParameters");
0543     gDirectory->cd("TrackParameters");
0544     {
0545       TString chi2Name = "Chi2Prim";
0546       for(int iPart=0; iPart < KFPartEfficiencies::nParticles; iPart++)
0547       {
0548         TString chi2NamePart = "Chi2Prim";
0549         chi2NamePart += "_";
0550         chi2NamePart += fParteff.partName[iPart].data();
0551         hTrackParameters[iPart] = new TH1F(chi2NamePart.Data(), (GetDirectoryPath()+chi2NamePart).Data(), 1000, 0, 100);
0552 
0553       }
0554       hTrackParameters[KFPartEfficiencies::nParticles  ] = new TH1F("Chi2Prim_total", (GetDirectoryPath()+TString("Chi2Prim_total")), 1000, 0, 100);
0555       hTrackParameters[KFPartEfficiencies::nParticles+1] = new TH1F("Chi2Prim_prim", (GetDirectoryPath()+TString("Chi2Prim_prim")), 1000, 0, 100);
0556       hTrackParameters[KFPartEfficiencies::nParticles+2] = new TH1F("Chi2Prim_sec", (GetDirectoryPath()+TString("Chi2Prim_sec")), 1000, 0, 100);
0557       hTrackParameters[KFPartEfficiencies::nParticles+3] = new TH1F("Chi2Prim_ghost", (GetDirectoryPath()+TString("Chi2Prim_ghost")), 1000, 0, 100);
0558       
0559       hTrackParameters[KFPartEfficiencies::nParticles+4] = new TH1F("ProbPrim_total", (GetDirectoryPath()+TString("ProbPrim_total")), 10000, 0, 1);
0560       hTrackParameters[KFPartEfficiencies::nParticles+5] = new TH1F("ProbPrim_prim", (GetDirectoryPath()+TString("ProbPrim_prim")), 10000, 0, 1);
0561       hTrackParameters[KFPartEfficiencies::nParticles+6] = new TH1F("ProbPrim_sec", (GetDirectoryPath()+TString("ProbPrim_sec")), 10000, 0, 1);
0562       hTrackParameters[KFPartEfficiencies::nParticles+7] = new TH1F("ProbPrim_ghost", (GetDirectoryPath()+TString("ProbPrim_ghost")), 10000, 0, 1);
0563     }
0564     gDirectory->cd(".."); //particle directory
0565     curdir->cd();    
0566   }
0567 }
0568 
0569 void KFParticlePerformanceBase::CreateFitHistograms(TH1F* histo[nFitQA], int iPart)
0570 {
0571   /** Creates 1D histograms with fit QA for decay with "iPart" number.
0572    ** \param[in,out] histo - array with pointers, for which the memory is allocated
0573    ** \param[in] iPart - number of the decay in the KF Particle Finder reconstruction scheme
0574    **/
0575   TString res = "res";
0576   TString pull = "pull";
0577   
0578   TString AxisNameResidual[nFitQA/2];
0579   TString AxisNamePull[nFitQA/2];
0580 
0581   AxisNameResidual[0] = "Residual (x^{reco} - x^{mc}) [cm]";
0582   AxisNameResidual[1] = "Residual (y^{reco} - y^{mc}) [cm]";
0583   AxisNameResidual[2] = "Residual (z^{reco} - z^{mc}) [cm]";
0584   AxisNameResidual[3] = "Residual (P_{x}^{reco} - P_{x}^{mc}) [GeV/c]";
0585   AxisNameResidual[4] = "Residual (P_{y}^{reco} - P_{y}^{mc}) [GeV/c]";
0586   AxisNameResidual[5] = "Residual (P_{z}^{reco} - P_{z}^{mc}) [GeV/c]";
0587   AxisNameResidual[6] = "Residual (E^{reco} - E^{mc}) [GeV/c^{2}]";
0588   AxisNameResidual[7] = "Residual (M^{reco} - M^{mc}) [GeV/c^{2}]"; 
0589 
0590   AxisNamePull[0] = "Pull x";
0591   AxisNamePull[1] = "Pull y";
0592   AxisNamePull[2] = "Pull z";
0593   AxisNamePull[3] = "Pull P_{x}";
0594   AxisNamePull[4] = "Pull P_{y}";
0595   AxisNamePull[5] = "Pull P_{z}";
0596   AxisNamePull[6] = "Pull E";
0597   AxisNamePull[7] = "Pull M";
0598   
0599   gDirectory->mkdir("FitQA");
0600   gDirectory->cd("FitQA");
0601   {
0602     TString parName[nFitQA/2] = {"X","Y","Z","Px","Py","Pz","E","M"};
0603     int nBins = 50;
0604     float xMax[nFitQA/2] = {0.15,0.15,1.2,0.02,0.02,0.15,0.15,0.006};
0605     float mult[nFitQA/2]={1.f,1.f,1.f,1.f,1.f,1.f,1.f,1.f};
0606     if(iPart>63 && iPart<75)
0607       for(int iMult=3; iMult<nFitQA/2; iMult++)
0608         mult[iMult] = 3;
0609     if(iPart>45 && iPart<64)
0610     {
0611 #ifdef CBM
0612       for(int iMult=0; iMult<3; iMult++)
0613         mult[iMult] = 0.03;
0614       for(int iMult=3; iMult<nFitQA/2; iMult++)
0615         mult[iMult] = 3;
0616 #else
0617       mult[2] = 0.1;
0618       for(int iMult=3; iMult<nFitQA/2; iMult++)
0619         mult[iMult] = 10;
0620       mult[5] = 2;
0621       mult[6] = 2;
0622 #endif
0623     }
0624     if(iPart==44 || iPart==45)
0625     {
0626       mult[0] = 0.25;
0627       mult[1] = 0.5;
0628       mult[2] = 0.15;
0629       for(int iMult=3; iMult<nFitQA/2; iMult++)
0630         mult[iMult] = 4;
0631     }
0632     
0633     for( int iH=0; iH<nFitQA/2; iH++ )
0634     {
0635       histo[iH]   = new TH1F((res+parName[iH]).Data(),
0636                              (GetDirectoryPath()+res+parName[iH]).Data(), 
0637                              nBins, -mult[iH]*xMax[iH],mult[iH]*xMax[iH]);
0638       histo[iH]->GetXaxis()->SetTitle(AxisNameResidual[iH].Data());
0639       histo[iH+8] = new TH1F((pull+parName[iH]).Data(),
0640                              (GetDirectoryPath()+pull+parName[iH]).Data(), 
0641                              nBins, -6,6);
0642       histo[iH+8]->GetXaxis()->SetTitle(AxisNamePull[iH+8].Data());
0643     }
0644   }
0645   gDirectory->cd("..");
0646 }
0647 
0648 void KFParticlePerformanceBase::CreateEfficiencyHistograms(TProfile* histo[3][nPartEfficiency], TProfile2D* histo2[3][nPartEfficiency2D])
0649 {
0650   /** Creates efficiency plots in the current ROOT folder.
0651    ** \param[in,out] histo - 1D efficiency plots
0652    ** \param[in,out] histo2 - 2D efficiency plots
0653    **/
0654   gDirectory->mkdir("Efficiency");
0655   gDirectory->cd("Efficiency");
0656   {//vs p, pt, y, z, c*tau, decay length, l, r
0657     TString partNameEff[nPartEfficiency] = {"EffVsP","EffVsPt","EffVsY","EffVsZ","EffVsCT","EffVsDL","EffVsL","EffVsR","EffVsMt" };
0658     TString partAxisNameEff[nPartEfficiency] = {"p [GeV/c]","p_{t} [GeV/c]",
0659                                                 "y", "z [cm]", "Life time c#tau [cm]", "Decay length [cm]", 
0660                                                 "L [cm]", "Rxy [cm]", "m_{t} [GeV/c^{2}]"};
0661 #ifdef CBM
0662     int nBinsEff[nPartEfficiency]  = { 100 , 100 , 100 ,  360 ,  100 ,  100 , 200 , 200 , 100 };
0663     float xMinEff[nPartEfficiency] = {   0.,   0.,  0.,  -10.,    0.,    0.,    0.,    0. , 0.};
0664     float xMaxEff[nPartEfficiency] = {  20.,  5.,   6.,   80.,  100.,  100.,  100.,  50. , 4.};
0665 #else
0666     int nBinsEff[nPartEfficiency]  = { 100 , 100 ,  30  ,   100 ,  100 ,  100 ,  100 ,  100 , 100  };
0667     float xMinEff[nPartEfficiency] = {   0.,   0.,  -1.5,   -10.,    0.,    0.,    0.,    0.,   0. };
0668     float xMaxEff[nPartEfficiency] = {  10.,  10.,   1.5,    10.,   30.,    5.,    1.,    1.,  10. };
0669 #endif
0670     TString effTypeName[3] = {"All particles",
0671                               "Reconstructable daughters",
0672                               "Reconstructed daughters"};
0673     
0674     for(int iEff=0; iEff<3; iEff++)
0675     {
0676       gDirectory->mkdir(effTypeName[iEff].Data());
0677       gDirectory->cd(effTypeName[iEff].Data());
0678       {
0679         for(int iH=0; iH<nPartEfficiency; iH++)
0680         {
0681           histo[iEff][iH] = new TProfile( partNameEff[iH].Data(), (GetDirectoryPath()+partAxisNameEff[iH]).Data(), nBinsEff[iH], xMinEff[iH], xMaxEff[iH]);
0682           histo[iEff][iH]->GetYaxis()->SetTitle("Efficiency");                  
0683           histo[iEff][iH]->GetYaxis()->SetTitleOffset(1.0);  
0684           histo[iEff][iH]->GetXaxis()->SetTitle(partAxisNameEff[iH].Data());
0685         }
0686         
0687         histo2[iEff][0] = new TProfile2D( "EffVsPtVsY", (GetDirectoryPath()+partAxisNameEff[2]+partAxisNameEff[1]).Data(), 
0688                                           nBinsEff[2], xMinEff[2], xMaxEff[2], nBinsEff[1], xMinEff[1], xMaxEff[1]);
0689         histo2[iEff][0]->GetZaxis()->SetTitle("Efficiency");
0690         histo2[iEff][0]->GetXaxis()->SetTitle(partAxisNameEff[2].Data());
0691         histo2[iEff][0]->GetYaxis()->SetTitle(partAxisNameEff[1].Data());
0692         histo2[iEff][0]->GetYaxis()->SetTitleOffset(1.0);
0693         
0694         histo2[iEff][1] = new TProfile2D( "EffVsMtVsY", (GetDirectoryPath()+partAxisNameEff[2]+partAxisNameEff[8]).Data(), 
0695                                           nBinsEff[2], xMinEff[2], xMaxEff[2], nBinsEff[8], xMinEff[8], xMaxEff[8]);
0696         histo2[iEff][1]->GetZaxis()->SetTitle("Efficiency");
0697         histo2[iEff][1]->GetXaxis()->SetTitle(partAxisNameEff[2].Data());
0698         histo2[iEff][1]->GetYaxis()->SetTitle(partAxisNameEff[8].Data());
0699         histo2[iEff][1]->GetYaxis()->SetTitleOffset(1.0);
0700       }
0701       gDirectory->cd("..");// particle directory / Efficiency
0702     }
0703   }
0704   gDirectory->cd("..");// particle directory 
0705 }
0706 
0707 void KFParticlePerformanceBase::CreateParameterHistograms(TH1F* histoParameters[KFPartEfficiencies::nParticles][nHistoPartParam],
0708                                                           TH2F *histoParameters2D[KFPartEfficiencies::nParticles][nHistoPartParam2D],
0709                                                           TH3F *histoParameters3D[KFPartEfficiencies::nParticles][nHistoPartParam3D],
0710                                                           int iPart, bool drawZR)
0711 {
0712   /** Creates histograms with parameter distributions for decay with "iPart" number.
0713    ** \param[in,out] histoParameters - 1D histograms
0714    ** \param[in,out] histoParameters2D - 2D histograms
0715    ** \param[in,out] histoParameters3D - 3D histograms
0716    ** \param[in] iPart - number of the decay in the KF Particle Finder reconstruction scheme
0717    ** \param[in] drawZR - flag showing if Z-R histogram should be created
0718    **/
0719   TString parName[nHistoPartParam] = {"M","p","p_{t}","y","DecayL","c#tau","chi2ndf","prob","#theta","phi","X","Y","Z","R", "L", "l/dl","m_{t}","Multiplicity"};
0720   TString parTitle[nHistoPartParam];
0721   TString parName2D[nHistoPartParam2D] = {"y-p_{t}", "Z-R", "Armenteros", "y-m_{t}"};
0722   TString parTitle2D[nHistoPartParam2D];
0723   TString parName3D[nHistoPartParam3D] = {"y-p_{t}-M", "y-m_{t}-M", "centrality-pt-M", "centrality-y-M", "centrality-mt-M", "ct-pt-M"};
0724   TString parTitle3D[nHistoPartParam3D];
0725   for(int iParam=0; iParam<nHistoPartParam; iParam++)
0726   {
0727     TString path = GetDirectoryPath();
0728     parTitle[iParam] = path + parName[iParam];
0729     if(iParam<nHistoPartParam2D)
0730       parTitle2D[iParam] = path + parName2D[iParam];
0731     if(iParam<nHistoPartParam3D)
0732       parTitle3D[iParam] = path + parName3D[iParam];
0733   }
0734   
0735   TString parAxisName[nHistoPartParam] = {"m [GeV/c^{2}]","p [GeV/c]","p_{t} [GeV/c]",
0736                                           "y","Decay length [cm]","Life time c#tau [cm]",
0737                                           "chi2/ndf","prob","#theta [rad]",
0738                                           "phi [rad]","x [cm]","y [cm]","z [cm]","Rxy [cm]", "L [cm]", "L/dL","m_{t} [GeV/c^{2}]","Multiplicity"};
0739 #ifdef CBM
0740   int nBins[nHistoPartParam] =  {1000, // M
0741                                   100, // p
0742                                   100, // pt
0743                                    40, // y
0744                                    60, // DecayL
0745                                    60, // ctau
0746                                   100, // chi2/ndf
0747                                   100, // prob
0748                                   100, // theta
0749                                   100, // phi
0750                                   200, // X
0751                                   200, // Y
0752                                   360, // Z
0753                                    60, // R
0754                                   140, // L
0755                                   200, // L/dL
0756                                   100, // Mt
0757                                   fParteff.partMaxMult[iPart]+1};
0758   float xMin[nHistoPartParam] = { fParteff.partMHistoMin[iPart], // M
0759                                   0.f, // p
0760                                   0.f, // pt
0761                                   0.f, // y
0762                                 -10.f, // DecayL
0763                                 -10.f, // ctau
0764                                   0.f, // chi2/ndf
0765                                   0.f, // prob
0766                                  -2.f, // theta
0767                                  -2.f, // phi
0768                                 -50.f, // X
0769                                 -50.f, // Y
0770                                 -10.f, // Z
0771                                   0.f, // R
0772                                   0.f, // L
0773                                  -1.f, // L/dL
0774                                   0.f, // Mt
0775                                  -0.5f };
0776   float xMax[nHistoPartParam] = { fParteff.partMHistoMax[iPart], // M
0777                                   20.f, // p
0778                                    5.f, // pt
0779                                    4.f, // y
0780                                   50.f, // DecayL
0781                                   50.f, // ctau
0782                                   20.f, // chi2/ndf
0783                                    1.f, // prob
0784                                    2.f, // theta
0785                                    2.f, // phi
0786                                   50.f, // X
0787                                   50.f, // Y
0788                                   80.f, // Z
0789                                   30.f, // R
0790                                   70.f, // L
0791                                   35.f, // L/dL
0792                                   4.f, // Mt
0793                                   float(fParteff.partMaxMult[iPart])+0.5f};
0794 #else
0795   int nBins[nHistoPartParam] = {1000, // M
0796                                  100, // p
0797                                  100, // pt
0798                                   30, // y
0799                                   60, // DecayL
0800                                   60, // ctau
0801                                  100, // chi2/ndf
0802                                  100, // prob
0803                                  100, // theta
0804                                  100, // phi
0805                                  100, // X
0806                                  100, // Y
0807                                  100, // Z
0808                                  100, // R
0809                                  100, // L
0810                                 1000, // L/dL
0811                                  100, // Mt
0812                                  fParteff.partMaxMult[iPart]+1};
0813   float xMin[nHistoPartParam] = { fParteff.partMHistoMin[iPart], // M
0814                                   0.f, // p
0815                                   0.f, // pt
0816                                 -1.5f, // y
0817                                 -10.f, // DecayL
0818                                 -10.f, // ctau
0819                                   0.f, // chi2/ndf
0820                                   0.f, // prob
0821                                   0.f, // theta
0822                              -3.1416f, // phi
0823                                 -10.f, // X
0824                                 -10.f, // Y
0825                                 -30.f, // Z
0826                                   0.f, // R
0827                                   0.f, // L
0828                                  -1.f, // L/dL
0829                                   0.f, // Mt
0830                                  -0.5f };
0831   float xMax[nHistoPartParam] = { fParteff.partMHistoMax[iPart], // M
0832                                   10.f, // p
0833                                   10.f, // pt
0834                                   1.5f, // y
0835                                   50.f, // DecayL
0836                                   50.f, // ctau
0837                                   20.f, // chi2/ndf
0838                                    1.f, // prob
0839                                3.1416f, // theta
0840                                3.1416f, // phi
0841                                   10.f, // X
0842                                   10.f, // Y
0843                                   30.f, // Z
0844                                   50.f, // R
0845                                   50.f, // L
0846                                   35.f, // L/dL
0847                                   10.f, // Mt
0848                                   float(fParteff.partMaxMult[iPart])+0.5f};
0849   if(iPart < 9)
0850   {
0851     xMin[10] =-50; xMin[11] =-50; xMin[12] =-100;
0852     xMax[10] = 50; xMax[11] = 50; xMax[12] = 100; xMax[13] = 50; xMax[14] = 50; 
0853   }
0854 #endif
0855   for(int iH=0; iH<nHistoPartParam; iH++)
0856   {
0857     histoParameters[iPart][iH] = new TH1F(parName[iH].Data(),parTitle[iH].Data(),
0858                                           nBins[iH],xMin[iH],xMax[iH]);
0859     histoParameters[iPart][iH]->GetXaxis()->SetTitle(parAxisName[iH].Data());
0860   }
0861 
0862   histoParameters2D[iPart][0] = new TH2F(parName2D[0].Data(),parTitle2D[0].Data(),
0863                                     nBins[3],xMin[3],xMax[3],
0864                                     nBins[2],xMin[2],xMax[2]);
0865   histoParameters2D[iPart][0]->GetXaxis()->SetTitle("y");
0866   histoParameters2D[iPart][0]->GetYaxis()->SetTitle("p_{t} [GeV/c]");
0867   histoParameters2D[iPart][0]->GetYaxis()->SetTitleOffset(1.0);
0868 
0869   if(drawZR)
0870   {
0871     histoParameters2D[iPart][1] = new TH2F(parName2D[1].Data(),parTitle2D[1].Data(),
0872                                       nBins[12],xMin[12],xMax[12],
0873                                       nBins[13],xMin[13],xMax[13]);
0874     histoParameters2D[iPart][1]->GetXaxis()->SetTitle("Z [cm]");
0875     histoParameters2D[iPart][1]->GetYaxis()->SetTitle("R [cm]");
0876     histoParameters2D[iPart][1]->GetYaxis()->SetTitleOffset(1.0);
0877   }
0878   else
0879     histoParameters2D[iPart][1] = NULL;
0880   
0881   //create armenteros plot
0882   if(IsCollectArmenteros(iPart))
0883   {
0884     histoParameters2D[iPart][2] = new TH2F(parName2D[2].Data(),parTitle2D[2].Data(),
0885                                            50, -1.f, 1.f,
0886                                           150,  0.f, 1.f);
0887     histoParameters2D[iPart][2]->GetXaxis()->SetTitle("#alpha (p_{L}^{+}-p_{L}^{-})/(p_{L}^{+}+p_{L}^{-})");
0888     histoParameters2D[iPart][2]->GetYaxis()->SetTitle("q_{t} [GeV/c]");
0889     histoParameters2D[iPart][2]->GetYaxis()->SetTitleOffset(1.0);
0890   }
0891   else
0892     histoParameters2D[iPart][2] = NULL;
0893   //create y-mt plot
0894   histoParameters2D[iPart][3] = new TH2F(parName2D[3].Data(),parTitle2D[3].Data(),
0895                                          nBins[3],xMin[3], xMax[3],     //y
0896                                          nBins[16],xMin[16],xMax[16]); //Mt
0897   histoParameters2D[iPart][3]->GetXaxis()->SetTitle("y");
0898   histoParameters2D[iPart][3]->GetYaxis()->SetTitle("m_{t} [GeV/c]");
0899   histoParameters2D[iPart][3]->GetYaxis()->SetTitleOffset(1.0);
0900   
0901   
0902   if( histoParameters3D && IsCollect3DHistogram(iPart) )
0903   {
0904     histoParameters3D[iPart][0] = new TH3F(parName3D[0].Data(),parTitle3D[0].Data(),
0905                                       nBins[3],xMin[3],xMax[3],
0906                                       nBins[2],xMin[2],xMax[2],
0907                                       nBins[0],xMin[0],xMax[0]);
0908     histoParameters3D[iPart][0]->GetXaxis()->SetTitle("y");
0909     histoParameters3D[iPart][0]->GetYaxis()->SetTitle("p_{t} [GeV/c]");
0910     histoParameters3D[iPart][0]->GetYaxis()->SetTitleOffset(1.0);
0911     histoParameters3D[iPart][0]->GetZaxis()->SetTitle("M");
0912     
0913     histoParameters3D[iPart][1] = new TH3F(parName3D[1].Data(),parTitle3D[1].Data(),
0914                                            nBins[3],xMin[3],xMax[3],
0915                                            nBins[16],xMin[16],xMax[16],
0916                                            nBins[0],xMin[0],xMax[0]);
0917     histoParameters3D[iPart][1]->GetXaxis()->SetTitle("y");
0918     histoParameters3D[iPart][1]->GetYaxis()->SetTitle("m_{t} [GeV/c]");
0919     histoParameters3D[iPart][1]->GetYaxis()->SetTitleOffset(1.0);
0920     histoParameters3D[iPart][1]->GetZaxis()->SetTitle("M");
0921     
0922     int centralityHisto[3] = {2,3,16};
0923     for(int iCH = 0; iCH<3; iCH++)
0924     {
0925       histoParameters3D[iPart][2+iCH] = new TH3F(parName3D[2+iCH].Data(),parTitle3D[2+iCH].Data(),
0926                                                  10,0.,10.,
0927                                                  nBins[centralityHisto[iCH]],xMin[centralityHisto[iCH]],xMax[centralityHisto[iCH]],
0928                                                  nBins[0],xMin[0],xMax[0]);
0929       histoParameters3D[iPart][2+iCH]->GetXaxis()->SetTitle("centrality bin");
0930       histoParameters3D[iPart][2+iCH]->GetYaxis()->SetTitle(parAxisName[centralityHisto[iCH]]);
0931       histoParameters3D[iPart][2+iCH]->GetYaxis()->SetTitleOffset(1.0);
0932       histoParameters3D[iPart][2+iCH]->GetZaxis()->SetTitle("M");
0933     }
0934     
0935     histoParameters3D[iPart][5] = new TH3F(parName3D[5].Data(),parTitle3D[5].Data(),
0936                                            nBins[5],xMin[5],xMax[5],
0937                                            nBins[2],xMin[2],xMax[2],
0938                                            nBins[0],xMin[0],xMax[0]);
0939     histoParameters3D[iPart][5]->GetXaxis()->SetTitle("c#tau [cm]");
0940     histoParameters3D[iPart][5]->GetYaxis()->SetTitle("p_{t} [GeV/c]");
0941     histoParameters3D[iPart][5]->GetYaxis()->SetTitleOffset(1.0);
0942     histoParameters3D[iPart][5]->GetZaxis()->SetTitle("M");
0943   }
0944   else if(histoParameters3D)
0945   {
0946     histoParameters3D[iPart][0] = NULL;
0947     histoParameters3D[iPart][1] = NULL;
0948     for(int iCH = 0; iCH<3; iCH++)
0949       histoParameters3D[iPart][2+iCH] = NULL;
0950     histoParameters3D[iPart][5] = NULL;
0951   }
0952 }
0953 
0954 bool KFParticlePerformanceBase::IsCollectZRHistogram(int iParticle) const
0955 {
0956   /** Checks if Z-R histogram for decay "iParticle" should be created. */
0957   return (abs(fParteff.partPDG[iParticle]) == 310 ||
0958           abs(fParteff.partPDG[iParticle]) == 3122 ||
0959           abs(fParteff.partPDG[iParticle]) == 3312 ||
0960           abs(fParteff.partPDG[iParticle]) == 3334 ||
0961           abs(fParteff.partPDG[iParticle]) == 22) && fStoreMCHistograms && fStoreZRHistograms;
0962 }
0963 
0964 bool KFParticlePerformanceBase::IsCollect3DHistogram(int iParticle) const
0965 {
0966   /** Checks if 3D histograms for decay "iParticle" should be created. */
0967   return abs(fParteff.partPDG[iParticle]) == 310 ||
0968          abs(fParteff.partPDG[iParticle]) == 3122 ||
0969          abs(fParteff.partPDG[iParticle]) == 3312 ||
0970          abs(fParteff.partPDG[iParticle]) == 3334 ||
0971          abs(fParteff.partPDG[iParticle]) == 3003 ||
0972          abs(fParteff.partPDG[iParticle]) == 3103 ||
0973          abs(fParteff.partPDG[iParticle]) == 3004 ||
0974          abs(fParteff.partPDG[iParticle]) == 3005 ||
0975 #ifdef CBM
0976          abs(fParteff.partPDG[iParticle]) == 7003112 ||
0977          abs(fParteff.partPDG[iParticle]) == 7003222 ||
0978          abs(fParteff.partPDG[iParticle]) == 7003312 ||
0979          abs(fParteff.partPDG[iParticle]) == 8003222 ||
0980          abs(fParteff.partPDG[iParticle]) == 9000321;
0981 #else
0982          abs(fParteff.partPDG[iParticle]) == 421 ||
0983          abs(fParteff.partPDG[iParticle]) == 429 ||
0984          abs(fParteff.partPDG[iParticle]) == 426 ||
0985          abs(fParteff.partPDG[iParticle]) == 411 ||
0986          abs(fParteff.partPDG[iParticle]) == 431 ||
0987          abs(fParteff.partPDG[iParticle]) == 4122 ||
0988          abs(fParteff.partPDG[iParticle]) == 521 ||
0989          abs(fParteff.partPDG[iParticle]) == 511;
0990 #endif
0991 }
0992 
0993 bool KFParticlePerformanceBase::IsCollectArmenteros(int iParticle) const
0994 {
0995   /** Checks if Armenteros-Podoliansky plot for decay "iParticle" should be created. */
0996   return abs(fParteff.partPDG[iParticle]) == 310 ||
0997          abs(fParteff.partPDG[iParticle]) == 3122 ||
0998          abs(fParteff.partPDG[iParticle]) == 3312 ||
0999          abs(fParteff.partPDG[iParticle]) == 3334 ||
1000          abs(fParteff.partPDG[iParticle]) == 22 ||
1001          abs(fParteff.partPDG[iParticle]) == 111 ||
1002          abs(fParteff.partPDG[iParticle]) == 3003 ||
1003          abs(fParteff.partPDG[iParticle]) == 3103 ||
1004          abs(fParteff.partPDG[iParticle]) == 3004 ||
1005          abs(fParteff.partPDG[iParticle]) == 3005 ||
1006          abs(fParteff.partPDG[iParticle]) == 3203 ||
1007          abs(fParteff.partPDG[iParticle]) == 3008 ||
1008          abs(fParteff.partPDG[iParticle]) == 3000 ||
1009          abs(fParteff.partPDG[iParticle]) == 333 ||
1010 #ifdef CBM
1011          abs(fParteff.partPDG[iParticle]) == 7003112 ||
1012          abs(fParteff.partPDG[iParticle]) == 7003222 ||
1013          abs(fParteff.partPDG[iParticle]) == 7003312 ||
1014          abs(fParteff.partPDG[iParticle]) == 8003222 ||
1015          abs(fParteff.partPDG[iParticle]) == 9000321;
1016 #else
1017          abs(fParteff.partPDG[iParticle]) == 421 ||
1018          abs(fParteff.partPDG[iParticle]) == 420 ||
1019          abs(fParteff.partPDG[iParticle]) == 426 ||
1020          abs(fParteff.partPDG[iParticle]) == 521 ||
1021          abs(fParteff.partPDG[iParticle]) == 511;        
1022 #endif
1023 }
1024 
1025 void KFParticlePerformanceBase::CreateParameterSubfolder(TString folderName, 
1026                                                          TH1F* histoParameters[nParametersSet][KFPartEfficiencies::nParticles][nHistoPartParam],
1027                                                          TH2F* histoParameters2D[nParametersSet][KFPartEfficiencies::nParticles][nHistoPartParam2D],
1028                                                          TH1F* histoFit[KFPartEfficiencies::nParticles][nFitQA], int iPart, bool withWrongPVHypothesis)
1029 {
1030   /** Creates all subfolders in the current ROOT directory for the current decay channel. */
1031   gDirectory->mkdir(folderName.Data());
1032   gDirectory->cd(folderName.Data());
1033   {
1034     gDirectory->mkdir("Signal");
1035     gDirectory->cd("Signal");
1036     {
1037       CreateParameterHistograms(histoParameters[1], histoParameters2D[1], 0, iPart);
1038     }
1039     gDirectory->cd("..");
1040     if(withWrongPVHypothesis)
1041     {
1042       gDirectory->mkdir("WrongPVHypothesis");
1043       gDirectory->cd("WrongPVHypothesis");
1044       {
1045         CreateParameterHistograms(histoParameters[4], histoParameters2D[4], 0, iPart);
1046       }
1047       gDirectory->cd("..");
1048     }
1049     gDirectory->mkdir("Background");
1050     gDirectory->cd("Background");
1051     {
1052       CreateParameterHistograms(histoParameters[2], histoParameters2D[2], 0, iPart);
1053     }
1054     gDirectory->cd("..");
1055     gDirectory->mkdir("Ghost");
1056     gDirectory->cd("Ghost");
1057     {
1058       CreateParameterHistograms(histoParameters[3], histoParameters2D[3], 0, iPart);
1059     }
1060     gDirectory->cd("..");
1061     
1062     CreateParameterHistograms(histoParameters[0], histoParameters2D[0], 0, iPart);
1063     if(histoFit!=0)
1064       CreateFitHistograms(histoFit[iPart], iPart);
1065   }
1066   gDirectory->cd("..");
1067 }
1068 
1069 TString KFParticlePerformanceBase::GetDirectoryPath()
1070 {
1071   /** Returns the path to the current folder. It is used as an addition to the histogram name. */
1072   TString path = gDirectory->GetPath();
1073   int fileNamePosition = path.Index("Finder/");
1074   path.Remove(0, fileNamePosition+7);
1075   path.ReplaceAll("Particles/", "");
1076   path.ReplaceAll("/Parameters", "");
1077   path+=" ";
1078   return path;
1079 }
1080 
1081 #endif //DO_TPCCATRACKER_EFF_PERFORMANCE
1082