]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGUD/selectors/trigger/AliTriggerTask.cxx
Compilation with Root6: TH1::GetXaxis returns now const TAxis*
[u/mrichter/AliRoot.git] / PWGUD / selectors / trigger / AliTriggerTask.cxx
1 /* $Id: AliTriggerTask.cxx 35782 2009-10-22 11:54:31Z jgrosseo $ */
2
3 #include "AliTriggerTask.h"
4
5 #include <TCanvas.h>
6 #include <TFile.h>
7 #include <TChain.h>
8 #include <TH1F.h>
9 #include <TH2F.h>
10 #include <TH3F.h>
11
12 #include <AliLog.h>
13 #include <AliESDEvent.h>
14 #include <AliHeader.h>
15 #include <AliAnalysisManager.h>
16 #include <AliESDInputHandler.h>
17 #include <AliESDHeader.h>
18 #include <AliPhysicsSelection.h>
19
20 ClassImp(AliTriggerTask)
21
22 AliTriggerTask::AliTriggerTask(const char* opt) :
23   AliAnalysisTask("AliTriggerTask", ""),
24   fESD(0),
25   fOutput(0),
26   fOption(opt),
27   fStartTime(0),
28   fEndTime(0),
29   fUseOrbits(kFALSE),
30   fFirstOrbit(0),
31   fLastOrbit(0),
32   fNTriggers(0),
33   fTriggerList(0),
34   fNTriggerClasses(0),
35   fTriggerClassesList(0),
36   fStats(0),
37   fPhysicsSelection(0)
38 {
39   //
40   // Constructor. Initialization of pointers
41   //
42
43   // Define input and output slots here
44   DefineInput(0, TChain::Class());
45   DefineOutput(0, TList::Class());
46   
47   fNTriggers = 14;
48   
49   static AliTriggerAnalysis::Trigger triggerList[] = { AliTriggerAnalysis::kAcceptAll, AliTriggerAnalysis::kFPANY, AliTriggerAnalysis::kMB1, AliTriggerAnalysis::kMB2, AliTriggerAnalysis::kMB3, AliTriggerAnalysis::kSPDGFO, AliTriggerAnalysis::kSPDGFOBits, AliTriggerAnalysis::kV0A, AliTriggerAnalysis::kV0C, AliTriggerAnalysis::kZDC, AliTriggerAnalysis::kZDCA, AliTriggerAnalysis::kZDCC, AliTriggerAnalysis::kFMDA, AliTriggerAnalysis::kFMDC };
50   fTriggerList = triggerList;
51   
52   fNTriggerClasses = 4;
53   static const char* triggerClassesList[] = { "CINT1B-ABCE-NOPF-ALL", "CINT1C-ABCE-NOPF-ALL", "CINT1A-ABCE-NOPF-ALL", "CINT1-E-NOPF-ALL" };
54   fTriggerClassesList = triggerClassesList;
55   
56   fStats = new TH1**[fNTriggerClasses];
57   for (Int_t i=0; i<fNTriggerClasses; i++)
58     fStats[i] = new TH1*[fNTriggers];
59 }
60
61 AliTriggerTask::~AliTriggerTask()
62 {
63   //
64   // Destructor
65   //
66
67   // histograms are in the output list and deleted when the output
68   // list is deleted by the TSelector dtor
69
70   if (fOutput) {
71     delete fOutput;
72     fOutput = 0;
73   }
74 }
75
76 //________________________________________________________________________
77 void AliTriggerTask::ConnectInputData(Option_t *)
78 {
79   // Connect ESD
80   // Called once
81
82   Printf("AliTriggerTask::ConnectInputData called");
83
84   AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
85
86   if (!esdH) {
87     Printf("ERROR: Could not get ESDInputHandler");
88   } else {
89     fESD = esdH->GetEvent();
90     
91     TString branches("AliESDHeader Vertex AliMultiplicity ALIESDVZERO ALIESDZDC FMD");
92
93     // Enable only the needed branches
94     esdH->SetActiveBranches(branches);
95   }
96 }
97
98 void AliTriggerTask::CreateOutputObjects()
99 {
100   // create result objects and add to output list
101
102   Printf("AliTriggerTask::CreateOutputObjects");
103
104   fOutput = new TList;
105   fOutput->SetOwner();
106   
107   if (fStartTime == fEndTime)
108     AliWarning("Start and endtime not set. Automatic binning will be used. This does not work in parallel systems");
109
110   Int_t nBins = 1000;
111   if (fEndTime - fStartTime > 0)
112     nBins = fEndTime - fStartTime + 1;
113   if (nBins > 10000)
114     nBins = 10000;
115   
116   Int_t start = 0;
117   Int_t end = fEndTime - fStartTime;
118   
119   if (fUseOrbits)
120   {
121     start = fStartTime;
122     end = fEndTime;
123   }
124   
125   for (Int_t j=0; j<fNTriggerClasses; j++)
126   {
127     for (Int_t i=0; i<fNTriggers; i++)
128     {
129       fStats[j][i] = new TH1F(Form("fStats_%d_%d", j, i), Form("%s %s;%s;counts", fTriggerClassesList[j], AliTriggerAnalysis::GetTriggerName(fTriggerList[i]), (fUseOrbits) ? "orbit number" : "time"), nBins, start - 0.5, end + 0.5);
130       fOutput->Add(fStats[j][i]);
131     }
132   }
133     
134   fFirstOrbit = new TParameter<Long_t> ("fFirstOrbit", 0);
135   fLastOrbit = new TParameter<Long_t> ("fLastOrbit", 0);
136   fOutput->Add(fFirstOrbit);
137   fOutput->Add(fLastOrbit);
138   
139   fOutput->Add(fPhysicsSelection);
140   
141   AliLog::SetClassDebugLevel("AliPhysicsSelection", AliLog::kDebug);
142 }
143
144 void AliTriggerTask::Exec(Option_t*)
145 {
146   // process the event
147
148   // post the data already here
149   PostData(0, fOutput);
150
151   if (!fESD)
152   {
153     AliError("ESD branch not available");
154     return;
155   }
156   
157   // fill histograms
158   fPhysicsSelection->IsCollisionCandidate(fESD);
159   
160   // check event type (should be PHYSICS = 7)
161   AliESDHeader* esdHeader = fESD->GetHeader();
162   if (!esdHeader)
163   {
164     Printf("ERROR: esdHeader could not be retrieved");
165     return;
166   }
167   
168   UInt_t eventType = esdHeader->GetEventType();
169   if (eventType != 7)
170   {
171     Printf("Skipping event because it is of type %d", eventType);
172     return;
173   }
174   
175   // TODO select on hardware trigger for histograms...
176   Int_t triggerClass = 0;
177   while (triggerClass < fNTriggerClasses && !fESD->IsTriggerClassFired(fTriggerClassesList[triggerClass]))
178     triggerClass++;
179   if (triggerClass == fNTriggerClasses)
180   {
181     Printf("Unknown trigger class %s. Skipping event", fESD->GetFiredTriggerClasses().Data());
182     return;
183   }
184   
185   Long64_t timeStamp = 0;
186   if (fUseOrbits)
187   {
188     timeStamp = fESD->GetBunchCrossNumber();
189     timeStamp += (Long64_t) 3564 * (fESD->GetOrbitNumber() + fESD->GetPeriodNumber() * 16777215);
190     timeStamp = (Long64_t) (25e-9 * timeStamp);
191     timeStamp -= fStartTime;
192   }
193   else
194     timeStamp = fESD->GetTimeStamp() - fStartTime;
195     
196   //Printf("%d", timeStamp);
197   
198   //Annalisa Time (s) = 1440*period + 88*10-6 * orbit + 25*10-9 *bc
199   
200   for (Int_t i = 0; i < fNTriggers; i++)
201   {
202     Bool_t triggered = fPhysicsSelection->GetTriggerAnalysis()->IsOfflineTriggerFired(fESD, fTriggerList[i]);
203     if (triggered)
204       fStats[triggerClass][i]->Fill(timeStamp);
205     //Printf("%s: %d", AliTriggerAnalysis::GetTriggerName(fTriggerList[i]), triggered);
206   }
207   
208   if (fFirstOrbit->GetVal() == 0)
209     fFirstOrbit->SetVal(timeStamp);
210   else
211     fFirstOrbit->SetVal(TMath::Min(fFirstOrbit->GetVal(), (Long_t) timeStamp));
212   
213   fLastOrbit->SetVal(TMath::Max(fLastOrbit->GetVal(), (Long_t) timeStamp));
214 }
215
216 void AliTriggerTask::Terminate(Option_t *)
217 {
218   // The Terminate() function is the last function to be called during
219   // a query. It always runs on the client, it can be used to present
220   // the results graphically or save the results to file.
221
222   fOutput = dynamic_cast<TList*> (GetOutputData(0));
223   if (!fOutput)
224     Printf("ERROR: fOutput not available");
225     
226   //fOutput->Print();
227
228   if (fOutput)
229   {
230     for (Int_t j=0; j<fNTriggerClasses; j++)
231       for (Int_t i=0; i<fNTriggers; i++)
232         fStats[j][i] = dynamic_cast<TH1*> (fOutput->FindObject(Form("fStats_%d_%d", j, i)));
233     fPhysicsSelection = dynamic_cast<AliPhysicsSelection*> (fOutput->FindObject("AliPhysicsSelection"));
234   }
235
236   TFile* fout = new TFile("trigger.root", "RECREATE");
237
238   for (Int_t j=0; j<fNTriggerClasses; j++)
239     for (Int_t i=0; i<fNTriggers; i++)
240       if (fStats[j][i])
241         fStats[j][i]->Write();
242         
243   if (fPhysicsSelection)
244   {
245     fPhysicsSelection->SaveHistograms("physics_selection");
246     fPhysicsSelection->Print();
247   }
248     
249   if (fFirstOrbit)
250     fFirstOrbit->Dump();
251   if (fLastOrbit)
252     fLastOrbit->Dump();
253     
254   fout->Write();
255   fout->Close();
256   
257   Int_t nX = (Int_t) TMath::Sqrt(fNTriggers);
258   Int_t nY = nX;
259   
260   while (nX * nY < fNTriggers)
261   {
262     if (nX == nY)
263       nX++;
264     else
265       nY++;
266   }
267   
268   TCanvas* c = new TCanvas("c", "c", 800, 800);
269   c->Divide(nX, nY);
270   
271   Printf("+++++++++ TRIGGER STATS:");
272
273   Int_t triggerClass = 0;
274
275   Int_t base = 1;
276   if (fStats[triggerClass][0])
277     base = (Int_t) fStats[triggerClass][0]->Integral();
278
279   Int_t length = fEndTime - fStartTime;
280   
281   for (Int_t i=0; i<fNTriggers; i++)
282     if (fStats[triggerClass][i])
283     {
284       c->cd(i+1);
285       fStats[triggerClass][i]->Draw();
286       Printf("%s: %d triggers | %f %% of all triggered | Rate: %f Hz", AliTriggerAnalysis::GetTriggerName(fTriggerList[i]), (UInt_t) fStats[triggerClass][i]->Integral(), 100.0 * fStats[triggerClass][i]->Integral() / base, (length > 0) ? (fStats[triggerClass][i]->Integral() / length) : -1);
287     }
288
289   Printf("Writting result to trigger.root");
290 }