]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliAnalysisTaskSE.cxx
coverity
[u/mrichter/AliRoot.git] / ANALYSIS / AliAnalysisTaskSE.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 /* $Id$ */
17  
18 #include <TROOT.h>
19 #include <TSystem.h>
20 #include <TInterpreter.h>
21 #include <TChain.h>
22 #include <TFile.h>
23 #include <TList.h>
24
25 #include "AliAnalysisTaskSE.h"
26 #include "AliAnalysisManager.h"
27 #include "AliAnalysisCuts.h"
28 #include "AliAnalysisDataSlot.h"
29 #include "AliAnalysisDataContainer.h"
30
31 #include "AliESDEvent.h"
32 #include "AliESDfriend.h"
33 #include "AliESD.h"
34 #include "AliAODEvent.h"
35 #include "AliAODHeader.h"
36 #include "AliAODTracklets.h"
37 #include "AliAODCaloCells.h"
38 #include "AliAODCaloTrigger.h"
39 #include "AliAODMCParticle.h"
40 #include "AliVEvent.h"
41 #include "AliAODHandler.h"
42 #include "AliAODInputHandler.h"
43 #include "AliMCEventHandler.h"
44 #include "AliInputEventHandler.h"
45 #include "AliMultiInputEventHandler.h"
46 #include "AliESDInputHandler.h"
47 #include "AliMCEvent.h"
48 #include "AliStack.h"
49 #include "AliLog.h"
50 #include "AliAODDimuon.h"
51
52
53 ClassImp(AliAnalysisTaskSE)
54
55 ////////////////////////////////////////////////////////////////////////
56 AliAODHeader*    AliAnalysisTaskSE::fgAODHeader         = NULL;
57 TClonesArray*    AliAnalysisTaskSE::fgAODTracks         = NULL;
58 TClonesArray*    AliAnalysisTaskSE::fgAODVertices       = NULL;
59 TClonesArray*    AliAnalysisTaskSE::fgAODV0s            = NULL;
60 TClonesArray*    AliAnalysisTaskSE::fgAODPMDClusters    = NULL;
61 TClonesArray*    AliAnalysisTaskSE::fgAODJets           = NULL;
62 TClonesArray*    AliAnalysisTaskSE::fgAODFMDClusters    = NULL;
63 TClonesArray*    AliAnalysisTaskSE::fgAODCaloClusters   = NULL;
64 TClonesArray*    AliAnalysisTaskSE::fgAODMCParticles    = NULL;
65 AliAODTracklets* AliAnalysisTaskSE::fgAODTracklets      = NULL;
66 AliAODCaloCells* AliAnalysisTaskSE::fgAODEmcalCells     = NULL;
67 AliAODCaloCells* AliAnalysisTaskSE::fgAODPhosCells      = NULL;
68 AliAODCaloTrigger* AliAnalysisTaskSE::fgAODEMCALTrigger = NULL;
69 AliAODCaloTrigger* AliAnalysisTaskSE::fgAODPHOSTrigger  = NULL;
70 TClonesArray*    AliAnalysisTaskSE::fgAODDimuons        = NULL;
71
72 AliAnalysisTaskSE::AliAnalysisTaskSE():
73     AliAnalysisTask(),
74     fDebug(0),
75     fEntry(0),
76     fInputEvent(0x0),
77     fESDfriend(0x0),
78     fInputHandler(0x0),
79     fOutputAOD(0x0),
80     fMCEvent(0x0),
81     fTreeA(0x0),
82     fCurrentRunNumber(-1),
83     fHistosQA(0x0),
84     fOfflineTriggerMask(0),
85     fMultiInputHandler(0),
86     fMCEventHandler(0)
87 {
88   // Default constructor
89 }
90
91 AliAnalysisTaskSE::AliAnalysisTaskSE(const char* name):
92     AliAnalysisTask(name, "AnalysisTaskSE"),
93     fDebug(0),
94     fEntry(0),
95     fInputEvent(0x0),
96     fESDfriend(0x0),
97     fInputHandler(0x0),
98     fOutputAOD(0x0),
99     fMCEvent(0x0),
100     fTreeA(0x0),
101     fCurrentRunNumber(-1),
102     fHistosQA(0x0),
103     fOfflineTriggerMask(0),
104     fMultiInputHandler(0),
105     fMCEventHandler(0)
106 {
107   // Default constructor
108     DefineInput (0, TChain::Class());
109     DefineOutput(0,  TTree::Class());
110 }
111
112 AliAnalysisTaskSE::AliAnalysisTaskSE(const AliAnalysisTaskSE& obj):
113     AliAnalysisTask(obj),
114     fDebug(0),
115     fEntry(0),
116     fInputEvent(0x0),
117     fESDfriend(0x0),
118     fInputHandler(0x0),
119     fOutputAOD(0x0),
120     fMCEvent(0x0),
121     fTreeA(0x0),
122     fCurrentRunNumber(-1),
123     fHistosQA(0x0),
124     fOfflineTriggerMask(0),
125     fMultiInputHandler(obj.fMultiInputHandler),
126     fMCEventHandler(obj.fMCEventHandler)
127 {
128 // Copy constructor
129     fDebug            = obj.fDebug;
130     fEntry            = obj.fEntry;
131     fInputEvent       = obj.fInputEvent;
132     fESDfriend        = obj.fESDfriend;
133     fInputHandler     = obj.fInputHandler;
134     fOutputAOD        = obj.fOutputAOD;
135     fMCEvent          = obj.fMCEvent;
136     fTreeA            = obj.fTreeA;    
137     fCurrentRunNumber = obj.fCurrentRunNumber;
138     fHistosQA         = obj.fHistosQA;
139
140 }
141
142
143 AliAnalysisTaskSE& AliAnalysisTaskSE::operator=(const AliAnalysisTaskSE& other)
144 {
145 // Assignment
146     AliAnalysisTask::operator=(other);
147     fDebug            = other.fDebug;
148     fEntry            = other.fEntry;
149     fInputEvent       = other.fInputEvent;
150     fESDfriend        = other.fESDfriend;
151     fInputHandler     = other.fInputHandler;
152     fOutputAOD        = other.fOutputAOD;
153     fMCEvent          = other.fMCEvent;
154     fTreeA            = other.fTreeA;    
155     fCurrentRunNumber = other.fCurrentRunNumber;
156     fHistosQA         = other.fHistosQA;
157     fOfflineTriggerMask = other.fOfflineTriggerMask;
158     fMultiInputHandler  = other.fMultiInputHandler;
159     fMCEventHandler     = other.fMCEventHandler;
160     return *this;
161 }
162
163 //______________________________________________________________________________
164 void AliAnalysisTaskSE::ConnectInputData(Option_t* /*option*/)
165 {
166 // Connect the input data
167     if (fDebug > 1) printf("AnalysisTaskSE::ConnectInputData() \n");
168
169    // Connect input handlers (multi input handler is handled)
170     ConnectMultiHandler();
171     
172     if (fInputHandler) {
173         if ((fInputHandler->GetTree())->GetBranch("ESDfriend."))
174             fESDfriend = ((AliESDInputHandler*)fInputHandler)->GetESDfriend();
175
176         fInputEvent = fInputHandler->GetEvent();
177     } else if( fMCEvent ) {
178          AliWarning("No Input Event Handler connected, only MC Truth Event Handler") ; 
179     } else {
180          AliError("No Input Event Handler connected") ; 
181          return ; 
182     }
183     // Disconnect multi handler
184     DisconnectMultiHandler();
185 }
186
187 void AliAnalysisTaskSE::CreateOutputObjects()
188 {
189 // Create the output container
190 //
191 //  Default AOD
192     if (fDebug > 1) printf("AnalysisTaskSE::CreateOutPutData() \n");
193
194     AliAODHandler* handler = dynamic_cast<AliAODHandler*> 
195          ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
196     
197     Bool_t merging = kFALSE;
198     AliAODInputHandler* aodIH = static_cast<AliAODInputHandler*>((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
199     if (aodIH) {
200         if (aodIH->GetMergeEvents()) merging = kTRUE;
201     }
202
203     // Check if AOD replication has been required
204     if (handler) {
205         fOutputAOD   = handler->GetAOD();
206         fTreeA = handler->GetTree();
207         if (fOutputAOD && !(handler->IsStandard())) {
208             if ((handler->NeedsHeaderReplication() || merging) && !(fgAODHeader)) 
209                 {
210                  if (fDebug > 1) AliInfo("Replicating header");
211                  fgAODHeader = new AliAODHeader;
212                  handler->AddBranch("AliAODHeader", &fgAODHeader);
213                 }
214             if ((handler->NeedsTracksBranchReplication() || merging) && !(fgAODTracks))      
215             {   
216                 if (fDebug > 1) AliInfo("Replicating track branch\n");
217                 fgAODTracks = new TClonesArray("AliAODTrack",500);
218                 fgAODTracks->SetName("tracks");
219                 handler->AddBranch("TClonesArray", &fgAODTracks);
220             }    
221             if ((handler->NeedsVerticesBranchReplication() || merging) && !(fgAODVertices))
222             {
223                 if (fDebug > 1) AliInfo("Replicating vertices branch\n");
224                 fgAODVertices = new TClonesArray("AliAODVertex",500);
225                 fgAODVertices->SetName("vertices");
226                 handler->AddBranch("TClonesArray", &fgAODVertices);
227             }   
228             if ((handler->NeedsV0sBranchReplication()) && !(fgAODV0s))    
229             {   
230                 if (fDebug > 1) AliInfo("Replicating V0s branch\n");
231                 fgAODV0s = new TClonesArray("AliAODv0",500);
232                 fgAODV0s->SetName("v0s");
233                 handler->AddBranch("TClonesArray", &fgAODV0s);
234             }
235             if ((handler->NeedsTrackletsBranchReplication()) && !(fgAODTracklets))        
236             {   
237                 if (fDebug > 1) AliInfo("Replicating Tracklets branch\n");
238                 fgAODTracklets = new AliAODTracklets("tracklets","tracklets");
239                 handler->AddBranch("AliAODTracklets", &fgAODTracklets);
240             }
241             if ((handler->NeedsPMDClustersBranchReplication()) && !(fgAODPMDClusters))    
242             {   
243                 if (fDebug > 1) AliInfo("Replicating PMDClusters branch\n");
244                 fgAODPMDClusters = new TClonesArray("AliAODPmdCluster",500);
245                 fgAODPMDClusters->SetName("pmdClusters");
246                 handler->AddBranch("TClonesArray", &fgAODPMDClusters);
247             }
248             if ((handler->NeedsJetsBranchReplication() || merging) && !(fgAODJets))       
249             {   
250                 if (fDebug > 1) AliInfo("Replicating Jets branch\n");
251                 fgAODJets = new TClonesArray("AliAODJet",500);
252                 fgAODJets->SetName("jets");
253                 handler->AddBranch("TClonesArray", &fgAODJets);
254             }
255             if ((handler->NeedsFMDClustersBranchReplication()) && !(fgAODFMDClusters))    
256             {   
257                 AliInfo("Replicating FMDClusters branch\n");
258                 fgAODFMDClusters = new TClonesArray("AliAODFmdCluster",500);
259                 fgAODFMDClusters->SetName("fmdClusters");
260                 handler->AddBranch("TClonesArray", &fgAODFMDClusters);
261             }
262             if ((handler->NeedsCaloClustersBranchReplication() || merging) && !(fgAODCaloClusters))       
263             {   
264                 if (fDebug > 1) AliInfo("Replicating CaloClusters branch\n");
265                 fgAODCaloClusters = new TClonesArray("AliAODCaloCluster",500);
266                 fgAODCaloClusters->SetName("caloClusters");
267                 handler->AddBranch("TClonesArray", &fgAODCaloClusters);
268
269                 fgAODEmcalCells = new AliAODCaloCells("emcalCells","emcalCells",AliVCaloCells::kEMCALCell);
270                 handler->AddBranch("AliAODCaloCells", &fgAODEmcalCells);
271
272                 fgAODPhosCells = new AliAODCaloCells("phosCells","phosCells",AliVCaloCells::kPHOSCell);
273                 handler->AddBranch("AliAODCaloCells", &fgAODPhosCells);
274                 }
275             if ((handler->NeedsCaloTriggerBranchReplication() || merging) && !(fgAODEMCALTrigger))        
276             {   
277                 if (fDebug > 1) AliInfo("Replicating EMCAL Calo Trigger branches\n");
278                 fgAODEMCALTrigger = new AliAODCaloTrigger("emcalTrigger","emcalTrigger");
279                 handler->AddBranch("AliAODCaloTrigger", &fgAODEMCALTrigger);
280                 }
281                 if ((handler->NeedsCaloTriggerBranchReplication() || merging) && !(fgAODPHOSTrigger))     
282                 {   
283                 if (fDebug > 1) AliInfo("Replicating PHOS Calo Trigger branches\n");
284                 fgAODPHOSTrigger = new AliAODCaloTrigger("phosTrigger","phosTrigger");
285                 handler->AddBranch("AliAODCaloTrigger", &fgAODPHOSTrigger);
286             }
287             if ((handler->NeedsMCParticlesBranchReplication() || merging) && !(fgAODMCParticles))         
288             {   
289                 if (fDebug > 1) AliInfo("Replicating MCParticles branch\n");
290                 fgAODMCParticles = new TClonesArray("AliAODMCParticle",500);
291                 fgAODMCParticles->SetName("mcparticles");
292                 handler->AddBranch("TClonesArray", &fgAODMCParticles);
293             }
294             if ((handler->NeedsDimuonsBranchReplication() || merging) && !(fgAODDimuons))      
295             {   
296                 if (fDebug > 1) AliInfo("Replicating dimuon branch\n");
297                 fgAODDimuons = new TClonesArray("AliAODDimuon",0);
298                 fgAODDimuons->SetName("dimuons");
299                 handler->AddBranch("TClonesArray", &fgAODDimuons);
300             }    
301
302             // cache the pointerd in the AODEvent
303             fOutputAOD->GetStdContent();
304         }
305     }
306     ConnectMultiHandler();
307     UserCreateOutputObjects();
308     DisconnectMultiHandler();
309 }
310
311 void AliAnalysisTaskSE::Exec(Option_t* option)
312 {
313 //
314 // Exec analysis of one event
315
316     ConnectMultiHandler();
317
318     if ( fDebug >= 10)
319       printf("Task is active %5d\n", IsActive());
320     
321     if (fDebug > 1) AliInfo("AliAnalysisTaskSE::Exec() \n");
322 //
323     AliAODHandler* handler = dynamic_cast<AliAODHandler*> 
324         ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
325
326     AliAODInputHandler* aodH = dynamic_cast<AliAODInputHandler*>(fInputHandler);
327 //
328 // Was event selected ? If no event selection mechanism, the event SHOULD be selected (AG)
329     UInt_t isSelected = AliVEvent::kAny;
330     if( fInputHandler && (fInputHandler->GetEventSelection() || aodH)) {
331       // Get the actual offline trigger mask for the event and AND it with the
332       // requested mask. If no mask requested select by default the event.
333       if (fOfflineTriggerMask)
334          isSelected = fOfflineTriggerMask & fInputHandler->IsEventSelected();
335     }
336 //  Functionality below moved in the filter tasks (AG)
337 //    if (handler) handler->SetFillAOD(isSelected);
338
339     if( fInputHandler ) {
340         fEntry = fInputHandler->GetReadEntry();
341         fESDfriend = ((AliESDInputHandler*)fInputHandler)->GetESDfriend();
342     }
343     
344
345 // Notify the change of run number
346     if (InputEvent() && (InputEvent()->GetRunNumber() != fCurrentRunNumber)) {
347         fCurrentRunNumber = InputEvent()->GetRunNumber();
348         NotifyRun();
349     }    
350            
351     else if( fMCEvent )
352        fEntry = fMCEvent->Header()->GetEvent(); 
353     if ( !((Entry()-1)%100) && fDebug > 0) 
354          AliInfo(Form("%s ----> Processing event # %lld", CurrentFileName(), Entry()));
355
356     
357     
358
359     if (handler && aodH) {
360         fMCEvent = aodH->MCEvent();
361         Bool_t merging = aodH->GetMergeEvents();
362       
363   // Do not analyze merged events if last embedded file has less events than normal event, 
364   // skip analysis after last embeded event 
365   if(merging){
366     if(aodH->GetReadEntry() + aodH->GetMergeOffset() >= aodH->GetTreeToMerge()->GetEntriesFast()){
367       //printf("Skip Entry %lld, Offset %d, Tree Entries %d\n",aodH->GetReadEntry(),aodH->GetMergeOffset(), aodH->GetTreeToMerge()->GetEntries());
368           
369       // Do I need to add the lines before the return?
370       // Added protection in case the derived task is not an AOD producer.
371       AliAnalysisDataSlot *out0 = GetOutputSlot(0);
372       if (out0 && out0->IsConnected()) PostData(0, fTreeA);    
373           
374       DisconnectMultiHandler();
375           
376       return;
377     }
378     //else   printf("MERGE Entry %lld, Offset %d, Tree Entries %d\n",aodH->GetReadEntry(),aodH->GetMergeOffset(), aodH->GetTreeToMerge()->GetEntries());
379   }
380       
381         AliAODEvent* aod = dynamic_cast<AliAODEvent*>(InputEvent());
382
383         if (aod && !(handler->IsStandard()) && !(handler->AODIsReplicated())) {
384             if ((handler->NeedsHeaderReplication() || merging) && (fgAODHeader))
385             {
386               // copy the contents by assigment
387               *fgAODHeader =  *(aod->GetHeader());
388             }
389             if ((handler->NeedsTracksBranchReplication() || (merging &&  aodH->GetMergeTracks())) && (fgAODTracks))
390             {
391                 TClonesArray* tracks = aod->GetTracks();
392                 new (fgAODTracks) TClonesArray(*tracks);
393             }
394             if ((handler->NeedsVerticesBranchReplication() || merging) && (fgAODVertices))
395             {
396                 TClonesArray* vertices = aod->GetVertices();
397                 new (fgAODVertices) TClonesArray(*vertices);
398             }
399             if ((handler->NeedsV0sBranchReplication()) && (fgAODV0s))
400             {
401                 TClonesArray* v0s = aod->GetV0s();
402                 new (fgAODV0s) TClonesArray(*v0s);
403             }
404             if ((handler->NeedsTrackletsBranchReplication()) && (fgAODTracklets))
405             {
406               *fgAODTracklets = *aod->GetTracklets();
407             }
408             if ((handler->NeedsPMDClustersBranchReplication()) && (fgAODPMDClusters))
409             {
410                 TClonesArray* pmdClusters = aod->GetPmdClusters();
411                 new (fgAODPMDClusters) TClonesArray(*pmdClusters);
412             }
413             if ((handler->NeedsJetsBranchReplication() || (merging &&aodH->GetMergeTracks())) && (fgAODJets))
414             {
415                 TClonesArray* jets = aod->GetJets();
416                 new (fgAODJets) TClonesArray(*jets);
417             }
418             if ((handler->NeedsFMDClustersBranchReplication()) && (fgAODFMDClusters))
419             {
420                 TClonesArray* fmdClusters = aod->GetFmdClusters();
421                 new (fgAODFMDClusters) TClonesArray(*fmdClusters);
422             }
423             if ((handler->NeedsCaloClustersBranchReplication() || 
424                  (merging && (aodH->GetMergeEMCALClusters() || aodH->GetMergePHOSClusters()))) 
425                 && (fgAODCaloClusters))
426             {
427                 TClonesArray* caloClusters = aod->GetCaloClusters();
428                 new (fgAODCaloClusters) TClonesArray(*caloClusters);
429             }
430
431             if ((handler->NeedsMCParticlesBranchReplication() || merging) && (fgAODMCParticles))
432             {
433                 TClonesArray* mcParticles = (TClonesArray*) (aod->FindListObject("mcparticles"));
434                 if( mcParticles )
435                   new (fgAODMCParticles) TClonesArray(*mcParticles);
436             }
437             
438             if ((handler->NeedsDimuonsBranchReplication() || (merging && aodH->GetMergeTracks())) && (fgAODDimuons))
439             {
440                 fgAODDimuons->Clear();
441                 TClonesArray& dimuons = *fgAODDimuons;
442                 TClonesArray& tracksnew = *fgAODTracks;
443                 
444                 Int_t nMuonTrack[100]; 
445                 for(Int_t imuon = 0; imuon < 100; imuon++) nMuonTrack[imuon] = 0;
446                 Int_t nMuons=0;
447                 for(Int_t ii=0; ii < fgAODTracks->GetEntries(); ii++){
448                     AliAODTrack *track = (AliAODTrack*) fgAODTracks->At(ii);
449                     if(track->IsMuonTrack()) {
450                         nMuonTrack[nMuons]= ii;
451                         nMuons++;
452                     }  
453                 }
454                 Int_t jDimuons=0;
455                 if(nMuons >= 2){
456                     for(Int_t i = 0; i < nMuons; i++){
457                         Int_t index0 = nMuonTrack[i];
458                         for(Int_t j = i+1; j < nMuons; j++){
459                             Int_t index1 = nMuonTrack[j];
460                             tracksnew.At(index0)->ResetBit(kIsReferenced);
461                             tracksnew.At(index0)->SetUniqueID(0); 
462                             tracksnew.At(index1)->ResetBit(kIsReferenced);
463                             tracksnew.At(index1)->SetUniqueID(0);
464                             new(dimuons[jDimuons++]) AliAODDimuon(tracksnew.At(index0),tracksnew.At(index1));
465                         }
466                     }    
467                 }
468             }
469             // Additional merging if needed
470             if (merging) {
471               Int_t nc;
472
473                 // mcParticles
474                 TClonesArray* mcparticles = (TClonesArray*) ((aodH->GetEventToMerge())->FindListObject("mcparticles"));
475                 if( mcparticles ){
476                   Int_t npart = mcparticles->GetEntries();
477                   nc = fgAODMCParticles->GetEntries();
478                   
479                   for (Int_t i = 0; i < npart; i++) {
480                     AliAODMCParticle* particle = (AliAODMCParticle*) mcparticles->At(i);
481                     new((*fgAODMCParticles)[nc++]) AliAODMCParticle(*particle);
482                   }
483                 }
484
485                 // tracks
486                 TClonesArray* tracks = aodH->GetEventToMerge()->GetTracks();
487                 if(tracks && aodH->GetMergeTracks()){
488                   Int_t ntr = tracks->GetEntries();
489                   nc  = fgAODTracks->GetEntries();      
490                   for (Int_t i = 0; i < ntr; i++) {
491                     AliAODTrack*    track = (AliAODTrack*) tracks->At(i);
492                     AliAODTrack* newtrack = new((*fgAODTracks)[nc++]) AliAODTrack(*track);
493                     newtrack->SetLabel(newtrack->GetLabel() + fgAODMCParticles->GetEntries());
494                   }
495                   
496                   for (Int_t i = 0; i < nc; i++) 
497                     {
498                       AliAODTrack* track = (AliAODTrack*) fgAODTracks->At(i);
499                       track->ResetBit(kIsReferenced);
500                       track->SetUniqueID(0);
501                     }
502                 }
503                 
504                 // clusters
505                 TClonesArray* clusters = aodH->GetEventToMerge()->GetCaloClusters();
506                 if( clusters  && (aodH->GetMergeEMCALClusters() || aodH->GetMergePHOSClusters())) {
507                   Int_t ncl  = clusters->GetEntries();
508                   nc         =  fgAODCaloClusters->GetEntries();
509                   for (Int_t i = 0; i < ncl; i++) {
510                     AliAODCaloCluster*    cluster = (AliAODCaloCluster*) clusters->At(i);
511                     if(cluster->IsEMCAL() && !aodH->GetMergeEMCALClusters() ) continue;
512                     if(cluster->IsPHOS()  && !aodH->GetMergePHOSClusters()  ) continue;   
513                     new((*fgAODCaloClusters)[nc++]) AliAODCaloCluster(*cluster);
514                   }
515                 }
516
517                 // EMCAL cells
518                 //*fgAODEmcalCells =  *(aod->GetEMCALCells()); // This will be valid after 10.Mar.2011.
519                 if(aodH->GetMergeEMCALCells()) {
520                     AliAODCaloCells* copycells = aod->GetEMCALCells();
521                     fgAODEmcalCells->CreateContainer(copycells->GetNumberOfCells());
522                     nc  = copycells->GetNumberOfCells();
523                     while( nc-- ){ fgAODEmcalCells->SetCell(nc,copycells->GetCellNumber(nc),copycells->GetAmplitude(nc)); }
524
525                     AliAODCaloCells* cellsA = aodH->GetEventToMerge()->GetEMCALCells();
526                     if( cellsA ){
527                         Int_t ncells  = cellsA->GetNumberOfCells();
528                         nc = fgAODEmcalCells->GetNumberOfCells();
529                         for (Int_t i  = 0; i < ncells; i++) {
530                             Int_t cn  = cellsA->GetCellNumber(i);
531                             Int_t pos = fgAODEmcalCells->GetCellPosition(cn);
532                             if (pos >= 0) {
533                                 Double_t amp = cellsA->GetAmplitude(i) + fgAODEmcalCells->GetAmplitude(pos);
534                                 fgAODEmcalCells->SetCell(pos, cn, amp);
535                             } else {
536                                 AliAODCaloCells* copycells1 = new AliAODCaloCells(*fgAODEmcalCells);
537                                 fgAODEmcalCells->CreateContainer(nc+1);
538                                 Int_t nn = copycells1->GetNumberOfCells();
539                                 while( nn-- ){ fgAODEmcalCells->SetCell(nn,copycells1->GetCellNumber(nn),copycells1->GetAmplitude(nn)); }
540                                 fgAODEmcalCells->SetCell(nc++,cn,cellsA->GetAmplitude(i));
541                                 delete copycells1;
542                             }
543                         }
544                         fgAODEmcalCells->Sort();
545                     }
546                 } // merge emcal cells
547                 
548                 
549                 // PHOS cells
550                 //*fgAODPhosCells =  *(aod->GetPHOSCells()); // This will be valid after 10.Mar.2011.
551                 if(aodH->GetMergePHOSCells()) {
552                     AliAODCaloCells* copycells = aod->GetPHOSCells();
553                     fgAODPhosCells->CreateContainer(copycells->GetNumberOfCells());
554                     nc  = copycells->GetNumberOfCells();
555                     while( nc-- ){ fgAODPhosCells->SetCell(nc,copycells->GetCellNumber(nc),copycells->GetAmplitude(nc)); }
556                     AliAODCaloCells* cellsP = aodH->GetEventToMerge()->GetPHOSCells();
557                     if( cellsP ){
558                         Int_t ncellsP  = cellsP->GetNumberOfCells();
559                         nc = fgAODPhosCells->GetNumberOfCells();
560                         
561                         for (Int_t i  = 0; i < ncellsP; i++) {
562                             Int_t cn  = cellsP->GetCellNumber(i);
563                             Int_t pos = fgAODPhosCells->GetCellPosition(cn);
564                             if (pos >= 0) {
565                                 Double_t amp = cellsP->GetAmplitude(i) + fgAODPhosCells->GetAmplitude(pos);
566                                 fgAODPhosCells->SetCell(pos, cn, amp);
567                             } else {
568                                 AliAODCaloCells* copycells1 = new AliAODCaloCells(*fgAODPhosCells);
569                                 fgAODPhosCells->CreateContainer(nc+1);
570                                 Int_t nn = copycells1->GetNumberOfCells();
571                                 while( nn-- ){ fgAODPhosCells->SetCell(nn,copycells1->GetCellNumber(nn),copycells1->GetAmplitude(nn)); }
572                                 fgAODPhosCells->SetCell(nc++,cn,cellsP->GetAmplitude(i));
573                                 delete copycells1;
574                             }
575                         }
576                         fgAODPhosCells->Sort();
577                     }
578                 } // Merge PHOS Cells
579                 
580                 if (aodH->GetMergeEMCALTrigger()) 
581                 {
582                         Int_t   EMCALts[48][64], px, py, ts;
583                         Float_t EMCALfo[48][64], am;
584                         
585                         for (Int_t i = 0; i < 48; i++) for (Int_t j = 0; j < 64; j++) 
586                         {
587                                 EMCALts[i][j] = 0;
588                                 EMCALfo[i][j] = 0.;
589                         }
590                         
591                         AliAODCaloTrigger& trg0 = *(aod->GetCaloTrigger("EMCAL"));
592                         
593                         trg0.Reset();
594                         while (trg0.Next())
595                         {
596                                 trg0.GetPosition(px, py);
597                                 
598                                 if (px > -1 && py > -1) 
599                                 {
600                                         trg0.GetL1TimeSum(ts);
601                                         if (ts > -1) EMCALts[px][py] += ts;
602                                         
603                                         trg0.GetAmplitude(am);
604                                         if (am > -1) EMCALfo[px][py] += am;
605                                 }
606                         }
607                         
608                         AliAODCaloTrigger& trg1 = *((aodH->GetEventToMerge())->GetCaloTrigger("EMCAL"));
609                         
610                         trg1.Reset();
611                         while (trg1.Next())
612                         {
613                                 trg1.GetPosition(px, py);
614                                 
615                                 if (px > -1 && py > -1) 
616                                 {
617                                         trg1.GetL1TimeSum(ts);
618                                         if (ts > -1) EMCALts[px][py] += ts;
619                                         
620                                         trg1.GetAmplitude(am);
621                                         if (am > -1) EMCALfo[px][py] += am;
622                                 }
623                         }
624                         
625                         int nEntries = 0;
626                         for (Int_t i = 0; i < 48; i++) 
627                                 for (Int_t j = 0; j < 64; j++) 
628                                         if (EMCALts[i][j] || EMCALfo[i][j]) nEntries++;
629                 
630                         fgAODEMCALTrigger->Allocate(nEntries);
631                         Int_t L0times[10]; for (int i = 0; i < 10; i++) L0times[i] = -1;
632                 
633                         for (Int_t i = 0; i < 48; i++) 
634                                 for (Int_t j = 0; j < 64; j++) 
635                                         if (EMCALts[i][j] || EMCALfo[i][j]) 
636                                                 fgAODEMCALTrigger->Add(i, j, EMCALfo[i][j], -1., L0times, 0, EMCALts[i][j], 0);
637                 }
638                 
639                 if (aodH->GetMergePHOSTrigger()) 
640                 {
641                         // To be implemented by PHOS
642                 }
643             } // merging
644             
645             handler->SetAODIsReplicated();
646         }
647     }
648
649
650 // Call the user analysis    
651     if (!fMCEventHandler) {
652         if (isSelected) 
653             UserExec(option);
654     } else {
655         if (isSelected && (fMCEventHandler->InitOk())) 
656             UserExec(option);
657     }
658     
659 // Added protection in case the derived task is not an AOD producer.
660     AliAnalysisDataSlot *out0 = GetOutputSlot(0);
661     if (out0 && out0->IsConnected()) PostData(0, fTreeA);    
662
663     DisconnectMultiHandler();
664 }
665
666 const char* AliAnalysisTaskSE::CurrentFileName()
667 {
668 // Returns the current file name    
669     if( fInputHandler )
670       return fInputHandler->GetTree()->GetCurrentFile()->GetName();
671     else if( fMCEvent )
672       return ((AliMCEventHandler*) ((AliAnalysisManager::GetAnalysisManager())->GetMCtruthEventHandler()))->TreeK()->GetCurrentFile()->GetName();
673     else return "";
674 }
675
676 void AliAnalysisTaskSE::AddAODBranch(const char* cname, void* addobj, const char *fname)
677 {
678     // Add a new branch to the aod tree
679     AliAODHandler* handler = dynamic_cast<AliAODHandler*> 
680         ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
681     if (handler) {
682         handler->AddBranch(cname, addobj, fname);
683     }
684 }
685
686 Bool_t AliAnalysisTaskSE::IsStandardAOD() const
687 {
688 // Check if the output AOD handler is configured for standard or delta AOD.
689 // Users should first check that AODEvent() returns non-null.
690     AliAODHandler* handler = dynamic_cast<AliAODHandler*> 
691          ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
692     if (!handler) {
693        Error("IsStandardAOD", "No AOD handler. Please use AODEvent() to check this first");
694        return kTRUE;
695     }
696     return handler->IsStandard();   
697 }
698
699 Bool_t AliAnalysisTaskSE::Notify()
700 {
701     return (UserNotify());
702 }
703
704 const AliEventTag *AliAnalysisTaskSE::EventTag() const
705 {
706 // Returns tag for the current event, if any. The return value should always be checked by the user.
707    if (!fInputHandler) {
708       Error("EventTag", "Input handler not yet available. Call this in UserExec");
709       return NULL;
710    }
711    return fInputHandler->GetEventTag();
712 }
713
714 void AliAnalysisTaskSE::LoadBranches() const
715 {
716 // Load all branches declared in fBranchNames data member of the parent class.
717 // Should be called in UserExec.
718   if (!fInputHandler) {
719      Error("LoadBranches", "Input handler not available yet. Call this in UserExec");
720      return;
721   }
722   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
723   if (mgr->GetAutoBranchLoading()) return;
724   TString taskbranches;
725   GetBranches(fInputHandler->GetDataType(), taskbranches);
726   if (taskbranches.IsNull()) return;
727   TObjArray *arr = taskbranches.Tokenize(",");
728   TIter next(arr);
729   TObject *obj;
730   while ((obj=next())) mgr->LoadBranch(obj->GetName());
731 }
732
733
734 //_________________________________________________________________________________________________
735 void AliAnalysisTaskSE::ConnectMultiHandler()
736 {
737    //
738    // Connect MultiHandler
739    //
740    fInputHandler = (AliInputEventHandler *)((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
741    fMultiInputHandler = dynamic_cast<AliMultiInputEventHandler *>(fInputHandler);
742    if (fMultiInputHandler) {
743       fInputHandler = dynamic_cast<AliInputEventHandler *>(fMultiInputHandler->GetFirstInputEventHandler());
744       fMCEventHandler = dynamic_cast<AliMCEventHandler *>(fMultiInputHandler->GetFirstMCEventHandler());
745    } else { 
746       fMCEventHandler = dynamic_cast<AliMCEventHandler *>((AliAnalysisManager::GetAnalysisManager())->GetMCtruthEventHandler());
747    }
748    if (fMCEventHandler) fMCEvent = fMCEventHandler->MCEvent();
749 }
750
751 //_________________________________________________________________________________________________
752 void AliAnalysisTaskSE::DisconnectMultiHandler()
753 {
754    //
755    // Disconnect MultiHandler
756    //
757    if (fMultiInputHandler) fInputHandler = fMultiInputHandler;
758 }