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