ROOT logo
//____________________________________________________________________
//
// $Id: DrawSDigits.C 20907 2007-09-25 08:44:03Z cholm $
//
// Script that contains a class to draw eloss from hits, versus ADC
// counts from digits, using the AliFMDInputHits class in the util library. 
//
// It draws the energy loss versus the p/(mq^2).  It can be overlayed
// with the Bethe-Bloc curve to show how the simulation behaves
// relative to the expected. 
//
// Use the script `Compile.C' to compile this class using ACLic. 
//
#include <TH1D.h>
#include <TH2D.h>
#include <TProfile2D.h>
#include <TCanvas.h>
#include <TMath.h>
#include <AliFMDSDigit.h>
#include <AliFMDInput.h>
#include <AliFMDEdepMap.h>
#include <AliFMDGeometry.h>
#include <iostream>
#include <TStyle.h>
#include <TLatex.h>
#include <TArrayF.h>
#include <AliLog.h>

/** @class DrawSDigits
    @brief Draw hit energy loss versus digit ADC
    @code 
    Root> .L Compile.C
    Root> Compile("DrawSDigits.C")
    Root> DrawSDigits c
    Root> c.Run();
    @endcode
    @ingroup FMD_script
 */
class DrawSDigits : public AliFMDInput
{
private:
  TH1D* fAdc; // Histogram 
  TProfile2D* fPrimRatio[5]; // Histograms
public:
  void Idx2Det(UShort_t idx, UShort_t& d, Char_t& r) const
  {
    d = 0;
    r = '\0';
    switch (idx) { 
    case 0: d = 1; r = 'I'; break;
    case 1: d = 2; r = 'I'; break;
    case 2: d = 2; r = 'O'; break;
    case 3: d = 3; r = 'I'; break;
    case 4: d = 3; r = 'O'; break;
    }
  }
  Short_t Det2Idx(UShort_t d, Char_t r) const
  { 
    Short_t idx = -1;
    switch (d) { 
    case 1: idx = 0; break;
    case 2: idx = 1; break;
    case 3: idx = 3; break;
    }
    return (idx + ((r == 'O' || r == 'o') ? 1 : 0));
  }
  
    
  //__________________________________________________________________
  DrawSDigits(Int_t m=1100, Double_t amin=-0.5, Double_t amax=1023.5) 
    : AliFMDInput("galice.root")
  { 
    AddLoad(kSDigits);
    AddLoad(kGeometry);
    fAdc = new TH1D("adc", "ADC", m, amin, amax);
    fAdc->SetXTitle("ADC value");
    fAdc->Sumw2();

    Int_t   eN   = 50;
    Float_t eMin = -4;
    Float_t eMax = 6;
    Int_t   pN   = 40;
    Float_t pMin = -.5;
    Float_t pMax = 359.5;
    for (UShort_t i = 0; i < 5; i++) { 
      UShort_t d;
      Char_t   r;
      Idx2Det(i, d, r);
      fPrimRatio[i] = new TProfile2D(Form("primRatio_fmd%d%c", d, r), 
				     Form("Primary/Total - FMD%d%c", d, r), 
				     eN, eMin, eMax, pN, pMin, pMax);
      fPrimRatio[i]->SetXTitle("#eta");
      fPrimRatio[i]->SetYTitle("#varphi");
      fPrimRatio[i]->SetZTitle("M_{ch,primary}/M_{ch,total}");
      fPrimRatio[i]->GetXaxis()->SetTitleFont(132);
      fPrimRatio[i]->GetYaxis()->SetTitleFont(132);
      fPrimRatio[i]->GetZaxis()->SetTitleFont(132);
      fPrimRatio[i]->GetXaxis()->SetLabelFont(132);
      fPrimRatio[i]->GetYaxis()->SetLabelFont(132);
      fPrimRatio[i]->GetZaxis()->SetLabelFont(132);
      fPrimRatio[i]->GetXaxis()->SetTitleSize(.06);
      fPrimRatio[i]->GetYaxis()->SetTitleSize(.06);
      fPrimRatio[i]->GetZaxis()->SetTitleSize(.06);
      fPrimRatio[i]->GetXaxis()->SetLabelSize(.06);
      fPrimRatio[i]->GetYaxis()->SetLabelSize(.06);
      fPrimRatio[i]->GetZaxis()->SetLabelSize(.06);
      fPrimRatio[i]->GetXaxis()->SetNdivisions(10);
      fPrimRatio[i]->GetYaxis()->SetNdivisions(10);
      fPrimRatio[i]->GetZaxis()->SetNdivisions(10);
    }
  }
  Bool_t Init() 
  {
    Bool_t ret = AliFMDInput::Init();
    AliFMDGeometry::Instance()->Init();
    AliFMDGeometry::Instance()->InitTransformations();
    return ret;
  }
    
  //__________________________________________________________________
  Bool_t ProcessSDigit(AliFMDSDigit* digit)
  {
    if (!digit) return kTRUE;
    fAdc->Fill(digit->Counts());
    digit->Print("lp");
    if (digit->NParticles() == 0) return kTRUE;
    

    Short_t primIdx = Det2Idx(digit->Detector(), digit->Ring());
    if (primIdx < 0) return kTRUE;
    
    AliFMDGeometry* geom = AliFMDGeometry::Instance();
    Double_t x, y, z;
    geom->Detector2XYZ(digit->Detector(), 
		       digit->Ring(), 
		       digit->Sector(), 
		       digit->Strip(), 
		       x, y, z);
    // We should get the primary vertex and do 
    // z     += fCurrentVertex;
    Double_t phi   =  TMath::ATan2(y, x) * 180. / TMath::Pi();
    Double_t r     =  TMath::Sqrt(y * y + x * x);
    Double_t theta =  TMath::ATan2(r, z);
    Double_t eta   = -TMath::Log(TMath::Tan(theta / 2));
    Double_t ratio = digit->NPrimaries() / digit->NParticles();
    if (phi < 0) phi += 360;
    fPrimRatio[primIdx]->Fill(eta, phi, ratio);
    return kTRUE;
  }
  //__________________________________________________________________
  Bool_t Finish()
  {
    gStyle->SetPalette(1);
    gStyle->SetOptTitle(0);
    gStyle->SetCanvasColor(0);
    gStyle->SetCanvasBorderSize(0);
    gStyle->SetPadColor(0);
    gStyle->SetPadBorderSize(0);
    
    TCanvas* c1 = new TCanvas("fmdSdigitSpectra", "FMD SDigit spectra");
    fAdc->SetStats(kFALSE);
    if (fAdc->GetEntries() != 0) 
      fAdc->Scale(1. / fAdc->GetEntries());
    fAdc->Draw();
    c1->Modified();
    c1->Update();
    c1->cd();

    TCanvas* c2 = new TCanvas("fmdSdigitPrims", "FMD SDigit # prim/# ntotal",
			      800, 800);
    c2->SetTopMargin(0);
    c2->SetRightMargin(0);
    c2->Divide(2, 3, 0, 0);
    for (UShort_t i = 0; i < 5; i++) { 
      Int_t        np = i+1+(i == 0 ? 0 : 1);
      TVirtualPad* p   = c2->cd(np);
      p->SetFillColor(0);
      p->SetTopMargin((np <= 2) ? 0.05 : 0);
      p->SetLeftMargin((np % 2 == 1) ? 0.20 : 0);
      p->SetRightMargin((np % 2 == 1) ? 0 : 0.20);
      p->SetBottomMargin((np >= 5) ? 0.15 : 0);
      fPrimRatio[i]->SetStats(kFALSE);
      fPrimRatio[i]->Draw(Form("COL%c", (np % 2 == 1) ? ' ' : 'Z'));
      UShort_t d;
      Char_t   r;
      Idx2Det(i, d, r);
      TLatex* text = new TLatex(0, 180, Form("FMD%d%c", d, r));
      text->SetTextFont(132);
      text->SetTextAlign(22);
      text->SetTextSize(0.08);
      text->Draw();
    }
    c2->Modified();
    c2->Update();
    c2->cd();

    return kTRUE;
  }

  ClassDef(DrawSDigits,0);
};

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