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