]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/AliRsnVAnalysisTask.cxx
Added a new cut on comparison of daughters' momenta (needed for Deltas)
[u/mrichter/AliRoot.git] / PWG2 / RESONANCES / AliRsnVAnalysisTask.cxx
1 //
2 // Class AliRsnVAnalysisTask
3 //
4 // Virtual Class derivated from AliAnalysisTaskSE which will be base class
5 // for all RSN SE tasks
6 //
7 // authors: Martin Vala (martin.vala@cern.ch)
8 //          Alberto Pulvirenti (alberto.pulvirenti@ct.infn.it)
9 //
10
11 #include <Riostream.h>
12
13 #include "AliESDEvent.h"
14 #include "AliMCEvent.h"
15 #include "AliAODEvent.h"
16 #include "AliAnalysisManager.h"
17 #include "AliMCEventHandler.h"
18 #include "AliMultiInputEventHandler.h"
19 #include "AliMixInputEventHandler.h"
20
21
22 #include "AliRsnEvent.h"
23 #include "AliRsnTarget.h"
24
25 #include "AliRsnVAnalysisTask.h"
26
27 ClassImp(AliRsnVAnalysisTask)
28
29 //_____________________________________________________________________________
30 AliRsnVAnalysisTask::AliRsnVAnalysisTask
31 (const char *name, Bool_t mcOnly) :
32    AliAnalysisTaskSE(name),
33    fLogType(AliLog::kInfo),
34    fLogClassesString(""),
35    fIsMixing(kFALSE),
36    fMCOnly(mcOnly),
37    fInfoList(0x0),
38    fTaskInfo(name),
39    fMixedEH(0),
40    fUseMixingRange(kTRUE)
41 {
42 //
43 // Default constructor.
44 // Define the output slot for histograms.
45 //
46
47    DefineOutput(1, TList::Class());
48    DefineOutput(2, TList::Class());
49    for (Int_t i = 0; i < 2; i++) {
50       fESDEvent[i] = 0;
51       fMCEvent[i] = 0;
52       fAODEventIn[i] = 0;
53       fAODEventOut[i] = 0;
54    }
55 }
56
57 //_____________________________________________________________________________
58 AliRsnVAnalysisTask::AliRsnVAnalysisTask(const AliRsnVAnalysisTask& copy) :
59    AliAnalysisTaskSE(copy),
60    fLogType(copy.fLogType),
61    fLogClassesString(copy.fLogClassesString),
62    fIsMixing(copy.fIsMixing),
63    fMCOnly(copy.fMCOnly),
64    fInfoList(0x0),
65    fTaskInfo(copy.fTaskInfo),
66    fMixedEH(copy.fMixedEH),
67    fUseMixingRange(copy.fUseMixingRange)
68 {
69 //
70 // Copy constructor.
71 // Defined for coding conventions compliance but never used.
72 //
73    AliDebug(AliLog::kDebug + 2, "<-");
74    for (Int_t i = 0; i < 2; i++) {
75       fESDEvent[i] = copy.fESDEvent[i];
76       fMCEvent[i] = copy.fMCEvent[i];
77       fAODEventIn[i] = copy.fAODEventIn[i];
78       fAODEventOut[i] = copy.fAODEventOut[i];
79       fRsnEvent[i] = copy.fRsnEvent[i];
80    }
81    AliDebug(AliLog::kDebug + 2, "->");
82
83 }
84
85 //_____________________________________________________________________________
86 void AliRsnVAnalysisTask::LocalInit()
87 {
88 //
89 // Local initialization.
90 // Defines the debug message level and calls the mother class LocalInit().
91 //
92
93    AliAnalysisTaskSE::LocalInit();
94    SetDebugForAllClasses();
95 }
96
97 //_____________________________________________________________________________
98 Bool_t AliRsnVAnalysisTask::UserNotify()
99 {
100 //
101 // Calls the mother class Notify()
102 //
103
104    return AliAnalysisTaskSE::UserNotify();
105 }
106
107 //_____________________________________________________________________________
108 void AliRsnVAnalysisTask::ConnectInputData(Option_t *opt)
109 {
110 //
111 // Connect input data, which consist in initializing properly
112 // the pointer to the input event, which is dynamically casted
113 // to all available types, and this allows to know its type.
114 // Calls also the mother class omonyme method.
115 //
116
117    AliAnalysisTaskSE::ConnectInputData(opt);
118
119    // get AliESDEvent and, if successful
120    // retrieve the corresponding MC if exists
121    fESDEvent[0] = dynamic_cast<AliESDEvent *>(fInputEvent);
122    if (fESDEvent[0]) {
123       fMCEvent[0] = (AliMCEvent*) MCEvent();
124       AliInfo(Form("Input event is of type ESD   (%p)", fESDEvent[0]));
125       if (fMCEvent[0]) AliInfo(Form("Input has an associated MC (%p)", fMCEvent[0]));
126    }
127
128    // get AliAODEvent from input and, if successful
129    // it will contain both the reconstructed and MC informations
130    fAODEventIn[0] = dynamic_cast<AliAODEvent *>(fInputEvent);
131    if (fAODEventIn[0]) {
132       AliInfo(Form("Input event if of type native AOD (%p)", fAODEventIn[0]));
133    }
134
135    // get AliAODEvent from output of previous task
136    fAODEventOut[0] = dynamic_cast<AliAODEvent *>(AODEvent());
137    if (fAODEventOut[0]) {
138       AliInfo(Form("Input event if of type produced AOD from previous step (%p)", fAODEventOut[0]));
139    }
140 }
141
142 //_____________________________________________________________________________
143 void AliRsnVAnalysisTask::UserCreateOutputObjects()
144 {
145 //
146 // Creates and links to task all output objects.
147 // Does explicitly the initialization for the event info class,
148 // and then calls the customized function which must be overloaded
149 // in the applications of this base class.
150 //
151
152    SetDebugForAllClasses();
153
154    // set event info outputs
155    fInfoList = new TList();
156    fInfoList->SetOwner();
157    fTaskInfo.GenerateInfoList(fInfoList);
158
159    // create customized outputs
160    RsnUserCreateOutputObjects();
161
162    PostData(1, fInfoList);
163 }
164
165 //_____________________________________________________________________________
166 void AliRsnVAnalysisTask::UserExec(Option_t* opt)
167 {
168 //
169 // Prepares for execution, setting the correct pointers of the
170 // RSN package event interface, which will point to the not NULL
171 // objects obtained from dynamic-casts called in ConnectInputData().
172 //
173    if (!IsMixing()) {
174
175
176       if (fMCOnly && fMCEvent) {
177          fRsnEvent[0].SetRef(fMCEvent[0]);
178          fRsnEvent[0].SetRefMC(fMCEvent[0]);
179       } else if (fESDEvent[0]) {
180          fRsnEvent[0].SetRef(fESDEvent[0]);
181          fRsnEvent[0].SetRefMC(fMCEvent[0]);
182       } else if (fAODEventOut) {
183          fRsnEvent[0].SetRef(fAODEventOut[0]);
184          fRsnEvent[0].SetRefMC(fAODEventOut[0]);
185       } else if (fAODEventIn) {
186          fRsnEvent[0].SetRef(fAODEventIn[0]);
187          fRsnEvent[0].SetRefMC(fAODEventIn[0]);
188       } else {
189          AliError("Unknown input event format. Skipping");
190          return;
191       }
192
193       // since this class is for single-event analysis
194       // both static pointers of AliRsnEvent class
195       // will point to the same unique datamember
196       AliRsnEvent::SetCurrentEvent1(&fRsnEvent[0], fEntry);
197       AliRsnEvent::SetCurrentEvent2(&fRsnEvent[0], fEntry);
198       AliRsnTarget::SwitchToFirst();
199
200       // call event preprocessing...
201       Bool_t preCheck = EventProcess();
202       // ...then fill the information object and print informations...
203       fTaskInfo.FillInfo();
204       fTaskInfo.PrintInfo(fTaskInfo.GetNumerOfEventsProcessed());
205       // ...and return if event did not pass selections
206       if (!preCheck) {
207          AliDebug(AliLog::kDebug, "Event preprocessing has failed. Skipping event");
208          return;
209       }
210
211
212       // call customized implementation for execution
213       RsnUserExec(opt);
214    }
215    // post outputs for the info object
216    // (eventually others are done in the derived classes)
217    PostData(1, fInfoList);
218 }
219
220 void AliRsnVAnalysisTask::UserExecMix(Option_t* option)
221 {
222    AliDebug(AliLog::kDebug + 2, "<-");
223
224
225    if (!IsMixing()) return;
226
227    SetupMixingEvents();
228
229    if (!fMixedEH) return;
230
231    if (fMCOnly && fMCEvent[0]) {
232       fRsnEvent[0].SetRef(fMCEvent[0]);
233       fRsnEvent[0].SetRefMC(fMCEvent[0]);
234       fRsnEvent[1].SetRef(fMCEvent[1]);
235       fRsnEvent[1].SetRefMC(fMCEvent[1]);
236    } else if (fESDEvent[0]) {
237       fRsnEvent[0].SetRef(fESDEvent[0]);
238       fRsnEvent[0].SetRefMC(fMCEvent[0]);
239       fRsnEvent[1].SetRef(fESDEvent[1]);
240       fRsnEvent[1].SetRefMC(fMCEvent[1]);
241    } else if (fAODEventOut) {
242       fRsnEvent[0].SetRef(fAODEventOut[0]);
243       fRsnEvent[0].SetRefMC(fAODEventOut[0]);
244       fRsnEvent[1].SetRef(fAODEventOut[1]);
245       fRsnEvent[1].SetRefMC(fAODEventOut[1]);
246    } else if (fAODEventIn) {
247       fRsnEvent[0].SetRef(fAODEventIn[0]);
248       fRsnEvent[0].SetRefMC(fAODEventIn[0]);
249       fRsnEvent[1].SetRef(fAODEventIn[1]);
250       fRsnEvent[1].SetRefMC(fAODEventIn[1]);
251    } else {
252       AliError("NO ESD or AOD object!!! Skipping ...");
253       return;
254    }
255    // since this class is for single-event analysis
256    // both static pointers of AliRsnEvent class
257    // will point to the same unique datamember
258
259
260    AliRsnEvent::SetCurrentEvent1(&fRsnEvent[0], fMixedEH->CurrentEntryMain());
261    AliRsnEvent::SetCurrentEvent2(&fRsnEvent[1], fMixedEH->CurrentEntryMix());
262    AliRsnTarget::SwitchToFirst();
263
264    // call event preprocessing...
265    Bool_t preCheck = EventProcess();
266    // ...then fill the information object and print informations...
267    fTaskInfo.FillInfo();
268    fTaskInfo.PrintInfo(fTaskInfo.GetNumerOfEventsProcessed());
269    // ...and return if event did not pass selections
270    if (!preCheck) {
271       AliDebug(AliLog::kDebug, "Event preprocessing has failed. Skipping event");
272       return;
273    }
274
275    RsnUserExecMix(option);
276    AliDebug(AliLog::kDebug + 2, "->");
277 }
278
279
280 //_____________________________________________________________________________
281 void AliRsnVAnalysisTask::Terminate(Option_t* opt)
282 {
283 //
284 // Termination routines.
285 // Stores all histograms (after checking they exist)
286 // and includes to the TList all task informations.
287 //
288
289    AliAnalysisTask::Terminate();
290
291    TList* list  = dynamic_cast<TList*>(GetOutputData(1));
292    if (!list) {
293       AliError(Form("At end of analysis, fOutList is %p", list));
294       return;
295    }
296
297    RsnTerminate(opt);
298
299    TH1I *hEventInfo = (TH1I*) list->FindObject(fTaskInfo.GetEventHistogramName());
300    if (!hEventInfo) {
301       AliError(Form("hEventInfo is %p", hEventInfo));
302       return;
303    }
304    AliInfo(Form("=== %s ==================", GetName()));
305    AliInfo(Form("Number Of Events Processed : %10lld", (Long64_t)hEventInfo->Integral()));
306    AliInfo(Form("Number Of Events Accepted  : %10lld", (Long64_t)hEventInfo->GetBinContent(2)));
307    AliInfo(Form("Number Of Events Skipped   : %10lld", (Long64_t)hEventInfo->GetBinContent(1)));
308    AliInfo(Form("=== end %s ==============", GetName()));
309
310    AliDebug(AliLog::kDebug + 2, "->");
311 }
312
313 //_____________________________________________________________________________
314 void AliRsnVAnalysisTask::RsnUserCreateOutputObjects()
315 {
316 //
317 // Define here all instructions to create output objects.
318 // This method will be called inside the "UserCreateOutputObjects"
319 // in the used task.
320 //
321
322    AliWarning("Implement this in derived classes");
323 }
324
325 //_____________________________________________________________________________
326 void AliRsnVAnalysisTask::RsnUserExec(Option_t*)
327 {
328 //
329 //
330 //
331
332    AliWarning("Implement this in derived classes");
333 }
334
335 void AliRsnVAnalysisTask::RsnUserExecMix(Option_t*)
336 {
337    //
338    //
339    //
340
341    AliWarning("Implement this in derived classes");
342 }
343
344 //_____________________________________________________________________________
345 void AliRsnVAnalysisTask::RsnTerminate(Option_t*)
346 {
347 //
348 // Overload this to add additional termination operations
349 //
350
351    AliWarning("Implement this in derived classes");
352 }
353
354 //_____________________________________________________________________________
355 Bool_t AliRsnVAnalysisTask::EventProcess()
356 {
357 //
358 // Performs some pre-processing of current event,
359 // which is useful for all the operations which
360 // need to be done only once for each event.
361 //
362
363    // if not using mixing cuts return kTRUE
364    if (!IsUsingMixingRange()) return kTRUE;
365
366    // cut if event was in range of mixing cuts
367    AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
368    AliMultiInputEventHandler *inEvHMain = dynamic_cast<AliMultiInputEventHandler *>(mgr->GetInputEventHandler());
369    if (inEvHMain) {
370       fMixedEH = dynamic_cast<AliMixInputEventHandler *>(inEvHMain->GetFirstMultiInputHandler());
371       if (fMixedEH) {
372          if (fMixedEH->CurrentBinIndex() < 0) return kFALSE;
373       }
374
375    }
376
377    // in this case, return always a success
378    return kTRUE;
379 }
380
381 void AliRsnVAnalysisTask::SetupMixingEvents()
382 {
383    // mixing setup
384
385    AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
386    AliMultiInputEventHandler *inEvHMain = dynamic_cast<AliMultiInputEventHandler *>(mgr->GetInputEventHandler());
387    if (inEvHMain) {
388       fMixedEH = dynamic_cast<AliMixInputEventHandler *>(inEvHMain->GetFirstMultiInputHandler());
389       if (fMixedEH) {
390          AliMultiInputEventHandler *inEvHMainTmpMix = dynamic_cast<AliMultiInputEventHandler *>(fMixedEH->InputEventHandler(0));
391          if (!inEvHMainTmpMix) return;
392          AliInputEventHandler *inEvHMixTmp = dynamic_cast<AliInputEventHandler *>(inEvHMainTmpMix->GetFirstInputEventHandler());
393          AliMCEventHandler *inEvHMCMixTmp = dynamic_cast<AliMCEventHandler *>(inEvHMainTmpMix->GetFirstMCEventHandler());
394          if (!inEvHMixTmp) return;
395          fESDEvent[1] = dynamic_cast<AliESDEvent *>(inEvHMixTmp->GetEvent());
396          if (fESDEvent[1]) AliDebug(AliLog::kDebug, Form("Input is ESD (%p) MIXED", fESDEvent[1]));
397          // getting AliAODEvent from input
398          fAODEventIn[1] = dynamic_cast<AliAODEvent *>(inEvHMixTmp->GetEvent());
399          if (fAODEventIn[1]) AliDebug(AliLog::kDebug, Form("Input is AOD (%p) MIXED", fAODEventIn[1]));
400          // getting AliAODEvent if it is output from previous task (not supported)
401          fAODEventOut[1] = 0;
402
403          if (inEvHMCMixTmp) {
404             fMCEvent[1] = inEvHMCMixTmp->MCEvent();
405             if (fMCEvent[1]) AliDebug(AliLog::kDebug, Form("Input is ESDMC (%p) MIXED", fMCEvent[1]));
406          }
407
408       }
409
410    }
411 }
412
413 //_____________________________________________________________________________
414 void AliRsnVAnalysisTask::SetDebugForAllClasses()
415 {
416 //
417 // Set debug level for all classes for which it is required
418 //
419
420    TObjArray  *array = fLogClassesString.Tokenize(":");
421    TObjString *objStr;
422    TString     str;
423    Int_t       i, n = array->GetEntriesFast();
424
425    for (i = 0; i < n; i++) {
426       objStr = (TObjString*)array->At(i);
427       str    = objStr->GetString();
428       AliLog::SetClassDebugLevel(str.Data(), fLogType);
429       AliInfo(Form("Setting Debug to %s", str.Data()));
430    }
431 }
432