ROOT logo
//____________________________________________________________________
//
// $Id$
//
// 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 <TH2D.h>
#include <AliFMDHit.h>
#include <AliFMDDigit.h>
#include <AliFMDRecPoint.h>
#include <AliFMDInput.h>
#include <AliFMDEdepMap.h>
#include <AliFMDFloatMap.h>
#include <iostream>
#include <TStyle.h>
#include <TArrayF.h>
#include <TCanvas.h>
#include <TParticle.h>
#include <TClonesArray.h>
#include <TTree.h>
#include <AliStack.h>
#include <AliLog.h>
#include <TF1.h>

//____________________________________________________________________
/** @class DrawHitsRecs
    @brief Draw hit energy loss versus rec point mult
    @code 
    Root> .L Compile.C
    Root> Compile("DrawHitsRecs.C")
    Root> DrawHitsRecs c
    Root> c.Run();
    @endcode
    @ingroup FMD_script
 */
class DrawHitsRecs : public AliFMDInput
{
private:
  TH2D* fHitEvsAdc;
  TH2D* fHitEvsRecM;  // Histogram 
  TH2D* fHitEvsRecE;  // Histogram 
  TH1D* fDiffE;       // Histogram 
  TH2D* fHitsVsRecM;  // Histogram 
  TH2D* fDiffM;       // Histogram 
  TH1*  fHitEloss;
  TH1*  fRecEloss;
  AliFMDEdepMap  fMap;
  AliFMDFloatMap fEta;
  AliFMDFloatMap fPhi;
  AliFMDFloatMap fMult;
  Bool_t fPrimary;
public:
  //__________________________________________________________________
  DrawHitsRecs(Bool_t primary=kFALSE,
	       Int_t n=900, Double_t emin=1e-3, Double_t emax=10, 
	       Int_t m=21, Double_t mmin=-0.5, Double_t mmax=20.5) 
  { 
    fPrimary = primary;
    AddLoad(kRecPoints);
    AddLoad(kHits);
    AddLoad(kDigits);
    if (fPrimary) AddLoad(kKinematics);
    TArrayF eloss(MakeLogScale(n, emin, emax));
    TArrayF mults(m+1);
    mults[0] = mmin;
    for (Int_t i = 1; i < m+1; i++) mults[i] = mults[i-1] + (mmax-mmin)/m;

    fHitEvsAdc  = new TH2D("hitEvsAdc", "#Delta E_{sim} vs. ADC",
			   n, emin, emax, 1025, -.5, 1024.5);
    fHitEvsAdc->SetXTitle("#Delta E_{sim} [MeV]");
    fHitEvsAdc->SetYTitle("ADC");
    fHitEvsAdc->Sumw2();
    
    fHitEvsRecM = new TH2D("hitEvsRecM", "#Delta E_{sim} vs. M_{rec}", 
			   eloss.fN-1, eloss.fArray, mults.fN-1, mults.fArray);
    fHitEvsRecM->SetXTitle("#Delta E_{sim} [MeV]");
    fHitEvsRecM->SetYTitle("M_{rec}");
    fHitEvsRecM->Sumw2();

    fHitEvsRecE = new TH2D("hitEvsRecE", "#Delta E_{sim} vs. #Delta E_{rec}", 
			    n, emin, emax, n, emin, emax);
    fHitEvsRecE->SetXTitle("#Delta E_{sim} [MeV]");
    fHitEvsRecE->SetYTitle("#Delta E_{rec} [MeV]");
    fHitEvsRecE->Sumw2();


    fDiffE = new TH1D("diffE", 
		      "#frac{#Delta E_{sim}-#Delta E_{rec}}{#Delta E_{sim}}", 
		      1100, -1, 1.1);
    fDiffE->SetXTitle("#frac{#Delta E_{sim}-#Delta E_{rec}}{#Delta E_{sim}}");
    fDiffE->Sumw2();
    
    Double_t omin = mmin; // -.5;
    Double_t omax = mmax; // 7.5;
    Int_t    o    = m;    // 8;
    fHitsVsRecM = new TH2D("hitsVsStrip", "# of Hits vs. M_{rec}",
			   o, omin, omax, m, mmin, mmax);
    fHitsVsRecM->SetXTitle("# of Hits");
    fHitsVsRecM->SetYTitle("M_{rec}");
    fHitsVsRecM->Sumw2();

    fDiffM = new TH2D("diffM", "M_{sim} - M_{rec}", 
		      41, -20.5, 20.5, 70, 1.5, 5);
    // 36, -TMath::Pi(),TMath::Pi());
    fDiffM->SetXTitle("M_{sim} - M_{rec}");
    fDiffM->SetYTitle("|#eta|");
    fDiffM->Sumw2();
    // fDiffM->SetYTitle("Detector");

    fHitEloss = new TH1D("hitEloss","#frac{#Delta E_{sim}}{#Delta x} (MeV/cm)", 
			 100, 0, 10);
    fHitEloss->SetFillColor(2);
    fHitEloss->SetFillStyle(3001);
    fHitEloss->Sumw2();

    fRecEloss = new TH1D("recEloss","#frac{#Delta E_{rec}}{#Delta x} (MeV/cm)", 
			 100, 0, 10);
    fRecEloss->SetFillColor(4);
    fRecEloss->SetFillStyle(3001);
    fRecEloss->Sumw2();
  }
  //__________________________________________________________________
  /** Begining of event
      @param ev Event number
      @return @c false on error */
  Bool_t Begin(Int_t ev) 
  {
    fMap.Reset();
    return AliFMDInput::Begin(ev);
  }
  //__________________________________________________________________
  Bool_t ProcessHit(AliFMDHit* hit, TParticle*) 
  {
    // Cache the energy loss 
    if (!hit) return kFALSE;
    UShort_t det = hit->Detector();
    Char_t   rng = hit->Ring();
    UShort_t sec = hit->Sector();
    UShort_t str = hit->Strip();
    if (str > 511) {
      AliWarning(Form("Bad strip number %d in hit", str));
      return kTRUE;
    }
    if (fPrimary) {
      TParticle* kine = fStack->Particle(hit->Track());
      if (!kine) return kTRUE;
      if (!kine->IsPrimary()) return kTRUE;
    }
    
    if (hit->Edep()/hit->Length() > 0.1) 
      fHitEloss->Fill(hit->Edep() / hit->Length());
    fMap(det, rng, sec, str).fEdep += hit->Edep();
    fMap(det, rng, sec, str).fN++;
    return kTRUE;
  }

  //__________________________________________________________________
  Bool_t ProcessDigit(AliFMDDigit* digit)
  {
    if (!digit) return kTRUE;

    UShort_t det = digit->Detector();
    Char_t   rng = digit->Ring();
    UShort_t sec = digit->Sector();
    UShort_t str = digit->Strip();
    if (str > 511) {
      AliWarning(Form("Bad strip number %d in digit", str));
      return kTRUE;
    }
    Double_t edep = fMap(det, rng, sec, str).fEdep;
    if (edep > 0) fHitEvsAdc->Fill(edep, digit->Counts());
      
    return kTRUE;
  }

  //__________________________________________________________________
  Bool_t ProcessRecPoint(AliFMDRecPoint* single) 
  {
    if (!single) return kTRUE;
    UShort_t det = single->Detector();
    Char_t   rng = single->Ring();
    UShort_t sec = single->Sector();
    UShort_t str = single->Strip();
    if (str > 511) {
      AliWarning(Form("Bad strip number %d in single", str));
      return kTRUE;
    }
    Double_t edep = fMap(det, rng, sec, str).fEdep;
    Int_t    nhit = fMap(det, rng, sec, str).fN;
    if (edep > 0) {
      fHitEvsRecM->Fill(edep, single->Particles());
      fHitEvsRecE->Fill(edep, single->Edep());
      fDiffE->Fill((single->Edep() - edep) / edep);
    }
    if (nhit > 0) fHitsVsRecM->Fill(nhit, single->Particles());
    fDiffM->Fill(nhit - single->Particles(), TMath::Abs(single->Eta()));
    if (single->Edep()/.03 > 0.1) fRecEloss->Fill(single->Edep() / 0.0300);
    return kTRUE;
  }
  //__________________________________________________________________
  Bool_t Finish()
  {
    gStyle->SetPalette(1);
    gStyle->SetOptTitle(0);
    gStyle->SetCanvasColor(0);
    gStyle->SetCanvasBorderSize(0);
    gStyle->SetPadColor(0);
    gStyle->SetPadBorderSize(0);
    TCanvas* c = 0;

    c = new TCanvas("c0", fHitEvsAdc->GetTitle());
    fHitEvsAdc->SetStats(kFALSE);
    fHitEvsAdc->Draw("COLZ");

    c = new TCanvas("c1", fHitEvsRecM->GetTitle());
    fHitEvsRecM->SetStats(kFALSE);
    fHitEvsRecM->Draw("COLZ");

    c = new TCanvas("c2", fHitEvsRecE->GetTitle());
    fHitEvsRecE->SetStats(kFALSE);
    fHitEvsRecE->Draw("COLZ");

    c = new TCanvas("c3", fDiffE->GetTitle());
    c->SetLogz();
    fDiffE->Draw();

    c = new TCanvas("c4", fHitsVsRecM->GetTitle());
    c->SetLogz();
    fHitsVsRecM->SetStats(kFALSE);
    fHitsVsRecM->Draw("COLZ");

    c = new TCanvas("c5", fDiffM->GetTitle());
    fDiffM->SetFillColor(2);
    fDiffM->SetFillStyle(3001);
    c->SetLogz();
    fDiffM->Draw("colz");

    c = new TCanvas("c6", "Hit Eloss, Reco Eloss");
    fRecEloss->Scale(1./fRecEloss->GetEntries());
    fRecEloss->Draw("hist e");
    fRecEloss->Fit("landau", "", "SAME", 2, 4);
    TF1* recResp = new TF1(*fRecEloss->GetFunction("landau"));
    fHitEloss->Scale(1./fHitEloss->GetEntries());
    fHitEloss->Draw("same hist e");
    fHitEloss->Fit("landau", "", "SAME", 2, 10);
    TF1* hitResp = new TF1(*fHitEloss->GetFunction("landau"));
    std::cout << "From hits: MPV=" 
	      << hitResp->GetParameter(1) << "+/-"
	      << hitResp->GetParError(1) << "  Width=" 
	      << hitResp->GetParameter(2) << "+/-"
	      << hitResp->GetParError(2) << "\nFrom recs: MPV=" 
	      << recResp->GetParameter(1) << "+/-"
	      << recResp->GetParError(1) << "  Width=" 
	      << recResp->GetParameter(2) << "+/-"
	      << recResp->GetParError(2) << "\nRatio: MPV(hit/rec)="
	      << hitResp->GetParameter(1) / recResp->GetParameter(1) 
	      << std::endl;
    c->SetLogy();

    return kTRUE;
  }

  ClassDef(DrawHitsRecs,0);
};

//____________________________________________________________________
//
// EOF
//
 DrawHitsRecs.C:1
 DrawHitsRecs.C:2
 DrawHitsRecs.C:3
 DrawHitsRecs.C:4
 DrawHitsRecs.C:5
 DrawHitsRecs.C:6
 DrawHitsRecs.C:7
 DrawHitsRecs.C:8
 DrawHitsRecs.C:9
 DrawHitsRecs.C:10
 DrawHitsRecs.C:11
 DrawHitsRecs.C:12
 DrawHitsRecs.C:13
 DrawHitsRecs.C:14
 DrawHitsRecs.C:15
 DrawHitsRecs.C:16
 DrawHitsRecs.C:17
 DrawHitsRecs.C:18
 DrawHitsRecs.C:19
 DrawHitsRecs.C:20
 DrawHitsRecs.C:21
 DrawHitsRecs.C:22
 DrawHitsRecs.C:23
 DrawHitsRecs.C:24
 DrawHitsRecs.C:25
 DrawHitsRecs.C:26
 DrawHitsRecs.C:27
 DrawHitsRecs.C:28
 DrawHitsRecs.C:29
 DrawHitsRecs.C:30
 DrawHitsRecs.C:31
 DrawHitsRecs.C:32
 DrawHitsRecs.C:33
 DrawHitsRecs.C:34
 DrawHitsRecs.C:35
 DrawHitsRecs.C:36
 DrawHitsRecs.C:37
 DrawHitsRecs.C:38
 DrawHitsRecs.C:39
 DrawHitsRecs.C:40
 DrawHitsRecs.C:41
 DrawHitsRecs.C:42
 DrawHitsRecs.C:43
 DrawHitsRecs.C:44
 DrawHitsRecs.C:45
 DrawHitsRecs.C:46
 DrawHitsRecs.C:47
 DrawHitsRecs.C:48
 DrawHitsRecs.C:49
 DrawHitsRecs.C:50
 DrawHitsRecs.C:51
 DrawHitsRecs.C:52
 DrawHitsRecs.C:53
 DrawHitsRecs.C:54
 DrawHitsRecs.C:55
 DrawHitsRecs.C:56
 DrawHitsRecs.C:57
 DrawHitsRecs.C:58
 DrawHitsRecs.C:59
 DrawHitsRecs.C:60
 DrawHitsRecs.C:61
 DrawHitsRecs.C:62
 DrawHitsRecs.C:63
 DrawHitsRecs.C:64
 DrawHitsRecs.C:65
 DrawHitsRecs.C:66
 DrawHitsRecs.C:67
 DrawHitsRecs.C:68
 DrawHitsRecs.C:69
 DrawHitsRecs.C:70
 DrawHitsRecs.C:71
 DrawHitsRecs.C:72
 DrawHitsRecs.C:73
 DrawHitsRecs.C:74
 DrawHitsRecs.C:75
 DrawHitsRecs.C:76
 DrawHitsRecs.C:77
 DrawHitsRecs.C:78
 DrawHitsRecs.C:79
 DrawHitsRecs.C:80
 DrawHitsRecs.C:81
 DrawHitsRecs.C:82
 DrawHitsRecs.C:83
 DrawHitsRecs.C:84
 DrawHitsRecs.C:85
 DrawHitsRecs.C:86
 DrawHitsRecs.C:87
 DrawHitsRecs.C:88
 DrawHitsRecs.C:89
 DrawHitsRecs.C:90
 DrawHitsRecs.C:91
 DrawHitsRecs.C:92
 DrawHitsRecs.C:93
 DrawHitsRecs.C:94
 DrawHitsRecs.C:95
 DrawHitsRecs.C:96
 DrawHitsRecs.C:97
 DrawHitsRecs.C:98
 DrawHitsRecs.C:99
 DrawHitsRecs.C:100
 DrawHitsRecs.C:101
 DrawHitsRecs.C:102
 DrawHitsRecs.C:103
 DrawHitsRecs.C:104
 DrawHitsRecs.C:105
 DrawHitsRecs.C:106
 DrawHitsRecs.C:107
 DrawHitsRecs.C:108
 DrawHitsRecs.C:109
 DrawHitsRecs.C:110
 DrawHitsRecs.C:111
 DrawHitsRecs.C:112
 DrawHitsRecs.C:113
 DrawHitsRecs.C:114
 DrawHitsRecs.C:115
 DrawHitsRecs.C:116
 DrawHitsRecs.C:117
 DrawHitsRecs.C:118
 DrawHitsRecs.C:119
 DrawHitsRecs.C:120
 DrawHitsRecs.C:121
 DrawHitsRecs.C:122
 DrawHitsRecs.C:123
 DrawHitsRecs.C:124
 DrawHitsRecs.C:125
 DrawHitsRecs.C:126
 DrawHitsRecs.C:127
 DrawHitsRecs.C:128
 DrawHitsRecs.C:129
 DrawHitsRecs.C:130
 DrawHitsRecs.C:131
 DrawHitsRecs.C:132
 DrawHitsRecs.C:133
 DrawHitsRecs.C:134
 DrawHitsRecs.C:135
 DrawHitsRecs.C:136
 DrawHitsRecs.C:137
 DrawHitsRecs.C:138
 DrawHitsRecs.C:139
 DrawHitsRecs.C:140
 DrawHitsRecs.C:141
 DrawHitsRecs.C:142
 DrawHitsRecs.C:143
 DrawHitsRecs.C:144
 DrawHitsRecs.C:145
 DrawHitsRecs.C:146
 DrawHitsRecs.C:147
 DrawHitsRecs.C:148
 DrawHitsRecs.C:149
 DrawHitsRecs.C:150
 DrawHitsRecs.C:151
 DrawHitsRecs.C:152
 DrawHitsRecs.C:153
 DrawHitsRecs.C:154
 DrawHitsRecs.C:155
 DrawHitsRecs.C:156
 DrawHitsRecs.C:157
 DrawHitsRecs.C:158
 DrawHitsRecs.C:159
 DrawHitsRecs.C:160
 DrawHitsRecs.C:161
 DrawHitsRecs.C:162
 DrawHitsRecs.C:163
 DrawHitsRecs.C:164
 DrawHitsRecs.C:165
 DrawHitsRecs.C:166
 DrawHitsRecs.C:167
 DrawHitsRecs.C:168
 DrawHitsRecs.C:169
 DrawHitsRecs.C:170
 DrawHitsRecs.C:171
 DrawHitsRecs.C:172
 DrawHitsRecs.C:173
 DrawHitsRecs.C:174
 DrawHitsRecs.C:175
 DrawHitsRecs.C:176
 DrawHitsRecs.C:177
 DrawHitsRecs.C:178
 DrawHitsRecs.C:179
 DrawHitsRecs.C:180
 DrawHitsRecs.C:181
 DrawHitsRecs.C:182
 DrawHitsRecs.C:183
 DrawHitsRecs.C:184
 DrawHitsRecs.C:185
 DrawHitsRecs.C:186
 DrawHitsRecs.C:187
 DrawHitsRecs.C:188
 DrawHitsRecs.C:189
 DrawHitsRecs.C:190
 DrawHitsRecs.C:191
 DrawHitsRecs.C:192
 DrawHitsRecs.C:193
 DrawHitsRecs.C:194
 DrawHitsRecs.C:195
 DrawHitsRecs.C:196
 DrawHitsRecs.C:197
 DrawHitsRecs.C:198
 DrawHitsRecs.C:199
 DrawHitsRecs.C:200
 DrawHitsRecs.C:201
 DrawHitsRecs.C:202
 DrawHitsRecs.C:203
 DrawHitsRecs.C:204
 DrawHitsRecs.C:205
 DrawHitsRecs.C:206
 DrawHitsRecs.C:207
 DrawHitsRecs.C:208
 DrawHitsRecs.C:209
 DrawHitsRecs.C:210
 DrawHitsRecs.C:211
 DrawHitsRecs.C:212
 DrawHitsRecs.C:213
 DrawHitsRecs.C:214
 DrawHitsRecs.C:215
 DrawHitsRecs.C:216
 DrawHitsRecs.C:217
 DrawHitsRecs.C:218
 DrawHitsRecs.C:219
 DrawHitsRecs.C:220
 DrawHitsRecs.C:221
 DrawHitsRecs.C:222
 DrawHitsRecs.C:223
 DrawHitsRecs.C:224
 DrawHitsRecs.C:225
 DrawHitsRecs.C:226
 DrawHitsRecs.C:227
 DrawHitsRecs.C:228
 DrawHitsRecs.C:229
 DrawHitsRecs.C:230
 DrawHitsRecs.C:231
 DrawHitsRecs.C:232
 DrawHitsRecs.C:233
 DrawHitsRecs.C:234
 DrawHitsRecs.C:235
 DrawHitsRecs.C:236
 DrawHitsRecs.C:237
 DrawHitsRecs.C:238
 DrawHitsRecs.C:239
 DrawHitsRecs.C:240
 DrawHitsRecs.C:241
 DrawHitsRecs.C:242
 DrawHitsRecs.C:243
 DrawHitsRecs.C:244
 DrawHitsRecs.C:245
 DrawHitsRecs.C:246
 DrawHitsRecs.C:247
 DrawHitsRecs.C:248
 DrawHitsRecs.C:249
 DrawHitsRecs.C:250
 DrawHitsRecs.C:251
 DrawHitsRecs.C:252
 DrawHitsRecs.C:253
 DrawHitsRecs.C:254
 DrawHitsRecs.C:255
 DrawHitsRecs.C:256
 DrawHitsRecs.C:257
 DrawHitsRecs.C:258
 DrawHitsRecs.C:259
 DrawHitsRecs.C:260
 DrawHitsRecs.C:261
 DrawHitsRecs.C:262
 DrawHitsRecs.C:263
 DrawHitsRecs.C:264
 DrawHitsRecs.C:265
 DrawHitsRecs.C:266
 DrawHitsRecs.C:267
 DrawHitsRecs.C:268
 DrawHitsRecs.C:269
 DrawHitsRecs.C:270
 DrawHitsRecs.C:271
 DrawHitsRecs.C:272
 DrawHitsRecs.C:273
 DrawHitsRecs.C:274
 DrawHitsRecs.C:275
 DrawHitsRecs.C:276
 DrawHitsRecs.C:277