]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/PartCorrBase/AliAnalysisTaskParticleCorrelation.cxx
Problem with overwriting memory corrected.
[u/mrichter/AliRoot.git] / PWG4 / PartCorrBase / AliAnalysisTaskParticleCorrelation.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 /* $Id: $ */
16
17 //_________________________________________________________________________
18 // Analysis task that executes the analysis classes
19 // that depend on the PartCorr frame, frame for Particle identification and correlations.
20 // Specially designed for calorimeters but also can be used for charged tracks
21 // Input of this task is a configuration file that contains all the settings of the analyis
22 //
23 // -- Author: Gustavo Conesa (INFN-LNF)
24
25 #include <cstdlib>
26
27 // --- Root ---
28 #include <TROOT.h>
29 #include <TInterpreter.h>
30 //#include <Riostream.h>
31
32 // --- Analysis ---
33 #include "AliAnalysisTaskParticleCorrelation.h"
34 #include "AliAnaPartCorrMaker.h"
35 #include "AliCaloTrackReader.h"
36 #include "AliPDG.h"
37 #include "AliAnalysisManager.h"
38
39 ClassImp(AliAnalysisTaskParticleCorrelation)
40
41 ////////////////////////////////////////////////////////////////////////
42 AliAnalysisTaskParticleCorrelation::AliAnalysisTaskParticleCorrelation():
43   AliAnalysisTaskSE(),
44   fAna(0x0),
45   fOutputContainer(0x0),
46   fConfigName("")
47 {
48   // Default constructor
49 }
50
51 //_____________________________________________________
52 AliAnalysisTaskParticleCorrelation::AliAnalysisTaskParticleCorrelation(const char* name):
53   AliAnalysisTaskSE(name),
54   fAna(0x0),
55   fOutputContainer(0x0),
56   fConfigName("")
57 {
58   // Default constructor
59   
60   DefineOutput(1, TList::Class());
61
62 }
63
64 //_____________________________________________________
65 AliAnalysisTaskParticleCorrelation::~AliAnalysisTaskParticleCorrelation() 
66 {
67   // Remove all pointers
68   
69   if(fOutputContainer){
70     fOutputContainer->Clear() ; 
71     delete fOutputContainer ;
72   }
73
74 }
75
76 //_____________________________________________________
77 void AliAnalysisTaskParticleCorrelation::UserCreateOutputObjects()
78 {
79   // Create the output container
80   if (DebugLevel() > 1) printf("AliAnalysisTaskParticleCorrelation::UserCreateOutputObjects() - Begin\n");
81   
82   //Get list of aod arrays, add each aod array to analysis frame 
83   TClonesArray * array = 0;
84   TList * list = fAna->GetAODBranchList();
85   TString deltaAODName = (fAna->GetReader())->GetDeltaAODFileName();
86   for(Int_t iaod = 0; iaod < list->GetEntries(); iaod++){
87     array = (TClonesArray*) list->At(iaod);
88         if(deltaAODName!="") AddAODBranch("TClonesArray", &array, deltaAODName);//Put it in DeltaAOD file
89         else AddAODBranch("TClonesArray", &array);//Put it in standard AOD file
90   } 
91   
92   //Histograms container
93   OpenFile(1);
94   fOutputContainer = fAna->GetOutputContainer();
95   if (DebugLevel() > 1) printf("AliAnalysisTaskParticleCorrelation::UserCreateOutputObjects() - End\n");
96 }
97
98 //_____________________________________________________
99 void AliAnalysisTaskParticleCorrelation::Init()
100 {
101   // Initialization
102  
103   if (DebugLevel() > 1) printf("AliAnalysisTaskParticleCorrelation::Init() - Begin\n");
104   
105   // Call configuration file if specified
106   
107   if (fConfigName.Length()) {
108     printf("AliAnalysisTaskParticleCorrelation::Init() - ### Configuration file is %s.C ###\n", fConfigName.Data());
109         gROOT->LoadMacro(fConfigName+".C");
110         fAna = (AliAnaPartCorrMaker*) gInterpreter->ProcessLine("ConfigAnalysis()");
111   }
112   
113   if(!fAna) {
114         printf("AliAnalysisTaskParticleCorrelation::Init() - Analysis maker pointer not initialized, no analysis specified, STOP !\n");
115         abort();
116   }
117   
118   // Add different generator particles to PDG Data Base 
119   // to avoid problems when reading MC generator particles
120   AliPDG::AddParticlesToPdgDataBase();
121
122   // Initialise analysis
123   fAna->Init();
124
125   //Delta AOD
126   if((fAna->GetReader())->GetDeltaAODFileName()!="")
127           AliAnalysisManager::GetAnalysisManager()->RegisterExtraFile((fAna->GetReader())->GetDeltaAODFileName());
128
129   if (DebugLevel() > 1) printf("AliAnalysisTaskParticleCorrelation::Init() - End\n");
130   
131 }
132
133
134 //_____________________________________________________
135 void AliAnalysisTaskParticleCorrelation::UserExec(Option_t */*option*/)
136 {
137   // Execute analysis for current event
138   //
139   if (DebugLevel() > 1) printf("AliAnalysisTaskParticleCorrelation::UserExec() - Begin\n");
140
141    //Get the type of data, check if type is correct
142   Int_t  datatype = fAna->GetReader()->GetDataType();
143   if(datatype != AliCaloTrackReader::kESD && datatype != AliCaloTrackReader::kAOD &&
144      datatype != AliCaloTrackReader::kMC){
145     printf("AliAnalysisTaskParticleCorrelation::UserExec() - Wrong type of data\n");
146     return ;
147   }
148   
149   fAna->GetReader()->SetInputOutputMCEvent(InputEvent(), AODEvent(), MCEvent());
150
151   //Process event
152   fAna->ProcessEvent((Int_t) Entry(), CurrentFileName());
153   //printf("AliAnalysisTaskParticleCorrelation::Current Event %d; Current File Name : %s\n",(Int_t) Entry(), CurrentFileName());
154   if (DebugLevel() > 1) printf("AliAnalysisTaskParticleCorrelation::UserExec() - End\n");
155   
156   PostData(1, fOutputContainer);
157   
158 }
159
160 //_____________________________________________________
161 void AliAnalysisTaskParticleCorrelation::Terminate(Option_t */*option*/)
162 {
163   // Terminate analysis
164   // Do some plots
165
166   // Get merged histograms from the output container
167   TList *outputList = (TList*)GetOutputData(1);
168   // Propagate histagrams to maker
169   fAna->Terminate(outputList);
170         
171 }
172