]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGDQ/dielectron/AliAnalysisTaskDielectronFilter.cxx
- coverity fix
[u/mrichter/AliRoot.git] / PWGDQ / dielectron / AliAnalysisTaskDielectronFilter.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 #include <AliLog.h>
25 #include <AliAODHandler.h>
26 #include <AliAODInputHandler.h>
27 #include <AliAnalysisManager.h>
28 #include <AliVEvent.h>
29 #include <AliTriggerAnalysis.h>
30 #include <AliInputEventHandler.h>
31 #include <AliAODInputHandler.h>
32 #include <AliESDInputHandler.h>
33
34 #include "AliDielectron.h"
35 #include "AliDielectronMC.h"
36 #include "AliDielectronHistos.h"
37 #include "AliDielectronVarManager.h"
38 #include "AliAnalysisTaskDielectronFilter.h"
39 #include "AliAODCaloCluster.h"
40
41 ClassImp(AliAnalysisTaskDielectronFilter)
42
43 //_________________________________________________________________________________
44 AliAnalysisTaskDielectronFilter::AliAnalysisTaskDielectronFilter() :
45 AliAnalysisTaskSE(),
46 fDielectron(0),
47 fSelectPhysics(kTRUE),
48 fTriggerMask(AliVEvent::kMB),
49 fExcludeTriggerMask(0),
50 fTriggerOnV0AND(kFALSE),
51 fRejectPileup(kFALSE),
52 fEventStat(0x0),
53 fTriggerLogic(kAny),
54 fTriggerAnalysis(0x0),
55 fStoreLikeSign(kFALSE),
56 fStoreRotatedPairs(kFALSE),
57 fCreateNanoAOD(kFALSE),
58 fEventFilter(0x0)
59 {
60   //
61   // Constructor
62   //
63 }
64
65 //_________________________________________________________________________________
66 AliAnalysisTaskDielectronFilter::AliAnalysisTaskDielectronFilter(const char *name) :
67 AliAnalysisTaskSE(name),
68 fDielectron(0),
69 fSelectPhysics(kTRUE),
70 fTriggerMask(AliVEvent::kMB),
71 fExcludeTriggerMask(0),
72 fTriggerOnV0AND(kFALSE),
73 fRejectPileup(kFALSE),
74 fEventStat(0x0),
75 fTriggerLogic(kAny),
76 fTriggerAnalysis(0x0),
77 fStoreLikeSign(kFALSE),
78 fStoreRotatedPairs(kFALSE),
79 fCreateNanoAOD(kFALSE),
80 fEventFilter(0x0)
81 {
82   //
83   // Constructor
84   //
85   DefineInput(0,TChain::Class());
86   DefineOutput(1, THashList::Class());
87   DefineOutput(2, TH1D::Class());
88 }
89
90 //_________________________________________________________________________________
91 void AliAnalysisTaskDielectronFilter::Init()
92 {
93   // Initialization
94   if (fDebug > 1) AliInfo("Init() \n");
95   
96 // require AOD handler
97   AliAODHandler *aodH = (AliAODHandler*)((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
98   if (!aodH) AliFatal("No AOD handler. Halting.");
99     
100   aodH->AddFilteredAOD("AliAOD.Dielectron.root", "DielectronEvents");
101 //   AddAODBranch("AliDielectronCandidates",fDielectron->GetPairArraysPointer(),"deltaAOD.Dielectron.root");
102 }
103
104 //_________________________________________________________________________________
105 void AliAnalysisTaskDielectronFilter::UserCreateOutputObjects()
106 {
107   //
108   // Initilise histograms
109   //
110
111   //require dielectron framework
112   if (!fDielectron) {
113     AliFatal("Dielectron framework class required. Please create and instance with proper cuts and set it via 'SetDielectron' before executing this task!!!");
114     return;
115   }
116   if(fStoreRotatedPairs) fDielectron->SetStoreRotatedPairs(kTRUE);
117   fDielectron->SetDontClearArrays(); 
118   fDielectron->Init();
119
120   Int_t nbins=kNbinsEvent+2;
121   if (!fEventStat){
122     fEventStat=new TH1D("hEventStat","Event statistics",nbins,0,nbins);
123     fEventStat->GetXaxis()->SetBinLabel(1,"Before Phys. Sel.");
124     fEventStat->GetXaxis()->SetBinLabel(2,"After Phys. Sel.");
125
126     //default names
127     fEventStat->GetXaxis()->SetBinLabel(3,"Bin3 not used");
128     fEventStat->GetXaxis()->SetBinLabel(4,"Bin4 not used");
129     fEventStat->GetXaxis()->SetBinLabel(5,"Bin5 not used");
130
131     if(fTriggerOnV0AND) fEventStat->GetXaxis()->SetBinLabel(3,"V0and triggers");
132     if (fEventFilter) fEventStat->GetXaxis()->SetBinLabel(4,"After Event Filter");
133     if (fRejectPileup) fEventStat->GetXaxis()->SetBinLabel(5,"After Pileup rejection");
134
135     fEventStat->GetXaxis()->SetBinLabel((kNbinsEvent+1),Form("#splitline{1 candidate}{%s}",fDielectron->GetName()));
136     fEventStat->GetXaxis()->SetBinLabel((kNbinsEvent+2),Form("#splitline{With >1 candidate}{%s}",fDielectron->GetName()));
137    }
138
139 Bool_t isAOD=AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()->IsA()==AliAODInputHandler::Class(); 
140 if(fCreateNanoAOD && isAOD){
141   AliAODHandler *aodH = (AliAODHandler*)((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
142   AliAODExtension *extDielectron = aodH->GetFilteredAOD("AliAOD.Dielectron.root");
143    TClonesArray *nanoAODTracks = new TClonesArray("AliAODTrack",500);
144    nanoAODTracks->SetName("tracks");
145    extDielectron->AddBranch("TClonesArray", &nanoAODTracks);
146    TClonesArray *nanoAODVertices = new TClonesArray("AliAODVertex",500);
147    nanoAODVertices->SetName("vertices");
148    extDielectron->AddBranch("TClonesArray", &nanoAODVertices);
149    TClonesArray *nanoAODCaloCluster = new TClonesArray("AliAODCaloCluster",500);
150    nanoAODCaloCluster->SetName("caloClusters");
151    extDielectron->AddBranch("TClonesArray", &nanoAODCaloCluster);  
152    extDielectron->GetAOD()->GetStdContent();
153    }else if(fCreateNanoAOD && !isAOD){AliWarning("Filtered-Nano AODs creation works only on AODs ");  } 
154   
155   PostData(2,fEventStat);
156 }
157
158 //_________________________________________________________________________________
159 void AliAnalysisTaskDielectronFilter::UserExec(Option_t *)
160 {
161   //
162   // Main loop. Called for every event
163   //
164   
165   if (!fDielectron) return;
166   
167   AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
168   Bool_t isESD=man->GetInputEventHandler()->IsA()==AliESDInputHandler::Class();
169   Bool_t isAOD=man->GetInputEventHandler()->IsA()==AliAODInputHandler::Class();
170
171   
172   AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
173   if (!inputHandler) return;
174
175   if ( inputHandler->GetPIDResponse() ){
176     AliDielectronVarManager::SetPIDResponse( inputHandler->GetPIDResponse() );
177   } else {
178     AliFatal("This task needs the PID response attached to the input event handler!");
179   }
180   
181   // Was event selected ?
182   ULong64_t isSelected = AliVEvent::kAny;
183   Bool_t isRejected = kFALSE;
184   if( fSelectPhysics && inputHandler){
185     if((isESD && inputHandler->GetEventSelection()) || isAOD){
186       isSelected = inputHandler->IsEventSelected();
187       if (fExcludeTriggerMask && (isSelected&fExcludeTriggerMask)) isRejected=kTRUE;
188       if (fTriggerLogic==kAny) isSelected&=fTriggerMask;
189       else if (fTriggerLogic==kExact) isSelected=((isSelected&fTriggerMask)==fTriggerMask);
190     }
191    }
192  
193   //before physics selection
194   fEventStat->Fill(kAllEvents);
195   if (isSelected==0||isRejected) {
196     PostData(3,fEventStat);
197     return;
198   }
199   //after physics selection
200   fEventStat->Fill(kSelectedEvents);
201
202   //V0and
203   if(fTriggerOnV0AND){
204   if(isESD){if (!fTriggerAnalysis->IsOfflineTriggerFired(static_cast<AliESDEvent*>(InputEvent()), AliTriggerAnalysis::kV0AND))
205             return;}
206   if(isAOD){if(!((static_cast<AliAODEvent*>(InputEvent()))->GetVZEROData()->GetV0ADecision() == AliVVZERO::kV0BB &&
207             (static_cast<AliAODEvent*>(InputEvent()))->GetVZEROData()->GetV0CDecision() == AliVVZERO::kV0BB) )
208             return;}
209    }
210
211    fEventStat->Fill(kV0andEvents);
212
213   //Fill Event histograms before the event filter
214   Double_t values[AliDielectronVarManager::kNMaxValues]={0};
215   Double_t valuesMC[AliDielectronVarManager::kNMaxValues]={0};
216   AliDielectronVarManager::SetEvent(InputEvent());
217   AliDielectronVarManager::Fill(InputEvent(),values);
218   AliDielectronVarManager::Fill(InputEvent(),valuesMC);
219
220   Bool_t hasMC=AliDielectronMC::Instance()->HasMC();
221   if (hasMC) {
222     if (AliDielectronMC::Instance()->ConnectMCEvent())
223       AliDielectronVarManager::Fill(AliDielectronMC::Instance()->GetMCEvent(),valuesMC);
224   }
225
226   AliDielectronHistos *h=fDielectron->GetHistoManager();
227     if (h){
228       if (h->GetHistogramList()->FindObject("Event_noCuts"))
229         h->FillClass("Event_noCuts",AliDielectronVarManager::kNMaxValues,values);
230       if (hasMC && h->GetHistogramList()->FindObject("MCEvent_noCuts"))
231         h->FillClass("Event_noCuts",AliDielectronVarManager::kNMaxValues,valuesMC);
232     }
233
234   //event filter
235   if (fEventFilter) {
236   if (!fEventFilter->IsSelected(InputEvent())) return;
237   }
238   fEventStat->Fill(kFilteredEvents);
239
240   //pileup
241   if (fRejectPileup){
242   if (InputEvent()->IsPileupFromSPD(3,0.8,3.,2.,5.)) return;
243   }
244   fEventStat->Fill(kPileupEvents);
245
246   //bz for AliKF
247   Double_t bz = InputEvent()->GetMagneticField();
248   AliKFParticle::SetField( bz );
249
250   AliDielectronPID::SetCorrVal((Double_t)InputEvent()->GetRunNumber());
251   
252   fDielectron->Process(InputEvent());
253   
254   Bool_t hasCand = kFALSE;
255   if(fStoreLikeSign) hasCand = (fDielectron->HasCandidates() || fDielectron->HasCandidatesLikeSign());
256   else hasCand = (fDielectron->HasCandidates());
257
258   if(fStoreRotatedPairs) hasCand = (hasCand || fDielectron->HasCandidatesTR());
259   
260   if(hasCand){
261     AliAODHandler *aodH=(AliAODHandler*)((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
262     AliAODEvent *aod = aodH->GetAOD();
263     
264     // reset bit for all tracks
265     if(isAOD){
266     for(Int_t it=0;it<aod->GetNumberOfTracks();it++){
267         aod->GetTrack(it)->ResetBit(kIsReferenced);  aod->GetTrack(it)->SetUniqueID(0);
268         }
269     }
270
271     //replace the references of the legs with the AOD references
272     TObjArray *obj = 0x0;
273     for(Int_t i=0; i < 11; i++ ){
274       obj = (TObjArray*)((*(fDielectron->GetPairArraysPointer()))->UncheckedAt(i));
275       if(!obj) continue;
276       for(int j=0;j<obj->GetEntriesFast();j++){
277         AliDielectronPair *pairObj = (AliDielectronPair*)obj->UncheckedAt(j);
278         Int_t id1 = ((AliVTrack*)pairObj->GetFirstDaughter())->GetID();
279         Int_t id2 = ((AliVTrack*)pairObj->GetSecondDaughter())->GetID();
280         
281         for(Int_t it=0;it<aod->GetNumberOfTracks();it++){
282           if(aod->GetTrack(it)->GetID() == id1) pairObj->SetRefFirstDaughter(aod->GetTrack(it)); 
283           if(aod->GetTrack(it)->GetID() == id2) pairObj->SetRefSecondDaughter(aod->GetTrack(it));
284         }
285       }
286     }
287     
288     AliAODExtension *extDielectron = aodH->GetFilteredAOD("AliAOD.Dielectron.root");
289     extDielectron->SelectEvent();
290     Int_t ncandidates=fDielectron->GetPairArray(1)->GetEntriesFast();
291     if (ncandidates==1) fEventStat->Fill((kNbinsEvent));
292     else if (ncandidates>1) fEventStat->Fill((kNbinsEvent+1));
293  
294     //see if dielectron candidate branch exists, if not create is
295     TTree *t=extDielectron->GetTree();
296
297     if(!t->GetListOfBranches()->GetEntries() && isAOD)
298       t->Branch(aod->GetList());
299     
300     if (!t->GetBranch("dielectrons"))
301       t->Bronch("dielectrons","TObjArray",fDielectron->GetPairArraysPointer());
302       
303       // store positive and negative tracks
304       if(fCreateNanoAOD && isAOD){
305       Int_t nTracks = (fDielectron->GetTrackArray(0))->GetEntries() + (fDielectron->GetTrackArray(1))->GetEntries();
306        AliAODEvent *nanoEv = extDielectron->GetAOD();
307        nanoEv->GetTracks()->Clear();      
308        nanoEv->GetVertices()->Clear();      
309        nanoEv->GetCaloClusters()->Clear();      
310        
311        AliAODVertex* tmp = ((static_cast<AliAODEvent*>(InputEvent()))->GetPrimaryVertex())->CloneWithoutRefs();
312        nanoEv->AddVertex(tmp);
313        AliAODVertex* tmpSpd = ((static_cast<AliAODEvent*>(InputEvent()))->GetPrimaryVertexSPD())->CloneWithoutRefs();
314        nanoEv->AddVertex(tmpSpd);
315        nanoEv->GetVertex(0)->SetNContributors((static_cast<AliAODEvent*>(InputEvent()))->GetPrimaryVertex()->GetNContributors());
316        nanoEv->GetVertex(1)->SetNContributors((static_cast<AliAODEvent*>(InputEvent()))->GetPrimaryVertexSPD()->GetNContributors());
317
318  
319          for(int kj=0; kj<(fDielectron->GetTrackArray(0))->GetEntries(); kj++){
320          Int_t posit = nanoEv->AddTrack((AliAODTrack*)fDielectron->GetTrackArray(0)->At(kj));
321          Int_t posVtx = nanoEv->AddVertex(((AliAODTrack*)fDielectron->GetTrackArray(0)->At(kj))->GetProdVertex()); 
322          nanoEv->GetVertex(posVtx)->ResetBit(kIsReferenced);  
323          nanoEv->GetVertex(posVtx)->SetUniqueID(0);
324          nanoEv->GetVertex(posVtx)->RemoveDaughters();
325          nanoEv->GetTrack(posit)->ResetBit(kIsReferenced);  
326          nanoEv->GetTrack(posit)->SetUniqueID(0); 
327          // calo cluster
328          Int_t caloIndex = ((AliAODTrack*)fDielectron->GetTrackArray(0)->At(kj))->GetEMCALcluster();
329          if(caloIndex > 0 && (static_cast<AliAODEvent*>(InputEvent()))->GetCaloCluster(caloIndex)){
330          Int_t posCaloCls = nanoEv->AddCaloCluster(static_cast<AliAODEvent*>(InputEvent())->GetCaloCluster(caloIndex));
331          nanoEv->GetTrack(posit)->SetEMCALcluster(posCaloCls);
332          AliAODCaloCluster *clCls = nanoEv->GetCaloCluster(posCaloCls);
333          for(int u=0; u<clCls->GetNTracksMatched(); u++) clCls->RemoveTrackMatched(clCls->GetTrackMatched(u)); 
334          nanoEv->GetCaloCluster(posCaloCls)->AddTrackMatched((AliAODTrack*)nanoEv->GetTrack(posit)); 
335          }
336          // set references for vtx
337          nanoEv->GetTrack(posit)->SetProdVertex(nanoEv->GetVertex(posVtx));
338          }
339
340          for(int kj=0; kj<(fDielectron->GetTrackArray(1))->GetEntries(); kj++){
341          Int_t negat = nanoEv->AddTrack((AliAODTrack*)fDielectron->GetTrackArray(1)->At(kj));
342          Int_t negVtx = nanoEv->AddVertex(((AliAODTrack*)fDielectron->GetTrackArray(1)->At(kj))->GetProdVertex());    
343          nanoEv->GetVertex(negVtx)->ResetBit(kIsReferenced);  
344          nanoEv->GetVertex(negVtx)->SetUniqueID(0);
345          nanoEv->GetVertex(negVtx)->RemoveDaughters();
346          nanoEv->GetTrack(negat)->ResetBit(kIsReferenced);  
347          nanoEv->GetTrack(negat)->SetUniqueID(0);
348          // calo cluster
349          Int_t caloIndex = ((AliAODTrack*)fDielectron->GetTrackArray(1)->At(kj))->GetEMCALcluster(); 
350          if(caloIndex > 0 && (static_cast<AliAODEvent*>(InputEvent()))->GetCaloCluster(caloIndex)){
351          Int_t negCaloCls = nanoEv->AddCaloCluster(static_cast<AliAODEvent*>(InputEvent())->GetCaloCluster(caloIndex));
352          nanoEv->GetTrack(negat)->SetEMCALcluster(negCaloCls);
353          AliAODCaloCluster *clCls = nanoEv->GetCaloCluster(negCaloCls); 
354          for(int u=0; u<clCls->GetNTracksMatched(); u++) clCls->RemoveTrackMatched(clCls->GetTrackMatched(u));
355          nanoEv->GetCaloCluster(negCaloCls)->AddTrackMatched((AliAODTrack*)nanoEv->GetTrack(negat));
356          }
357          nanoEv->GetTrack(negat)->SetProdVertex(nanoEv->GetVertex(negVtx)); 
358          }  
359         delete tmp; delete tmpSpd; 
360         nanoEv->GetTracks()->Expand(nTracks); 
361         nanoEv->GetVertices()->Expand(nTracks+2); 
362         nanoEv->GetCaloClusters()->Expand(nanoEv->GetNumberOfCaloClusters()); 
363    
364     }
365     if(isAOD) t->Fill();
366   }
367   
368   PostData(1, const_cast<THashList*>(fDielectron->GetHistogramList()));
369   PostData(2,fEventStat);
370   return;
371 }
372