]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/dNdEta/AlidNdEtaAnalysisMCSelector.cxx
o) adding log tags to all files
[u/mrichter/AliRoot.git] / PWG0 / dNdEta / AlidNdEtaAnalysisMCSelector.cxx
1 /* $Id$ */
2
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>
11 #include <TH3F.h>
12
13 #include <AliLog.h>
14 #include <AliGenEventHeader.h>
15
16 #include "dNdEtaAnalysis.h"
17
18
19 ClassImp(AlidNdEtaAnalysisMCSelector)
20
21 AlidNdEtaAnalysisMCSelector::AlidNdEtaAnalysisMCSelector() :
22   AlidNdEtaAnalysisSelector(),
23   fVertex(0)
24 {
25   //
26   // Constructor. Initialization of pointers
27   //
28 }
29
30 AlidNdEtaAnalysisMCSelector::~AlidNdEtaAnalysisMCSelector()
31 {
32   //
33   // Destructor
34   //
35 }
36
37 void AlidNdEtaAnalysisMCSelector::Init(TTree *tree)
38 {
39    AlidNdEtaAnalysisSelector::Init(tree);
40
41   tree->SetBranchStatus("ESD", 0);
42
43   fVertex = new TH3F("vertex", "vertex", 50, -50, 50, 50, -50, 50, 50, -50, 50);
44 }
45
46 Bool_t AlidNdEtaAnalysisMCSelector::Process(Long64_t entry)
47 {
48   //
49
50   if (AliSelector::Process(entry) == kFALSE)
51     return kFALSE;
52
53   TTree* particleTree = GetKinematics();
54   if (!particleTree)
55   {
56     AliDebug(AliLog::kError, "Kinematics not available");
57     return kFALSE;
58   }
59
60   AliHeader* header = GetHeader();
61   if (!header)
62   {
63     AliDebug(AliLog::kError, "Header not available");
64     return kFALSE;
65   }
66
67   // get the MC vertex
68   AliGenEventHeader* genHeader = header->GenEventHeader();
69
70   TArrayF vtxMC(3);
71   genHeader->PrimaryVertex(vtxMC);
72
73   particleTree->SetBranchStatus("*", 0);
74   particleTree->SetBranchStatus("fDaughter[2]", 1);
75   particleTree->SetBranchStatus("fPdgCode", 1);
76   particleTree->SetBranchStatus("fPx", 1);
77   particleTree->SetBranchStatus("fPy", 1);
78   particleTree->SetBranchStatus("fPz", 1);
79   particleTree->SetBranchStatus("fVx", 1);
80   particleTree->SetBranchStatus("fVy", 1);
81   particleTree->SetBranchStatus("fVz", 1);
82
83   TParticle* particle = 0;
84   particleTree->SetBranchAddress("Particles", &particle);
85
86   Int_t nPrim  = header->GetNprimary();
87   Int_t nTotal = header->GetNtrack();
88
89   for (Int_t i_mc = nTotal - nPrim; i_mc < nTotal; ++i_mc)
90   {
91     particleTree->GetEntry(i_mc);
92
93     if (!particle)
94       continue;
95
96     if (IsPrimaryCharged(particle, nPrim) == kFALSE)
97       continue;
98
99     AliDebug(AliLog::kDebug+1, Form("Accepted primary %d, unique ID: %d", i_mc, particle->GetUniqueID()));
100
101     fdNdEtaAnalysis->FillTrack(vtxMC[2], particle->Eta(), 1);
102     fVertex->Fill(particle->Vx(), particle->Vy(), particle->Vz());
103   }
104   fdNdEtaAnalysis->FillEvent(vtxMC[2]);
105
106   return kTRUE;
107 }
108
109 void AlidNdEtaAnalysisMCSelector::Terminate()
110 {
111   AlidNdEtaAnalysisSelector::Terminate();
112
113   new TCanvas;
114   fVertex->Draw();
115 }