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