]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/AliRsnAnalysisSE.cxx
Package revised - New AnalysisTask's - Added more functions
[u/mrichter/AliRoot.git] / PWG2 / RESONANCES / AliRsnAnalysisSE.cxx
1 //
2 // Class AliRsnAnalysisSE
3 //
4 // TODO
5 //
6 // authors: Martin Vala (martin.vala@cern.ch)
7 //          Alberto Pulvirenti (alberto.pulvirenti@ct.infn.it)
8 //
9
10 #include <TSystem.h>
11 #include <TFile.h>
12 #include <TFolder.h>
13
14 #include "AliLog.h"
15
16 #include "AliAnalysisManager.h"
17 #include "AliRsnPairMgr.h"
18 #include "AliRsnEventBuffer.h"
19
20 #include "AliMCEventHandler.h"
21 #include "AliESDEvent.h"
22
23 #include "AliRsnAnalysisSE.h"
24
25 ClassImp(AliRsnAnalysisSE)
26
27 //________________________________________________________________________
28 AliRsnAnalysisSE::AliRsnAnalysisSE(const char * name)
29     : AliRsnAnalysisTaskSEBase(name),
30     fPairMgrs(0),fOutList(0x0),fRsnEventBuffer(0x0),
31     fNumOfEventsInBuffer(1000)
32 {
33 //=========================================================
34 // Default constructor
35 //=========================================================
36
37   InitIOVars();
38   DefineOutput(1, TList::Class());
39 }
40
41 //________________________________________________________________________
42 AliRsnAnalysisSE::~AliRsnAnalysisSE()
43 {
44 //=========================================================
45 // Destructor
46 //=========================================================
47 }
48
49 //________________________________________________________________________
50 void AliRsnAnalysisSE::InitIOVars()
51 {
52 //=========================================================
53 // Init input output values
54 //=========================================================
55
56   AliDebug(AliLog::kDebug, "<-");
57   AliRsnAnalysisTaskSEBase::InitIOVars();
58
59   fRsnEventBuffer = 0;
60   fOutList = 0;
61
62   AliDebug(AliLog::kDebug, "->");
63 }
64
65 //________________________________________________________________________
66 void AliRsnAnalysisSE::UserCreateOutputObjects()
67 {
68 //=========================================================
69 // UserCreateOutputObjects() of AliAnalysisTaskSE
70 //=========================================================
71
72   AliDebug(AliLog::kDebug, "<-");
73 //   fPID.DumpPriors();
74   OpenFile ( 0 );
75   fOutList = new TList();
76   fOutList->SetOwner();
77   AliRsnPair *def=0;
78   AliRsnPairMgr *mgr=0;
79   TList *listTmp;
80   for (Int_t iMgr=0 ;iMgr< fPairMgrs.GetEntries();iMgr++)
81   {
82     mgr = (AliRsnPairMgr *) fPairMgrs.At(iMgr);
83     if (!mgr) continue;
84     listTmp = new TList();
85     listTmp->SetName(mgr->GetName());
86     for (Int_t i=0;i< mgr->GetPairs()->GetEntriesFast();i++)
87     {
88       def = (AliRsnPair *) mgr->GetPairs()->At(i);
89       if (def)
90       {
91         def->Init();
92         //listTmp->Add(def->GenerateHistograms(mgr->GetName()));
93         def->GenerateHistograms(mgr->GetName(), listTmp);
94         //def->Print();
95       }
96     }
97     fOutList->Add(listTmp);
98   }
99
100   fRsnEventBuffer = new AliRsnEventBuffer(100);
101 //   fRsnEventBuffer = new AliRsnEventBuffer ( 10000 ,kFALSE );
102   AliDebug(AliLog::kDebug, "->");
103
104 }
105
106 //________________________________________________________________________
107 void AliRsnAnalysisSE::UserExec(Option_t *)
108 {
109 //=========================================================
110 // UserExec() of AliAnalysisTaskSE
111 //=========================================================
112
113   if (fEntry++%1000==0)
114     AliInfo(Form("Event %d",-1));
115
116   AliRsnEvent *curEvent = GetRsnEventFromInputType();
117   if (!curEvent) return;
118
119   ProcessEventAnalysis(curEvent);
120   PostEventProcess();
121
122   PostData(1, fOutList);
123 }
124
125 //________________________________________________________________________
126 void AliRsnAnalysisSE::Terminate(Option_t *)
127 {
128 //=========================================================
129 // Terminate() of AliAnalysisTask
130 //=========================================================
131
132   fOutList = dynamic_cast<TList*>(GetOutputData(1));
133   if (!fOutList) { AliError(" fOutList not available"); return; }
134   fOutList->Print();
135 }
136
137 //________________________________________________________________________
138 void AliRsnAnalysisSE::ProcessEventAnalysis(AliRsnEvent *curEvent)
139 {
140 //=========================================================
141 // Process of one event
142 //=========================================================
143
144
145   // Adds event to Event Buffer
146   fRsnEventBuffer->AddEvent(curEvent);
147
148   // loop over all Pair managers
149   AliRsnPairMgr *mgr=0;
150   for (Int_t iMgr=0 ;iMgr< fPairMgrs.GetEntries();iMgr++)
151   {
152     mgr = (AliRsnPairMgr *) fPairMgrs.At(iMgr);
153     AliRsnPair *pair=0;
154     for (Int_t i=0;i< mgr->GetPairs()->GetEntriesFast();i++)
155     {
156       pair = (AliRsnPair *) mgr->GetPairs()->At(i);
157       pair->ProcessPair(fRsnEventBuffer);
158     }
159   }
160 }
161
162 //________________________________________________________________________
163 void AliRsnAnalysisSE::PostEventProcess(const Short_t & index)
164 {
165 //=========================================================
166 // Post process of one event
167 //=========================================================
168
169   switch (fInputType[index])
170   {
171     case kAOD:
172       break;
173     case kESD:
174       break;
175     case kESDMC:
176       break;
177     case kMC:
178       break;
179     case kRSN:
180     {
181       if (fRsnEventBuffer->GetDeleteBufferWhenReset() == kFALSE)
182       {
183         fRSN[index] = (AliRsnEvent*) fRsnEventBuffer->GetNextEvent();
184         SetBranchAddress(0 , "rsnEvents", &fRSN[index]);
185       }
186       break;
187     }
188     default:
189       break;
190   }
191
192 }
193
194 void  AliRsnAnalysisSE::AddPairMgr(AliRsnPairMgr * pairmgr)
195 {
196   fPairMgrs.Add(pairmgr);
197 }