]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG0/dNdEta/AlidNdEtaAnalysisMCSelector.cxx
moved AliSelector, AliSelectorRL to STEER
[u/mrichter/AliRoot.git] / PWG0 / dNdEta / AlidNdEtaAnalysisMCSelector.cxx
CommitLineData
dc740de4 1/* $Id$ */
2
4dd2ad81 3#include "AlidNdEtaAnalysisMCSelector.h"
4
5#include <TStyle.h>
6#include <TSystem.h>
7#include <TCanvas.h>
8#include <TParticle.h>
9#include <TParticlePDG.h>
10#include <TVector3.h>
5af55649 11#include <TH1F.h>
dc740de4 12#include <TH3F.h>
7029240a 13#include <TTree.h>
16e24ca3 14#include <TFile.h>
4dd2ad81 15
16#include <AliLog.h>
17#include <AliGenEventHeader.h>
7029240a 18#include <AliHeader.h>
0bd1f8a0 19#include <AliStack.h>
4dd2ad81 20
16e24ca3 21#include "dNdEta/dNdEtaAnalysis.h"
45e97e28 22#include "AliPWG0Helper.h"
4dd2ad81 23
dc740de4 24
4dd2ad81 25ClassImp(AlidNdEtaAnalysisMCSelector)
26
dc740de4 27AlidNdEtaAnalysisMCSelector::AlidNdEtaAnalysisMCSelector() :
16e24ca3 28 AliSelectorRL(),
29 fdNdEtaAnalysis(0),
5af55649 30 fVertex(0),
31 fPartEta(0),
4c351225 32 fPartPt(0),
5af55649 33 fEvents(0)
4dd2ad81 34{
35 //
36 // Constructor. Initialization of pointers
37 //
38}
39
40AlidNdEtaAnalysisMCSelector::~AlidNdEtaAnalysisMCSelector()
41{
42 //
43 // Destructor
44 //
45}
46
16e24ca3 47void AlidNdEtaAnalysisMCSelector::SlaveBegin(TTree * tree)
48{
49 // The SlaveBegin() function is called after the Begin() function.
50 // When running with PROOF SlaveBegin() is called on each slave server.
51 // The tree argument is deprecated (on PROOF 0 is passed).
52
53 AliSelectorRL::SlaveBegin(tree);
54
55 fdNdEtaAnalysis = new dNdEtaAnalysis("dndeta", "dndeta");
0bd1f8a0 56 fVertex = new TH3F("vertex_check", "vertex_check", 50, -50, 50, 50, -50, 50, 50, -50, 50);
57 fPartEta = new TH1F("dndeta_check", "dndeta_check", 120, -6, 6);
4c351225 58 fPartPt = new TH1F("dndeta_check_pt", "dndeta_check_pt", 1000, 0, 10);
0bd1f8a0 59 fPartEta->Sumw2();
16e24ca3 60}
61
944f0536 62void AlidNdEtaAnalysisMCSelector::Init(TTree *tree)
63{
16e24ca3 64 AliSelectorRL::Init(tree);
944f0536 65
4c351225 66 tree->SetBranchStatus("*", 0);
944f0536 67}
68
4dd2ad81 69Bool_t AlidNdEtaAnalysisMCSelector::Process(Long64_t entry)
70{
16e24ca3 71 // fill the dNdEtaAnalysis class from the monte carlo
4dd2ad81 72
16e24ca3 73 if (AliSelectorRL::Process(entry) == kFALSE)
4dd2ad81 74 return kFALSE;
75
0bd1f8a0 76 AliStack* stack = GetStack();
77 if (!stack)
dc740de4 78 {
0bd1f8a0 79 AliDebug(AliLog::kError, "Stack not available");
4dd2ad81 80 return kFALSE;
dc740de4 81 }
82
83 AliHeader* header = GetHeader();
84 if (!header)
85 {
86 AliDebug(AliLog::kError, "Header not available");
87 return kFALSE;
88 }
4dd2ad81 89
90 // get the MC vertex
dc740de4 91 AliGenEventHeader* genHeader = header->GenEventHeader();
4dd2ad81 92
93 TArrayF vtxMC(3);
94 genHeader->PrimaryVertex(vtxMC);
95
0bd1f8a0 96 // loop over mc particles
97 Int_t nPrim = stack->GetNprimary();
4dd2ad81 98
0bd1f8a0 99 for (Int_t iMc = 0; iMc < nPrim; ++iMc)
4dd2ad81 100 {
0bd1f8a0 101 TParticle* particle = stack->Particle(iMc);
4dd2ad81 102
103 if (!particle)
104 continue;
105
45e97e28 106 if (AliPWG0Helper::IsPrimaryCharged(particle, nPrim) == kFALSE)
4dd2ad81 107 continue;
108
0bd1f8a0 109 AliDebug(AliLog::kDebug+1, Form("Accepted primary %d, unique ID: %d", iMc, particle->GetUniqueID()));
dc740de4 110
45e97e28 111 fdNdEtaAnalysis->FillTrack(vtxMC[2], particle->Eta(), particle->Pt(), 1);
dc740de4 112 fVertex->Fill(particle->Vx(), particle->Vy(), particle->Vz());
5af55649 113
114 fPartEta->Fill(particle->Eta());
4c351225 115
116 if (TMath::Abs(particle->Eta()) < 0.8)
117 fPartPt->Fill(particle->Pt());
4dd2ad81 118 }
45e97e28 119 fdNdEtaAnalysis->FillEvent(vtxMC[2], 1);
4dd2ad81 120
5af55649 121 ++fEvents;
122
4dd2ad81 123 return kTRUE;
124}
dc740de4 125
16e24ca3 126void AlidNdEtaAnalysisMCSelector::SlaveTerminate()
127{
128 // The SlaveTerminate() function is called after all entries or objects
129 // have been processed. When running with PROOF SlaveTerminate() is called
130 // on each slave server.
131
132 AliSelectorRL::SlaveTerminate();
133
134 // Add the histograms to the output on each slave server
135 if (!fOutput)
136 {
137 AliDebug(AliLog::kError, Form("ERROR: Output list not initialized."));
138 return;
139 }
140
141 fOutput->Add(fdNdEtaAnalysis);
4c351225 142 fOutput->Add(fPartPt);
16e24ca3 143}
144
dc740de4 145void AlidNdEtaAnalysisMCSelector::Terminate()
146{
16e24ca3 147 //
148
149 AliSelectorRL::Terminate();
150
151 fdNdEtaAnalysis = dynamic_cast<dNdEtaAnalysis*> (fOutput->FindObject("dndeta"));
4c351225 152 fPartPt = dynamic_cast<TH1F*> (fOutput->FindObject("dndeta_check_pt"));
16e24ca3 153
4c351225 154 if (!fdNdEtaAnalysis || !fPartPt)
16e24ca3 155 {
4c351225 156 AliDebug(AliLog::kError, Form("ERROR: Histograms not available %p %p", (void*) fdNdEtaAnalysis, (void*) fPartPt));
16e24ca3 157 return;
158 }
159
d09fb536 160 fdNdEtaAnalysis->Finish(0, -1);
161
16e24ca3 162 TFile* fout = new TFile("analysis_mc.root","RECREATE");
163
164 fdNdEtaAnalysis->SaveHistograms();
165
166 fout->Write();
167 fout->Close();
dc740de4 168
0bd1f8a0 169 if (fPartEta)
170 {
171 fPartEta->Scale(1.0/fEvents);
172 fPartEta->Scale(1.0/fPartEta->GetBinWidth(1));
5af55649 173
0bd1f8a0 174 TCanvas* canvas = new TCanvas("control", "control", 900, 450);
175 canvas->Divide(2, 1);
5af55649 176
0bd1f8a0 177 canvas->cd(1);
178 fVertex->Draw();
5af55649 179
0bd1f8a0 180 canvas->cd(2);
181 fPartEta->Draw();
182 }
4c351225 183
184 if (fPartPt)
185 {
186 fPartPt->Scale(1.0/fEvents);
187 fPartPt->Scale(1.0/fPartPt->GetBinWidth(1));
188
189 new TCanvas("control2", "control2", 500, 500);
190 fPartPt->Draw();
191 }
dc740de4 192}