]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/FLOW/AliFlowTasks/AliAnalysisTaskCumulants.cxx
Bug fixes (Sadhana)
[u/mrichter/AliRoot.git] / PWG2 / FLOW / AliFlowTasks / AliAnalysisTaskCumulants.cxx
CommitLineData
924fafb7 1/*************************************************************************
2* Copyright(c) 1998-2008, 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
2188af53 16/**************************************
17 * analysis task for cumulant method *
18 * *
19 * authors: Naomi van der Kolk *
20 * (kolk@nikhef.nl) *
21 * Raimond Snellings *
22 * (snelling@nikhef.nl) *
23 * Ante Bilandzic *
24 * (anteb@nikhef.nl) *
25 * ***********************************/
26
2ed70edf 27class TFile;
28class TList;
29class AliAnalysisTaskSE;
30
2188af53 31#include "Riostream.h"
7183fe85 32#include "AliFlowEventSimple.h"
924fafb7 33#include "AliAnalysisTaskCumulants.h"
924fafb7 34#include "AliFlowAnalysisWithCumulants.h"
35
924fafb7 36ClassImp(AliAnalysisTaskCumulants)
37
2188af53 38//================================================================================================================
39
7183fe85 40AliAnalysisTaskCumulants::AliAnalysisTaskCumulants(const char *name, Bool_t useWeights):
2ed70edf 41AliAnalysisTaskSE(name),
42fEvent(NULL),
43fGFC(NULL),
44fListHistos(NULL),
45fUseWeights(useWeights),
46fUsePhiWeights(kFALSE),
47fUsePtWeights(kFALSE),
48fUseEtaWeights(kFALSE),
49fListWeights(NULL)
924fafb7 50{
2ed70edf 51 // Constructor
52 cout<<"AliAnalysisTaskCumulants::AliAnalysisTaskCumulants(const char *name, Bool_t useWeights)"<<endl;
fe488c8a 53
54 // Define input and output slots here
2ed70edf 55 // Input slot #0 works with an AliFlowEventSimple
7183fe85 56 DefineInput(0, AliFlowEventSimple::Class());
e5e75b58 57
2ed70edf 58 // Input slot #1 is needed for the weights input files
59 if(useWeights)
e5e75b58 60 {
61 DefineInput(1, TList::Class());
62 }
2ed70edf 63 // Output slot #0 is reserved
64 // Output slot #1 writes into a TList container
65 DefineOutput(1, TList::Class());
924fafb7 66}
2188af53 67
7183fe85 68AliAnalysisTaskCumulants::AliAnalysisTaskCumulants():
2ed70edf 69AliAnalysisTaskSE(),
70fEvent(NULL),
71fGFC(NULL),
72fListHistos(NULL),
73fUseWeights(kFALSE),
74fUsePhiWeights(kFALSE),
75fUsePtWeights(kFALSE),
76fUseEtaWeights(kFALSE),
77fListWeights(NULL)
aaebd73d 78{
2ed70edf 79 // Dummy constructor
2188af53 80 cout<<"AliAnalysisTaskCumulants::AliAnalysisTaskCumulants()"<<endl;
81}
82
83//================================================================================================================
924fafb7 84
2ed70edf 85void AliAnalysisTaskCumulants::UserCreateOutputObjects()
924fafb7 86{
2ed70edf 87 // Called at every worker node to initialize
88 cout<<"AliAnalysisTaskCumulants::UserCreateOutputObjects()"<<endl;
2188af53 89
2ed70edf 90 // Analyser:
91 fGFC = new AliFlowAnalysisWithCumulants();
e5e75b58 92
2ed70edf 93 // Weights:
e5e75b58 94 if(fUseWeights)
95 {
2ed70edf 96 // Pass the flags to class:
97 if(fUsePhiWeights) fGFC->SetUsePhiWeights(fUsePhiWeights);
98 if(fUsePtWeights) fGFC->SetUsePtWeights(fUsePtWeights);
99 if(fUseEtaWeights) fGFC->SetUseEtaWeights(fUseEtaWeights);
100 // Get data from input slot #1 which is used for weights:
e5e75b58 101 if(GetNinputs()==2)
102 {
103 fListWeights = (TList*)GetInputData(1);
104 }
2ed70edf 105 // Pass the list with weights to class:
106 if(fListWeights) fGFC->SetWeightsList(fListWeights);
e5e75b58 107 }
2188af53 108
2ed70edf 109 fGFC->Init();
110
111 if(fGFC->GetHistList())
2188af53 112 {
2ed70edf 113 fListHistos = fGFC->GetHistList();
2188af53 114 //fListHistos->Print();
2ed70edf 115 } else
116 {
117 Printf("ERROR: Could not retrieve histogram list (GFC, Task::UserCreateOutputObjects()) !!!!");
118 }
61e0c8c0 119
120 PostData(1,fListHistos);
2ed70edf 121
122} // end of void AliAnalysisTaskCumulants::UserCreateOutputObjects()
2188af53 123
124//================================================================================================================
924fafb7 125
2ed70edf 126void AliAnalysisTaskCumulants::UserExec(Option_t *)
924fafb7 127{
2ed70edf 128 // Main loop (called for each event)
7183fe85 129 fEvent = dynamic_cast<AliFlowEventSimple*>(GetInputData(0));
aaebd73d 130
2ed70edf 131 // Generating function cumulants (GFC):
7183fe85 132 if(fEvent)
133 {
2ed70edf 134 fGFC->Make(fEvent);
135 } else
136 {
137 cout<<"WARNING: No input data (GFC, Task::UserExec()) !!!!"<<endl;
138 }
7183fe85 139
2ed70edf 140 PostData(1,fListHistos);
141
142} // end of void AliAnalysisTaskCumulants::UserExec(Option_t *)
924fafb7 143
2188af53 144//================================================================================================================
145
924fafb7 146void AliAnalysisTaskCumulants::Terminate(Option_t *)
aaebd73d 147{
2ed70edf 148 // Accessing the output list which contains the merged 2D and 3D profiles from all worker nodes
149 fListHistos = (TList*)GetOutputData(1);
2188af53 150 //fListHistos->Print();
151
2ed70edf 152 fGFC = new AliFlowAnalysisWithCumulants();
fd46c3dd 153
2188af53 154 if(fListHistos)
52021ae2 155 {
2ed70edf 156 fGFC->GetOutputHistograms(fListHistos);
157 fGFC->Finish();
158 PostData(1,fListHistos);
fd46c3dd 159 } else
160 {
2ed70edf 161 cout<<"WARNING: histogram list pointer is empty (GFC, Task::Terminate()) !!!!"<<endl;
fd46c3dd 162 cout<<endl;
163 }
2ed70edf 164
165} // end of void AliAnalysisTaskCumulants::Terminate(Option_t *)
aaebd73d 166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
924fafb7 182
924fafb7 183
924fafb7 184
924fafb7 185
c75fdbdc 186