]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/AliSPDMCTrackDensity.cxx
A better way to specify the Nch axis for the MultDists task.
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / AliSPDMCTrackDensity.cxx
1 #include "AliSPDMCTrackDensity.h"
2 #include "AliMCEvent.h"
3 #include "AliTrackReference.h"
4 #include <TMath.h>
5 #include <AliLog.h>
6 #include <TROOT.h>
7 #include <TH2D.h>
8 #include <iostream>
9
10 //____________________________________________________________________
11 AliSPDMCTrackDensity::AliSPDMCTrackDensity()
12   : AliBaseMCTrackDensity(), 
13     fMinR(3.5), 
14     fMaxR(4.5),
15     fMinZ(-15), // -14.1), 
16     fMaxZ(+15), // +14.1)
17     fStored(0), 
18     fOutput(0)
19 {
20   // Default constructor 
21 }
22
23 //____________________________________________________________________
24 AliSPDMCTrackDensity::AliSPDMCTrackDensity(const char*)
25   : AliBaseMCTrackDensity("spdMCTrackDensity"), 
26     fMinR(3.5), 
27     fMaxR(4.5),
28     fMinZ(-14.1), 
29     fMaxZ(+14.1),
30     fStored(0), 
31     fOutput(0)
32 {
33   // Normal constructor constructor 
34 }
35
36 //____________________________________________________________________
37 AliSPDMCTrackDensity::AliSPDMCTrackDensity(const AliSPDMCTrackDensity& o)
38   : AliBaseMCTrackDensity(o),
39     fMinR(o.fMinR), 
40     fMaxR(o.fMaxR),
41     fMinZ(o.fMinZ), 
42     fMaxZ(o.fMaxZ), 
43     fStored(o.fStored), 
44     fOutput(o.fOutput)
45 {
46   // Normal constructor constructor 
47 }
48
49 //____________________________________________________________________
50 AliSPDMCTrackDensity&
51 AliSPDMCTrackDensity::operator=(const AliSPDMCTrackDensity& o)
52 {
53   // Assignment operator 
54   if (&o == this) return *this;
55   AliBaseMCTrackDensity::operator=(o);
56   fMinR             = o.fMinR;
57   fMaxR             = o.fMaxR;
58   fMinZ             = o.fMinZ;
59   fMaxZ             = o.fMaxZ;
60   fStored           = o.fStored;
61   fOutput           = o.fOutput;
62
63   return *this;
64 }
65
66 //____________________________________________________________________
67 Int_t
68 AliSPDMCTrackDensity::GetDetectorId() const
69 {
70   return AliTrackReference::kITS;
71 }
72
73
74 //____________________________________________________________________
75 void
76 AliSPDMCTrackDensity::BeginTrackRefs()
77 {
78   fStored = 0;
79 }
80
81 //____________________________________________________________________
82 Bool_t
83 AliSPDMCTrackDensity::CheckTrackRef(AliTrackReference*   ref) const
84 {
85   // Get radius and z where the track reference was made 
86   Double_t r = ref->R();
87   Double_t z = ref->Z();
88   if (r > fMaxR || r < fMinR) return false;
89   if (z > fMaxZ || z < fMinZ) return false;
90
91   return true;
92 }
93 //____________________________________________________________________
94 AliTrackReference*
95 AliSPDMCTrackDensity::ProcessRef(AliMCParticle*       /*particle*/,
96                                  const AliMCParticle* /*mother*/,
97                                  AliTrackReference*   ref)
98 {
99   if (fStored) return 0;
100
101   return fStored = ref;
102 }
103
104 //____________________________________________________________________
105 Double_t
106 AliSPDMCTrackDensity::StoreParticle(AliMCParticle* particle, 
107                                     const AliMCParticle* mother, 
108                                     AliTrackReference*   ref) const
109 {
110   Double_t w = AliBaseMCTrackDensity::StoreParticle(particle, mother, ref);
111   Double_t r = ref->R();
112   Double_t x = ref->X();
113   Double_t y = ref->Y();
114   Double_t z = ref->Z();
115
116   Double_t zr = z-fVz;
117   Double_t th = TMath::ATan2(r,zr);
118   if (th < 0) th += 2*TMath::Pi();
119   Double_t et = -TMath::Log(TMath::Tan(th/2));
120   Double_t ph = TMath::ATan2(y,x);
121   if (ph < 0) ph += 2*TMath::Pi();
122   fOutput->Fill(et,ph,w);
123
124   return w;
125 }
126
127
128 //____________________________________________________________________
129 Bool_t
130 AliSPDMCTrackDensity::Calculate(const AliMCEvent& event, 
131                                 Double_t          vz,
132                                 TH2D&             output, 
133                                 TH2D*             primary)
134 {
135   // 
136   // Filter the input kinematics and track references, using 
137   // some of the ESD information
138   // 
139   // Parameters:
140   //    input   Input ESD event
141   //    event   Input MC event
142   //    vz      Vertex position 
143   //    output  Output ESD-like object
144   //    primary Per-event histogram of primaries 
145   //
146   // Return:
147   //    True on succes, false otherwise 
148   //
149   fOutput = &output;
150
151   return ProcessTracks(event, vz, primary);
152 }
153
154 //____________________________________________________________________
155 void
156 AliSPDMCTrackDensity::Print(Option_t* option) const 
157 {
158   AliBaseMCTrackDensity::Print(option);
159   char ind[gROOT->GetDirLevel()+1];
160   for (Int_t i = 0; i < gROOT->GetDirLevel(); i++) ind[i] = ' ';
161   ind[gROOT->GetDirLevel()] = '\0';
162   std::cout << ind << " R range:                [" << fMinR << ',' << fMaxR
163             << "]\n"
164             << ind << " Z range:                [" << fMinZ << ',' << fMaxZ
165             << "]" << std::endl;
166   
167 }
168
169 //____________________________________________________________________
170 //
171 // EOF
172 //