]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/analysis/AliFMDAnalysisTaskBackgroundCorrection.cxx
Initialisations in constructor.
[u/mrichter/AliRoot.git] / FMD / analysis / AliFMDAnalysisTaskBackgroundCorrection.cxx
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>
9 #include "TH2F.h"
10 #include "AliFMDAnalysisTaskBackgroundCorrection.h"
11 #include "AliAnalysisManager.h"
12 #include "AliESDFMD.h"
13 #include "AliESDEvent.h"
14 #include "AliAODEvent.h"
15 #include "AliAODHandler.h"
16 #include "AliMCEventHandler.h"
17 #include "AliStack.h"
18 #include "AliLog.h"
19 #include "AliESDVertex.h"
20 #include "TMath.h"
21 #include "AliFMDAnaParameters.h"
22 //#include "AliFMDGeometry.h"
23
24 ClassImp(AliFMDAnalysisTaskBackgroundCorrection)
25
26
27 AliFMDAnalysisTaskBackgroundCorrection::AliFMDAnalysisTaskBackgroundCorrection()
28 : fDebug(0),
29   fOutputList(0),
30   fInputList(0),
31   fHitList(0),
32   fArray(0),
33   fInputArray(0),
34   fVertexString(0x0),
35   fNevents(),
36   fStandalone(kTRUE),
37   fOutputVertexString(0)
38 {
39   // Default constructor
40   DefineInput (0, TList::Class());
41   DefineOutput(0, TList::Class());
42   DefineOutput(1, TObjString::Class());
43 }
44 //_____________________________________________________________________
45 AliFMDAnalysisTaskBackgroundCorrection::AliFMDAnalysisTaskBackgroundCorrection(const char* name, Bool_t SE):
46     AliAnalysisTask(name, "Density"),
47     fDebug(0),
48     fOutputList(0),
49     fInputList(0),
50     fHitList(0),
51     fArray(),
52     fInputArray(0),
53     fVertexString(0x0),
54     fNevents(),
55     fStandalone(kTRUE),
56     fOutputVertexString(0)
57 {
58   fStandalone = SE;
59   if(fStandalone) {
60     DefineInput (0, TList::Class());
61     DefineOutput(0, TList::Class());
62     DefineOutput(1, TObjString::Class());
63   }
64 }
65 //_____________________________________________________________________
66 void AliFMDAnalysisTaskBackgroundCorrection::CreateOutputObjects()
67 {
68   AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
69   
70   fArray.SetName("FMD");
71   fArray.SetOwner();
72   
73   if(!fOutputList)
74     fOutputList = new TList();
75   fOutputList->SetName("BackgroundCorrectedPerEvent");
76   if(!fHitList)
77     fHitList = new TList();
78   fHitList->SetName("HitsList");
79   
80   
81   TH2F* hMult = 0;
82   TH2F* hHits = 0;
83   Int_t nVtxbins = pars->GetNvtxBins();
84   
85   for(Int_t det =1; det<=3;det++)
86     {
87       TObjArray* detArray = new TObjArray();
88       detArray->SetName(Form("FMD%d",det));
89       fArray.AddAtAndExpand(detArray,det);
90       Int_t nRings = (det==1 ? 1 : 2);
91       for(Int_t ring = 0;ring<nRings;ring++)
92         {
93           Char_t ringChar = (ring == 0 ? 'I' : 'O');
94           Int_t  nSec     = (ring == 0 ? 20 : 40);
95           
96           TObjArray* vtxArray = new TObjArray();
97           vtxArray->SetName(Form("FMD%d%c",det,ringChar));
98           detArray->AddAtAndExpand(vtxArray,ring);
99           for(Int_t i = 0; i< nVtxbins; i++) {
100             TH2F* hBg = pars->GetBackgroundCorrection(det, ringChar, i);
101             hMult  = new TH2F(Form("mult_FMD%d%c_vtxbin%d",det,ringChar,i),Form("mult_FMD%d%c_vtxbin%d",det,ringChar,i),
102                               hBg->GetNbinsX(),
103                               hBg->GetXaxis()->GetXmin(),
104                               hBg->GetXaxis()->GetXmax(),
105                               nSec, 0, 2*TMath::Pi());
106             hMult->Sumw2();
107             fOutputList->Add(hMult);
108             hHits  = new TH2F(Form("hits_FMD%d%c_vtxbin%d",det,ringChar,i),Form("hits_FMD%d%c_vtxbin%d",det,ringChar,i),
109                               hBg->GetNbinsX(),
110                               hBg->GetXaxis()->GetXmin(),
111                               hBg->GetXaxis()->GetXmax(),
112                               nSec, 0, 2*TMath::Pi());
113             
114             hHits->Sumw2();
115             fHitList->Add(hHits);
116             fOutputList->Add(hHits);
117             vtxArray->AddAtAndExpand(hMult,i);
118             
119           }
120         } 
121     }
122   if(fStandalone) {
123     fOutputVertexString = new TObjString();
124   }
125   fOutputList->Add(fOutputVertexString);
126   
127   
128   
129 }
130 //_____________________________________________________________________
131 void AliFMDAnalysisTaskBackgroundCorrection::ConnectInputData(Option_t */*option*/)
132 {
133   if(fStandalone) {
134     fInputList   = (TList*)GetInputData(0);
135     
136   }
137 }
138 //_____________________________________________________________________
139 void AliFMDAnalysisTaskBackgroundCorrection::Exec(Option_t */*option*/)
140 {
141   AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
142   
143   fInputArray   = (TObjArray*)fInputList->At(0);
144   fVertexString = (TObjString*)fInputList->At(1);
145   
146   
147   Int_t vtxbin   = fVertexString->GetString().Atoi();
148   fOutputVertexString->SetString(Form("%d",vtxbin));
149   
150   //fNevents.operator[](vtxbin)++;
151   fNevents.Fill(vtxbin);
152   
153   //Reset everything
154   for(UShort_t det=1;det<=3;det++) {
155     TObjArray* detArray = (TObjArray*)fArray.At(det);
156     Int_t nRings = (det==1 ? 1 : 2);
157     for (UShort_t ir = 0; ir < nRings; ir++) {
158       TObjArray* vtxArray = (TObjArray*)detArray->At(ir);
159       
160       TH2F* hMult   = (TH2F*)vtxArray->At(vtxbin); 
161       hMult->Reset();
162     }
163     
164   }
165   
166   
167   
168   for(UShort_t det=1;det<=3;det++) {
169     TObjArray* detInputArray = (TObjArray*)fInputArray->At(det);
170     TObjArray* detArray = (TObjArray*)fArray.At(det);
171     Int_t nRings = (det==1 ? 1 : 2);
172     for (UShort_t ir = 0; ir < nRings; ir++) {
173       Char_t ringChar = (ir == 0 ? 'I' : 'O');
174       TObjArray* vtxInputArray = (TObjArray*)detInputArray->At(ir);
175       TObjArray* vtxArray = (TObjArray*)detArray->At(ir);
176       TH2F* hMultTotal = (TH2F*)vtxArray->At(vtxbin);
177       TH2F* hMultInput = (TH2F*)vtxInputArray->At(vtxbin);
178       TH2F* hHits      = (TH2F*)fOutputList->FindObject(Form("hits_FMD%d%c_vtxbin%d",det,ringChar,vtxbin));
179       
180       hHits->Add(hMultInput);
181       TH2F* hBg        = pars->GetBackgroundCorrection(det, ringChar, vtxbin);
182       
183       //TH2F* hTmp       = (TH2F*)hMultInput->Clone("hMult_from_event");
184       hMultTotal->Add(hMultInput);
185       
186       hMultTotal->Divide(hBg);//,"B");
187       /*for(Int_t i = 1; i<=hTmp->GetNbinsX();i++) {
188         for(Int_t j = 1; j<=hTmp->GetNbinsY();j++) {
189           Float_t mult = hTmp->GetBinContent(i,j);
190           if(mult == 0) continue;
191           Float_t correction = hBg->GetBinContent(i,j);
192           
193           Float_t multcor = mult;
194           if(correction != 0)
195             multcor = multcor/correction;
196           else
197             std::cout<<"Warning! No correction for bin "<<i<<" , "<<j<<std::endl;
198           
199           hTmp->SetBinContent(i,j,multcor);
200         }
201         }*/
202       
203       
204       
205       // delete hTmp;
206       
207     }
208   }
209   if(fStandalone) {
210     PostData(0, fOutputList); 
211     PostData(1, fOutputVertexString);
212   }
213   
214 }
215 //_____________________________________________________________________
216 void AliFMDAnalysisTaskBackgroundCorrection::Terminate(Option_t */*option*/) {
217   
218   /*
219   AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
220   
221   Int_t nVtxbins = pars->GetNvtxBins();
222   
223   for(UShort_t det=1;det<=3;det++) {
224     TObjArray* detArray = (TObjArray*)fArray.At(det);
225     Int_t nRings = (det==1 ? 1 : 2);
226     for (UShort_t ir = 0; ir < nRings; ir++) {
227       TObjArray* vtxArray = (TObjArray*)detArray->At(ir);
228       for(Int_t i =0; i<nVtxbins; i++) {
229         TH2F* hMultTotal = (TH2F*)vtxArray->At(i);
230         if(fNevents.At(i))
231           hMultTotal->Scale(1/(Float_t)fNevents.At(i));
232       }
233     }
234   }
235   */
236 }
237 //_____________________________________________________________________
238 //
239 //
240 // EOF