]>
Commit | Line | Data |
---|---|---|
7c3e5162 | 1 | |
2 | #include <TROOT.h> | |
3 | #include <TSystem.h> | |
4 | #include <TInterpreter.h> | |
5 | #include <TChain.h> | |
6 | #include <TFile.h> | |
7 | #include <TList.h> | |
8 | #include <iostream> | |
bb8a464f | 9 | #include "TH1F.h" |
7c3e5162 | 10 | #include "TH2F.h" |
11 | #include "AliFMDAnalysisTaskDndeta.h" | |
12 | #include "AliAnalysisManager.h" | |
13 | #include "AliESDFMD.h" | |
14 | #include "AliESDEvent.h" | |
15 | #include "AliAODEvent.h" | |
16 | #include "AliAODHandler.h" | |
17 | #include "AliMCEventHandler.h" | |
18 | #include "AliStack.h" | |
19 | #include "AliLog.h" | |
20 | #include "AliESDVertex.h" | |
21 | #include "TMath.h" | |
22 | #include "AliFMDAnaParameters.h" | |
23 | #include "AliFMDGeometry.h" | |
bb8a464f | 24 | #include "AliGenEventHeader.h" |
7c3e5162 | 25 | |
26 | ClassImp(AliFMDAnalysisTaskDndeta) | |
27 | ||
28 | ||
29 | AliFMDAnalysisTaskDndeta::AliFMDAnalysisTaskDndeta() | |
30 | : fDebug(0), | |
31 | fOutputList(0), | |
32 | fInputList(0), | |
33 | fArray(0), | |
34 | fInputArray(0), | |
35 | fVertexString(0x0), | |
36 | fNevents(), | |
bb8a464f | 37 | fNMCevents(), |
38 | fStandalone(kTRUE), | |
39 | fMCevent(0) | |
7c3e5162 | 40 | { |
41 | // Default constructor | |
42 | DefineInput (0, TList::Class()); | |
43 | DefineOutput(0, TList::Class()); | |
44 | } | |
45 | //_____________________________________________________________________ | |
46 | AliFMDAnalysisTaskDndeta::AliFMDAnalysisTaskDndeta(const char* name, Bool_t SE): | |
47 | AliAnalysisTask(name, "Density"), | |
48 | fDebug(0), | |
49 | fOutputList(0), | |
50 | fInputList(0), | |
51 | fArray(), | |
52 | fInputArray(0), | |
53 | fVertexString(0x0), | |
54 | fNevents(), | |
bb8a464f | 55 | fNMCevents(), |
56 | fStandalone(kTRUE), | |
57 | fMCevent(0) | |
7c3e5162 | 58 | { |
59 | fStandalone = SE; | |
60 | if(fStandalone) { | |
61 | DefineInput (0, TList::Class()); | |
62 | DefineInput(1, TObjString::Class()); | |
63 | DefineOutput(0, TList::Class()); | |
64 | ||
65 | } | |
66 | } | |
67 | //_____________________________________________________________________ | |
68 | void AliFMDAnalysisTaskDndeta::CreateOutputObjects() | |
69 | { | |
70 | AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance(); | |
71 | ||
72 | fArray.SetName("FMD"); | |
73 | fArray.SetOwner(); | |
74 | ||
75 | if(!fOutputList) | |
76 | fOutputList = new TList(); | |
77 | fOutputList->SetName("BackgroundCorrected"); | |
78 | ||
79 | ||
80 | TH2F* hMult = 0; | |
bb8a464f | 81 | TH1F* hPrimVertexBin = 0; |
7c3e5162 | 82 | |
bb8a464f | 83 | |
84 | TH2F* hBg = pars->GetBackgroundCorrection(1, 'I', 0); | |
85 | TH1F* hPrimary = new TH1F("hMultvsEta","hMultvsEta", | |
86 | hBg->GetNbinsX(), | |
87 | hBg->GetXaxis()->GetXmin(), | |
88 | hBg->GetXaxis()->GetXmax()); | |
89 | hPrimary->Sumw2(); | |
90 | fOutputList->Add(hPrimary); | |
7c3e5162 | 91 | Int_t nVtxbins = pars->GetNvtxBins(); |
92 | ||
bb8a464f | 93 | |
7c3e5162 | 94 | for(Int_t det =1; det<=3;det++) |
95 | { | |
96 | TObjArray* detArray = new TObjArray(); | |
97 | detArray->SetName(Form("FMD%d",det)); | |
98 | fArray.AddAtAndExpand(detArray,det); | |
99 | Int_t nRings = (det==1 ? 1 : 2); | |
100 | for(Int_t ring = 0;ring<nRings;ring++) | |
101 | { | |
102 | Char_t ringChar = (ring == 0 ? 'I' : 'O'); | |
103 | Int_t nSec = (ring == 0 ? 20 : 40); | |
104 | ||
105 | TObjArray* vtxArray = new TObjArray(); | |
106 | vtxArray->SetName(Form("FMD%d%c",det,ringChar)); | |
107 | detArray->AddAtAndExpand(vtxArray,ring); | |
108 | for(Int_t i = 0; i< nVtxbins; i++) { | |
109 | TH2F* hBg = pars->GetBackgroundCorrection(det, ringChar, i); | |
110 | hMult = new TH2F(Form("dNdeta_FMD%d%c_vtxbin%d",det,ringChar,i),Form("dNdeta_FMD%d%c_vtxbin%d",det,ringChar,i), | |
111 | hBg->GetNbinsX(), | |
112 | hBg->GetXaxis()->GetXmin(), | |
113 | hBg->GetXaxis()->GetXmax(), | |
114 | nSec, 0, 2*TMath::Pi()); | |
115 | hMult->Sumw2(); | |
116 | fOutputList->Add(hMult); | |
117 | vtxArray->AddAtAndExpand(hMult,i); | |
118 | ||
119 | } | |
120 | } | |
121 | } | |
122 | ||
bb8a464f | 123 | for(Int_t i = 0; i< nVtxbins; i++) { |
124 | ||
125 | hPrimVertexBin = new TH1F(Form("primmult_vtxbin%d",i), | |
126 | Form("primmult_vtxbin%d",i), | |
127 | hBg->GetNbinsX(), | |
128 | hBg->GetXaxis()->GetXmin(), | |
129 | hBg->GetXaxis()->GetXmax()); | |
130 | hPrimVertexBin->Sumw2(); | |
131 | fOutputList->Add(hPrimVertexBin); | |
132 | ||
133 | } | |
134 | ||
7c3e5162 | 135 | fNevents.SetBins(nVtxbins,0,nVtxbins); |
136 | fNevents.SetName("nEvents"); | |
bb8a464f | 137 | fNMCevents.SetBins(nVtxbins,0,nVtxbins); |
138 | fNMCevents.SetName("nMCEvents"); | |
7c3e5162 | 139 | fOutputList->Add(&fNevents); |
bb8a464f | 140 | fOutputList->Add(&fNMCevents); |
7c3e5162 | 141 | |
142 | } | |
143 | //_____________________________________________________________________ | |
144 | void AliFMDAnalysisTaskDndeta::ConnectInputData(Option_t */*option*/) | |
145 | { | |
146 | if(fStandalone) { | |
147 | fInputList = (TList*)GetInputData(0); | |
148 | fVertexString = (TObjString*)GetInputData(1); | |
149 | } | |
150 | } | |
151 | //_____________________________________________________________________ | |
152 | void AliFMDAnalysisTaskDndeta::Exec(Option_t */*option*/) | |
153 | { | |
154 | Int_t vtxbin = fVertexString->GetString().Atoi(); | |
155 | fNevents.Fill(vtxbin); | |
156 | ||
157 | for(UShort_t det=1;det<=3;det++) { | |
158 | //TObjArray* detInputArray = (TObjArray*)fInputArray->At(det); | |
159 | TObjArray* detArray = (TObjArray*)fArray.At(det); | |
160 | Int_t nRings = (det==1 ? 1 : 2); | |
161 | for (UShort_t ir = 0; ir < nRings; ir++) { | |
162 | Char_t ringChar = (ir == 0 ? 'I' : 'O'); | |
163 | //TObjArray* vtxInputArray = (TObjArray*)detInputArray->At(ir); | |
164 | TObjArray* vtxArray = (TObjArray*)detArray->At(ir); | |
165 | TH2F* hMultTotal = (TH2F*)vtxArray->At(vtxbin); | |
166 | ||
167 | ||
168 | TH2F* hMultInput = (TH2F*)fInputList->FindObject(Form("mult_FMD%d%c_vtxbin%d",det,ringChar,vtxbin)); | |
1282ce49 | 169 | |
7c3e5162 | 170 | hMultTotal->Add(hMultInput); |
171 | ||
172 | ||
173 | } | |
174 | } | |
bb8a464f | 175 | |
176 | if(fMCevent) | |
177 | ProcessPrimary(); | |
178 | ||
7c3e5162 | 179 | if(fStandalone) { |
180 | PostData(0, fOutputList); | |
181 | } | |
182 | ||
183 | } | |
184 | //_____________________________________________________________________ | |
185 | void AliFMDAnalysisTaskDndeta::Terminate(Option_t */*option*/) { | |
186 | ||
1282ce49 | 187 | |
7c3e5162 | 188 | AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance(); |
189 | ||
190 | Int_t nVtxbins = pars->GetNvtxBins(); | |
191 | ||
192 | for(UShort_t det=1;det<=3;det++) { | |
193 | TObjArray* detArray = (TObjArray*)fArray.At(det); | |
194 | Int_t nRings = (det==1 ? 1 : 2); | |
195 | for (UShort_t ir = 0; ir < nRings; ir++) { | |
196 | TObjArray* vtxArray = (TObjArray*)detArray->At(ir); | |
1282ce49 | 197 | Char_t ringChar = (ir == 0 ? 'I' : 'O'); |
7c3e5162 | 198 | for(Int_t i =0; i<nVtxbins; i++) { |
199 | TH2F* hMultTotal = (TH2F*)vtxArray->At(i); | |
1282ce49 | 200 | TH1D* hMultProj = hMultTotal->ProjectionX(Form("dNdeta_FMD%d%c_vtxbin%d_proj",det,ringChar,i),1,hMultTotal->GetNbinsY()); |
201 | fOutputList->Add(hMultProj); | |
7c3e5162 | 202 | } |
203 | } | |
204 | } | |
1282ce49 | 205 | |
bb8a464f | 206 | } |
207 | //_____________________________________________________________________ | |
208 | void AliFMDAnalysisTaskDndeta::ProcessPrimary() { | |
209 | ||
210 | AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance(); | |
211 | ||
212 | AliMCParticle* particle = 0; | |
213 | AliStack* stack = fMCevent->Stack(); | |
214 | ||
215 | TH1F* hPrimary = (TH1F*)fOutputList->FindObject("hMultvsEta"); | |
216 | ||
217 | Bool_t firstTrack = kTRUE; | |
218 | Int_t nTracks = fMCevent->GetNumberOfTracks(); | |
219 | for(Int_t i = 0 ;i<nTracks;i++) { | |
220 | particle = fMCevent->GetTrack(i); | |
221 | if(!particle) | |
222 | continue; | |
223 | if(TMath::Abs(particle->Zv()) > pars->GetVtxCutZ()) | |
224 | continue; | |
225 | if(stack->IsPhysicalPrimary(i) && particle->Charge() != 0) { | |
226 | hPrimary->Fill(particle->Eta()); | |
227 | Double_t delta = 2*pars->GetVtxCutZ()/pars->GetNvtxBins(); | |
228 | Double_t vertexBinDouble = (particle->Zv() + pars->GetVtxCutZ()) / delta; | |
229 | Int_t vertexBin = (Int_t)vertexBinDouble; | |
230 | TH1F* hPrimVtxBin = (TH1F*)fOutputList->FindObject(Form("primmult_vtxbin%d",vertexBin)); | |
231 | hPrimVtxBin->Fill(particle->Eta()); | |
232 | if(firstTrack) { | |
233 | fNMCevents.Fill(vertexBin); | |
234 | firstTrack = kFALSE; | |
235 | } | |
236 | } | |
237 | ||
238 | } | |
239 | ||
240 | ||
241 | ||
242 | ||
243 | ||
7c3e5162 | 244 | } |
245 | //_____________________________________________________________________ | |
246 | // | |
247 | // | |
248 | // EOF |