ROOT logo
#include "TObject.h"
#include "AliRunLoader.h"
#include "AliRun.h"
#include "AliLoader.h"
#include "AliMFT.h"
#include "TClonesArray.h"
#include "AliMFTCluster.h"
#include "AliMFTSegmentation.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2D.h"
#include "AliLog.h"
#include "TString.h"

#include "AliMFTClusterQA.h"

//====================================================================================================================================================
//
// Class for the analysis of the MFT clusters (a.k.a. rec points). Few QA histograms are created
//
// Contact author: antonio.uras@cern.ch
//
//====================================================================================================================================================

ClassImp(AliMFTClusterQA)

//====================================================================================================================================================

AliMFTClusterQA::AliMFTClusterQA():
  TObject(),
  fMFTLoader(0),
  fRunLoader(0),
  fMFT(0),
  fNPlanes(0),
  fNEvents(0),
  fEv(0),
  fFileOut(0),
  fReadDir(0),
  fOutDir(0)
{
  
  // default constructor

  for (Int_t iPlane=0; iPlane<fNMaxPlanes; iPlane++) {
    fHistNClustersPerEvent[iPlane]     = 0;
    fHistNPixelsPerCluster[iPlane]     = 0;
    fHistClusterSizeX[iPlane]          = 0; 
    fHistClusterSizeY[iPlane]          = 0;
    fHistClusterRadialPosition[iPlane] = 0;
    fClusterScatterPlotXY[iPlane]      = 0;
  }

}

//====================================================================================================================================================

void AliMFTClusterQA::Init(Char_t *readDir, Char_t *outDir, Int_t nEventsToAnalyze) {

  fReadDir = readDir;
  fOutDir  = outDir;

  fRunLoader = AliRunLoader::Open(Form("%s/galice.root", fReadDir.Data()));
  gAlice = fRunLoader->GetAliRun();
  if (!gAlice) fRunLoader->LoadgAlice();
  fNEvents = fRunLoader->GetNumberOfEvents();
  if (nEventsToAnalyze>0) fNEvents = TMath::Min(fNEvents, nEventsToAnalyze);

  fMFT = (AliMFT*) gAlice->GetDetector("MFT"); 
  fNPlanes = fMFT->GetSegmentation()->GetNPlanes();

  BookHistos();

  fMFTLoader = fRunLoader->GetDetectorLoader("MFT");
  fMFTLoader -> LoadRecPoints("READ");

}

//====================================================================================================================================================

Bool_t AliMFTClusterQA::LoadNextEvent() {

  if (fEv>=fNEvents) return kFALSE;
  AliDebug(1, Form("event %5d",fEv));
  
  fRunLoader->GetEvent(fEv);
  fEv++;

  if (!fMFTLoader->TreeR()->GetEvent()) return kTRUE;

  for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
    Int_t nClusters = fMFT->GetRecPointsList(iPlane)->GetEntries();
    fHistNClustersPerEvent[iPlane] -> Fill(nClusters);
    fClusterScatterPlotXY[iPlane]  -> Fill(0., 0.);    // "scaler" bin
    AliDebug(1,Form("nClusters = %5d", nClusters));
    for (Int_t iCluster=0; iCluster<nClusters; iCluster++) {
      AliMFTCluster *cluster = (AliMFTCluster*) (fMFT->GetRecPointsList(iPlane))->At(iCluster);
      fHistNPixelsPerCluster[iPlane] -> Fill(cluster->GetSize());
      fHistClusterSizeX[iPlane]      -> Fill(cluster->GetErrX()*1.e4);   // converted in microns
      fHistClusterSizeY[iPlane]      -> Fill(cluster->GetErrY()*1.e4);   // converted in microns
      fHistClusterRadialPosition[iPlane] -> Fill(TMath::Sqrt(cluster->GetX()*cluster->GetX()+cluster->GetY()*cluster->GetY()));
      fClusterScatterPlotXY[iPlane]      -> Fill(cluster->GetX(), cluster->GetY());
    }
  }

  return kTRUE;

}  

//====================================================================================================================================================

void AliMFTClusterQA::BookHistos() {

  fFileOut = new TFile(Form("%s/MFT.RecPoints.QA.root",fOutDir.Data()), "recreate");

  for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {

    fHistNClustersPerEvent[iPlane] = new TH1D(Form("fHistNClustersPerEvent_Pl%02d",iPlane), 
					      Form("Number of clusters per event in Plane%02d",iPlane),
					      25000, -0.5, 24999.5);

    fHistNPixelsPerCluster[iPlane] = new TH1D(Form("fHistNPixelsPerCluster_Pl%02d",iPlane), 
					      Form("Number of pixels per cluster in Plane%02d",iPlane),
					      15, -0.5, 14.5);

    fHistClusterSizeX[iPlane]      = new TH1D(Form("fHistClusterSizeX_Pl%02d",iPlane), 
					      Form("#Deltax for clusters in Plane%02d",iPlane), 
					      100, 0., 100.);
    
    fHistClusterSizeY[iPlane]      = new TH1D(Form("fHistClusterSizeY_Pl%02d",iPlane), 
					      Form("#Deltay for clusters in Plane%02d",iPlane), 
					      100, 0., 100.);
    
    fHistNClustersPerEvent[iPlane] -> SetXTitle("N_{clusters} per Event");
    fHistNPixelsPerCluster[iPlane] -> SetXTitle("N_{pixels} per Cluster");
    fHistClusterSizeX[iPlane]      -> SetXTitle("#Deltax  [#mum]");
    fHistClusterSizeY[iPlane]      -> SetXTitle("#Deltay  [#mum]");

    fHistNClustersPerEvent[iPlane] -> Sumw2();
    fHistNPixelsPerCluster[iPlane] -> Sumw2();
    fHistClusterSizeX[iPlane]      -> Sumw2();
    fHistClusterSizeY[iPlane]      -> Sumw2();

    //------------------------------------------------------------

    Int_t rMax = Int_t(10.*(fMFT->GetSegmentation()->GetPlane(iPlane)->GetRMaxSupport()));  // expressed in mm

    fHistClusterRadialPosition[iPlane] = new TH1D(Form("fHistClusterRadialPosition_Pl%02d",iPlane),
						  Form("Cluster radial position (Plane%02d)",iPlane),
						  rMax, 0, Double_t(rMax)/10.);

    fClusterScatterPlotXY[iPlane] = new TH2D(Form("fClusterScatterPlotXY_Pl%02d",iPlane),
					     Form("Cluster scatter plot (Plane%02d)",iPlane),
					     2*rMax+1, (-rMax-0.5)/10., (rMax+0.5)/10., 2*rMax+1, (-rMax-0.5)/10., (rMax+0.5)/10.);
    
    fHistClusterRadialPosition[iPlane] -> SetXTitle("R  [cm]");
    fClusterScatterPlotXY[iPlane]      -> SetXTitle("X  [cm]");
    fClusterScatterPlotXY[iPlane]      -> SetYTitle("Y  [cm]");

    fHistClusterRadialPosition[iPlane] -> Sumw2();
    fClusterScatterPlotXY[iPlane]      -> Sumw2();
    
  }
  
}

//====================================================================================================================================================

void AliMFTClusterQA::Terminate() {

  AliInfo("Writing QA histos...");

  // ---- equalize radial clusters distribution ----------------------

  for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
    for (Int_t iBin=1; iBin<=fHistClusterRadialPosition[iPlane]->GetNbinsX(); iBin++) {
      Double_t rMin = fHistClusterRadialPosition[iPlane]->GetBinLowEdge(iBin);        // in cm
      Double_t rMax = fHistClusterRadialPosition[iPlane]->GetBinWidth(iBin) + rMin;   // in cm
      Double_t area = 100.*TMath::Pi()*(rMax*rMax - rMin*rMin);                       // in mm^2
      Double_t density = fHistClusterRadialPosition[iPlane]->GetBinContent(iBin)/area;
      fHistClusterRadialPosition[iPlane]->SetBinContent(iBin, density);
      fHistClusterRadialPosition[iPlane]->SetBinError(iBin, fHistClusterRadialPosition[iPlane]->GetBinError(iBin)/area);
    }
    fHistClusterRadialPosition[iPlane] -> SetBinContent(1, fEv);      // "scaler" bin
  }

  // -----------------------------------------------------------------

  fFileOut->cd();

  for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
    fHistNClustersPerEvent[iPlane]     -> Write();
    fHistNPixelsPerCluster[iPlane]     -> Write();
    fHistClusterSizeX[iPlane]          -> Write();
    fHistClusterSizeY[iPlane]          -> Write();
    fHistClusterRadialPosition[iPlane] -> Write();
    fClusterScatterPlotXY[iPlane]      -> Write();
  }

  fFileOut -> Close();

}

//====================================================================================================================================================
 AliMFTClusterQA.cxx:1
 AliMFTClusterQA.cxx:2
 AliMFTClusterQA.cxx:3
 AliMFTClusterQA.cxx:4
 AliMFTClusterQA.cxx:5
 AliMFTClusterQA.cxx:6
 AliMFTClusterQA.cxx:7
 AliMFTClusterQA.cxx:8
 AliMFTClusterQA.cxx:9
 AliMFTClusterQA.cxx:10
 AliMFTClusterQA.cxx:11
 AliMFTClusterQA.cxx:12
 AliMFTClusterQA.cxx:13
 AliMFTClusterQA.cxx:14
 AliMFTClusterQA.cxx:15
 AliMFTClusterQA.cxx:16
 AliMFTClusterQA.cxx:17
 AliMFTClusterQA.cxx:18
 AliMFTClusterQA.cxx:19
 AliMFTClusterQA.cxx:20
 AliMFTClusterQA.cxx:21
 AliMFTClusterQA.cxx:22
 AliMFTClusterQA.cxx:23
 AliMFTClusterQA.cxx:24
 AliMFTClusterQA.cxx:25
 AliMFTClusterQA.cxx:26
 AliMFTClusterQA.cxx:27
 AliMFTClusterQA.cxx:28
 AliMFTClusterQA.cxx:29
 AliMFTClusterQA.cxx:30
 AliMFTClusterQA.cxx:31
 AliMFTClusterQA.cxx:32
 AliMFTClusterQA.cxx:33
 AliMFTClusterQA.cxx:34
 AliMFTClusterQA.cxx:35
 AliMFTClusterQA.cxx:36
 AliMFTClusterQA.cxx:37
 AliMFTClusterQA.cxx:38
 AliMFTClusterQA.cxx:39
 AliMFTClusterQA.cxx:40
 AliMFTClusterQA.cxx:41
 AliMFTClusterQA.cxx:42
 AliMFTClusterQA.cxx:43
 AliMFTClusterQA.cxx:44
 AliMFTClusterQA.cxx:45
 AliMFTClusterQA.cxx:46
 AliMFTClusterQA.cxx:47
 AliMFTClusterQA.cxx:48
 AliMFTClusterQA.cxx:49
 AliMFTClusterQA.cxx:50
 AliMFTClusterQA.cxx:51
 AliMFTClusterQA.cxx:52
 AliMFTClusterQA.cxx:53
 AliMFTClusterQA.cxx:54
 AliMFTClusterQA.cxx:55
 AliMFTClusterQA.cxx:56
 AliMFTClusterQA.cxx:57
 AliMFTClusterQA.cxx:58
 AliMFTClusterQA.cxx:59
 AliMFTClusterQA.cxx:60
 AliMFTClusterQA.cxx:61
 AliMFTClusterQA.cxx:62
 AliMFTClusterQA.cxx:63
 AliMFTClusterQA.cxx:64
 AliMFTClusterQA.cxx:65
 AliMFTClusterQA.cxx:66
 AliMFTClusterQA.cxx:67
 AliMFTClusterQA.cxx:68
 AliMFTClusterQA.cxx:69
 AliMFTClusterQA.cxx:70
 AliMFTClusterQA.cxx:71
 AliMFTClusterQA.cxx:72
 AliMFTClusterQA.cxx:73
 AliMFTClusterQA.cxx:74
 AliMFTClusterQA.cxx:75
 AliMFTClusterQA.cxx:76
 AliMFTClusterQA.cxx:77
 AliMFTClusterQA.cxx:78
 AliMFTClusterQA.cxx:79
 AliMFTClusterQA.cxx:80
 AliMFTClusterQA.cxx:81
 AliMFTClusterQA.cxx:82
 AliMFTClusterQA.cxx:83
 AliMFTClusterQA.cxx:84
 AliMFTClusterQA.cxx:85
 AliMFTClusterQA.cxx:86
 AliMFTClusterQA.cxx:87
 AliMFTClusterQA.cxx:88
 AliMFTClusterQA.cxx:89
 AliMFTClusterQA.cxx:90
 AliMFTClusterQA.cxx:91
 AliMFTClusterQA.cxx:92
 AliMFTClusterQA.cxx:93
 AliMFTClusterQA.cxx:94
 AliMFTClusterQA.cxx:95
 AliMFTClusterQA.cxx:96
 AliMFTClusterQA.cxx:97
 AliMFTClusterQA.cxx:98
 AliMFTClusterQA.cxx:99
 AliMFTClusterQA.cxx:100
 AliMFTClusterQA.cxx:101
 AliMFTClusterQA.cxx:102
 AliMFTClusterQA.cxx:103
 AliMFTClusterQA.cxx:104
 AliMFTClusterQA.cxx:105
 AliMFTClusterQA.cxx:106
 AliMFTClusterQA.cxx:107
 AliMFTClusterQA.cxx:108
 AliMFTClusterQA.cxx:109
 AliMFTClusterQA.cxx:110
 AliMFTClusterQA.cxx:111
 AliMFTClusterQA.cxx:112
 AliMFTClusterQA.cxx:113
 AliMFTClusterQA.cxx:114
 AliMFTClusterQA.cxx:115
 AliMFTClusterQA.cxx:116
 AliMFTClusterQA.cxx:117
 AliMFTClusterQA.cxx:118
 AliMFTClusterQA.cxx:119
 AliMFTClusterQA.cxx:120
 AliMFTClusterQA.cxx:121
 AliMFTClusterQA.cxx:122
 AliMFTClusterQA.cxx:123
 AliMFTClusterQA.cxx:124
 AliMFTClusterQA.cxx:125
 AliMFTClusterQA.cxx:126
 AliMFTClusterQA.cxx:127
 AliMFTClusterQA.cxx:128
 AliMFTClusterQA.cxx:129
 AliMFTClusterQA.cxx:130
 AliMFTClusterQA.cxx:131
 AliMFTClusterQA.cxx:132
 AliMFTClusterQA.cxx:133
 AliMFTClusterQA.cxx:134
 AliMFTClusterQA.cxx:135
 AliMFTClusterQA.cxx:136
 AliMFTClusterQA.cxx:137
 AliMFTClusterQA.cxx:138
 AliMFTClusterQA.cxx:139
 AliMFTClusterQA.cxx:140
 AliMFTClusterQA.cxx:141
 AliMFTClusterQA.cxx:142
 AliMFTClusterQA.cxx:143
 AliMFTClusterQA.cxx:144
 AliMFTClusterQA.cxx:145
 AliMFTClusterQA.cxx:146
 AliMFTClusterQA.cxx:147
 AliMFTClusterQA.cxx:148
 AliMFTClusterQA.cxx:149
 AliMFTClusterQA.cxx:150
 AliMFTClusterQA.cxx:151
 AliMFTClusterQA.cxx:152
 AliMFTClusterQA.cxx:153
 AliMFTClusterQA.cxx:154
 AliMFTClusterQA.cxx:155
 AliMFTClusterQA.cxx:156
 AliMFTClusterQA.cxx:157
 AliMFTClusterQA.cxx:158
 AliMFTClusterQA.cxx:159
 AliMFTClusterQA.cxx:160
 AliMFTClusterQA.cxx:161
 AliMFTClusterQA.cxx:162
 AliMFTClusterQA.cxx:163
 AliMFTClusterQA.cxx:164
 AliMFTClusterQA.cxx:165
 AliMFTClusterQA.cxx:166
 AliMFTClusterQA.cxx:167
 AliMFTClusterQA.cxx:168
 AliMFTClusterQA.cxx:169
 AliMFTClusterQA.cxx:170
 AliMFTClusterQA.cxx:171
 AliMFTClusterQA.cxx:172
 AliMFTClusterQA.cxx:173
 AliMFTClusterQA.cxx:174
 AliMFTClusterQA.cxx:175
 AliMFTClusterQA.cxx:176
 AliMFTClusterQA.cxx:177
 AliMFTClusterQA.cxx:178
 AliMFTClusterQA.cxx:179
 AliMFTClusterQA.cxx:180
 AliMFTClusterQA.cxx:181
 AliMFTClusterQA.cxx:182
 AliMFTClusterQA.cxx:183
 AliMFTClusterQA.cxx:184
 AliMFTClusterQA.cxx:185
 AliMFTClusterQA.cxx:186
 AliMFTClusterQA.cxx:187
 AliMFTClusterQA.cxx:188
 AliMFTClusterQA.cxx:189
 AliMFTClusterQA.cxx:190
 AliMFTClusterQA.cxx:191
 AliMFTClusterQA.cxx:192
 AliMFTClusterQA.cxx:193
 AliMFTClusterQA.cxx:194
 AliMFTClusterQA.cxx:195
 AliMFTClusterQA.cxx:196
 AliMFTClusterQA.cxx:197
 AliMFTClusterQA.cxx:198
 AliMFTClusterQA.cxx:199
 AliMFTClusterQA.cxx:200
 AliMFTClusterQA.cxx:201
 AliMFTClusterQA.cxx:202
 AliMFTClusterQA.cxx:203