]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGDQ/dielectron/AliAnalysisTaskMultiDielectron.cxx
o update dielectron package
[u/mrichter/AliRoot.git] / PWGDQ / dielectron / AliAnalysisTaskMultiDielectron.cxx
1 /*************************************************************************
2 * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
3 *                                                                        *
4 * Author: The ALICE Off-line Project.                                    *
5 * Contributors are mentioned in the code where appropriate.              *
6 *                                                                        *
7 * Permission to use, copy, modify and distribute this software and its   *
8 * documentation strictly for non-commercial purposes is hereby granted   *
9 * without fee, provided that the above copyright notice appears in all   *
10 * copies and that both the copyright notice and this permission notice   *
11 * appear in the supporting documentation. The authors make no claims     *
12 * about the suitability of this software for any purpose. It is          *
13 * provided "as is" without express or implied warranty.                  *
14 **************************************************************************/
15
16 ///////////////////////////////////////////////////////////////////////////
17 //                                                                       //
18 //                        Basic Analysis Task                            //
19 //                                                                       //
20 ///////////////////////////////////////////////////////////////////////////
21
22 #include <TChain.h>
23 #include <TH1D.h>
24
25 #include <AliCFContainer.h>
26 #include <AliInputEventHandler.h>
27 #include <AliESDInputHandler.h>
28 #include <AliAODInputHandler.h>
29 #include <AliAnalysisManager.h>
30 #include <AliVEvent.h>
31 #include <AliTriggerAnalysis.h>
32 #include <AliPIDResponse.h>
33 #include <AliTPCPIDResponse.h>
34
35 #include "AliDielectron.h"
36 #include "AliDielectronHistos.h"
37 #include "AliDielectronCF.h"
38 #include "AliDielectronMC.h"
39 #include "AliDielectronMixingHandler.h"
40 #include "AliAnalysisTaskMultiDielectron.h"
41
42 ClassImp(AliAnalysisTaskMultiDielectron)
43
44 //_________________________________________________________________________________
45 AliAnalysisTaskMultiDielectron::AliAnalysisTaskMultiDielectron() :
46   AliAnalysisTaskSE(),
47   fListDielectron(),
48   fListHistos(),
49   fListCF(),
50   fSelectPhysics(kFALSE),
51   fTriggerMask(AliVEvent::kMB),
52   fExcludeTriggerMask(0),
53   fTriggerOnV0AND(kFALSE),
54   fRejectPileup(kFALSE),
55   fTriggerLogic(kAny),
56   fTriggerAnalysis(0x0),
57   fEventFilter(0x0),
58   fEventStat(0x0)
59 {
60   //
61   // Constructor
62   //
63 }
64
65 //_________________________________________________________________________________
66 AliAnalysisTaskMultiDielectron::AliAnalysisTaskMultiDielectron(const char *name) :
67   AliAnalysisTaskSE(name),
68   fListDielectron(),
69   fListHistos(),
70   fListCF(),
71   fSelectPhysics(kFALSE),
72   fTriggerMask(AliVEvent::kMB),
73   fExcludeTriggerMask(0),
74   fTriggerOnV0AND(kFALSE),
75   fRejectPileup(kFALSE),
76   fTriggerLogic(kAny),
77   fTriggerAnalysis(0x0),
78   fEventFilter(0x0),
79   fEventStat(0x0)
80 {
81   //
82   // Constructor
83   //
84   DefineInput(0,TChain::Class());
85   DefineOutput(1, TList::Class());
86   DefineOutput(2, TList::Class());
87   DefineOutput(3, TH1D::Class());
88   fListHistos.SetName("Dielectron_Histos_Multi");
89   fListCF.SetName("Dielectron_CF_Multi");
90   fListDielectron.SetOwner();
91   fListHistos.SetOwner();
92   fListCF.SetOwner();
93 }
94
95 //_________________________________________________________________________________
96 AliAnalysisTaskMultiDielectron::~AliAnalysisTaskMultiDielectron()
97 {
98   //
99   // Destructor
100   //
101
102   //histograms and CF are owned by the dielectron framework.
103   //however they are streamed to file, so in the first place the
104   //lists need to be owner...
105   fListHistos.SetOwner(kFALSE);
106   fListCF.SetOwner(kFALSE);
107   
108 }
109 //_________________________________________________________________________________
110 void AliAnalysisTaskMultiDielectron::UserCreateOutputObjects()
111 {
112   //
113   // Add all histogram manager histogram lists to the output TList
114   //
115
116   if (!fListHistos.IsEmpty()||!fListCF.IsEmpty()) return; //already initialised
117
118 //   AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
119 //   Bool_t isESD=man->GetInputEventHandler()->IsA()==AliESDInputHandler::Class();
120 //   Bool_t isAOD=man->GetInputEventHandler()->IsA()==AliAODInputHandler::Class();
121   
122   TIter nextDie(&fListDielectron);
123   AliDielectron *die=0;
124   while ( (die=static_cast<AliDielectron*>(nextDie())) ){
125     die->Init();
126     if (die->GetHistogramList()) fListHistos.Add(const_cast<THashList*>(die->GetHistogramList()));
127     if (die->GetCFManagerPair()) fListCF.Add(const_cast<AliCFContainer*>(die->GetCFManagerPair()->GetContainer()));
128   }
129
130   Int_t cuts=fListDielectron.GetEntries();
131   Int_t nbins=kNbinsEvent+2*cuts;
132   if (!fEventStat){
133     fEventStat=new TH1D("hEventStat","Event statistics",nbins,0,nbins);
134     fEventStat->GetXaxis()->SetBinLabel(1,"Before Phys. Sel.");
135     fEventStat->GetXaxis()->SetBinLabel(2,"After Phys. Sel.");
136
137     //default names
138     fEventStat->GetXaxis()->SetBinLabel(3,"Bin3 not used");
139     fEventStat->GetXaxis()->SetBinLabel(4,"Bin4 not used");
140     fEventStat->GetXaxis()->SetBinLabel(5,"Bin5 not used");
141     
142     if(fTriggerOnV0AND) fEventStat->GetXaxis()->SetBinLabel(3,"V0and triggers");
143     if (fEventFilter) fEventStat->GetXaxis()->SetBinLabel(4,"After Event Filter");
144     if (fRejectPileup) fEventStat->GetXaxis()->SetBinLabel(5,"After Pileup rejection");
145     
146     for (Int_t i=0; i<cuts; ++i){
147       fEventStat->GetXaxis()->SetBinLabel((kNbinsEvent+1)+2*i,Form("#splitline{1 candidate}{%s}",fListDielectron.At(i)->GetName()));
148       fEventStat->GetXaxis()->SetBinLabel((kNbinsEvent+2)+2*i,Form("#splitline{With >1 candidate}{%s}",fListDielectron.At(i)->GetName()));
149     }
150   }
151
152   if (!fTriggerAnalysis) fTriggerAnalysis=new AliTriggerAnalysis;
153   fTriggerAnalysis->EnableHistograms();
154   fTriggerAnalysis->SetAnalyzeMC(AliDielectronMC::Instance()->HasMC());
155   
156   PostData(1, &fListHistos);
157   PostData(2, &fListCF);
158   PostData(3, fEventStat);
159 }
160
161 //_________________________________________________________________________________
162 void AliAnalysisTaskMultiDielectron::UserExec(Option_t *)
163 {
164   //
165   // Main loop. Called for every event
166   //
167
168   if (fListHistos.IsEmpty()&&fListCF.IsEmpty()) return;
169
170   AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
171   Bool_t isESD=man->GetInputEventHandler()->IsA()==AliESDInputHandler::Class();
172   Bool_t isAOD=man->GetInputEventHandler()->IsA()==AliAODInputHandler::Class();
173   
174   AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
175   if (!inputHandler) return;
176   
177 //   AliPIDResponse *pidRes=inputHandler->GetPIDResponse();
178   if ( inputHandler->GetPIDResponse() ){
179     // for the 2.76 pass2 MC private train. Together with a sigma shift of -0.169
180 //    pidRes->GetTPCResponse().SetSigma(4.637e-3,2.41332105409873257e+04);
181     AliDielectronVarManager::SetPIDResponse( inputHandler->GetPIDResponse() );
182   } else {
183     AliFatal("This task needs the PID response attached to the input event handler!");
184   }
185   
186   // Was event selected ?
187   ULong64_t isSelected = AliVEvent::kAny;
188   Bool_t isRejected = kFALSE;
189   if( fSelectPhysics && inputHandler){
190     if((isESD && inputHandler->GetEventSelection()) || isAOD){
191       isSelected = inputHandler->IsEventSelected();
192       if (fExcludeTriggerMask && (isSelected&fExcludeTriggerMask)) isRejected=kTRUE;
193       if (fTriggerLogic==kAny) isSelected&=fTriggerMask;
194       else if (fTriggerLogic==kExact) isSelected=((isSelected&fTriggerMask)==fTriggerMask);
195     }
196    }
197  
198  
199   //Before physics selection
200   fEventStat->Fill(kAllEvents);
201   if (isSelected==0||isRejected) {
202     PostData(3,fEventStat);
203     return;
204   }
205   //after physics selection
206   fEventStat->Fill(kSelectedEvents);
207
208   //V0and
209   if(fTriggerOnV0AND){
210   if(isESD){if (!fTriggerAnalysis->IsOfflineTriggerFired(static_cast<AliESDEvent*>(InputEvent()), AliTriggerAnalysis::kV0AND))
211             return;}
212   if(isAOD){if(!((static_cast<AliAODEvent*>(InputEvent()))->GetVZEROData()->GetV0ADecision() == AliVVZERO::kV0BB &&
213             (static_cast<AliAODEvent*>(InputEvent()))->GetVZEROData()->GetV0CDecision() == AliVVZERO::kV0BB) )
214             return;}
215    }
216   
217
218   fEventStat->Fill(kV0andEvents);
219
220   //Fill Event histograms before the event filter
221   TIter nextDie(&fListDielectron);
222   AliDielectron *die=0;
223   Double_t values[AliDielectronVarManager::kNMaxValues]={0};
224   Double_t valuesMC[AliDielectronVarManager::kNMaxValues]={0};
225   AliDielectronVarManager::Fill(InputEvent(),values);
226   AliDielectronVarManager::Fill(InputEvent(),valuesMC);
227   Bool_t hasMC=AliDielectronMC::Instance()->HasMC();
228   if (hasMC) {
229     if (AliDielectronMC::Instance()->ConnectMCEvent())
230       AliDielectronVarManager::Fill(AliDielectronMC::Instance()->GetMCEvent(),valuesMC);
231   }
232
233   while ( (die=static_cast<AliDielectron*>(nextDie())) ){
234     AliDielectronHistos *h=die->GetHistoManager();
235     if (h){
236       if (h->GetHistogramList()->FindObject("Event_noCuts"))
237         h->FillClass("Event_noCuts",AliDielectronVarManager::kNMaxValues,values);
238       if (hasMC && h->GetHistogramList()->FindObject("MCEvent_noCuts"))
239         h->FillClass("Event_noCuts",AliDielectronVarManager::kNMaxValues,valuesMC);
240     }
241   }
242   nextDie.Reset();
243   
244   //event filter
245   if (fEventFilter) {
246     if (!fEventFilter->IsSelected(InputEvent())) return;
247   }
248   fEventStat->Fill(kFilteredEvents);
249   
250   //pileup
251   if (fRejectPileup){
252     if (InputEvent()->IsPileupFromSPD(3,0.8,3.,2.,5.)) return;
253   }
254   fEventStat->Fill(kPileupEvents);
255   
256   //bz for AliKF
257   Double_t bz = InputEvent()->GetMagneticField();
258   AliKFParticle::SetField( bz );
259
260   AliDielectronPID::SetCorrVal((Double_t)InputEvent()->GetRunNumber());
261   
262   //Process event in all AliDielectron instances
263   //   TIter nextDie(&fListDielectron);
264   //   AliDielectron *die=0;
265   Int_t idie=0;
266   while ( (die=static_cast<AliDielectron*>(nextDie())) ){
267     die->Process(InputEvent());
268     if (die->HasCandidates()){
269       Int_t ncandidates=die->GetPairArray(1)->GetEntriesFast();
270       if (ncandidates==1) fEventStat->Fill((kNbinsEvent)+2*idie);
271       else if (ncandidates>1) fEventStat->Fill((kNbinsEvent+1)+2*idie);
272     }
273     ++idie;
274   }
275   
276   PostData(1, &fListHistos);
277   PostData(2, &fListCF);
278   PostData(3,fEventStat);
279 }
280
281 //_________________________________________________________________________________
282 void AliAnalysisTaskMultiDielectron::FinishTaskOutput()
283 {
284   //
285   // Write debug tree
286   //
287   TIter nextDie(&fListDielectron);
288   AliDielectron *die=0;
289   while ( (die=static_cast<AliDielectron*>(nextDie())) ){
290     die->SaveDebugTree();
291     AliDielectronMixingHandler *mix=die->GetMixingHandler();
292 //    printf("\n\n\n===============\ncall mix in Terminate: %p (%p)\n=================\n\n",mix,die);
293     if (mix) mix->MixRemaining(die);
294   }
295   PostData(1, &fListHistos);
296   PostData(2, &fListCF);
297 }
298