]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/AliOfflineTrigger.cxx
Major update of the CMake compilation:
[u/mrichter/AliRoot.git] / PWG0 / AliOfflineTrigger.cxx
1 /* $Id: AliOfflineTrigger.cxx 35782 2009-10-22 11:54:31Z jgrosseo $ */
2
3 /**************************************************************************
4  * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
5  *                                                                        *
6  * Author: The ALICE Off-line Project.                                    *
7  * Contributors are mentioned in the code where appropriate.              *
8  *                                                                        *
9  * Permission to use, copy, modify and distribute this software and its   *
10  * documentation strictly for non-commercial purposes is hereby granted   *
11  * without fee, provided that the above copyright notice appears in all   *
12  * copies and that both the copyright notice and this permission notice   *
13  * appear in the supporting documentation. The authors make no claims     *
14  * about the suitability of this software for any purpose. It is          *
15  * provided "as is" without express or implied warranty.                  *
16  **************************************************************************/
17
18 //-------------------------------------------------------------------------
19 //                      Implementation of   Class AliOfflineTrigger
20 //   This class provides offline triggers from data in the ESD
21 //   Origin: Jan Fiete Grosse-Oetringhaus, CERN
22 //-------------------------------------------------------------------------
23
24 #include <TH1F.h>
25 #include <TList.h>
26 #include <TIterator.h>
27
28 #include <AliOfflineTrigger.h>
29
30 #include <AliLog.h>
31
32 #include <AliESDEvent.h>
33
34 #include <AliMultiplicity.h>
35 #include <AliESDVZERO.h>
36 #include <AliESDZDC.h>
37 #include <AliESDFMD.h>
38
39 ClassImp(AliOfflineTrigger)
40
41 AliOfflineTrigger::AliOfflineTrigger() :
42   fSPDGFOThreshold(1),
43   fV0AThreshold(1),
44   fV0CThreshold(1),
45   fFMDLowCut(0.7),
46   fFMDHitCut(1.2),
47   fHistSPD(0),
48   fHistV0A(0),       
49   fHistV0C(0),    
50   fHistZDC(0),    
51   fHistFMDA(0),    
52   fHistFMDC(0),   
53   fHistFMDSingle(0),
54   fHistFMDSum(0)
55 {
56 }
57
58 void AliOfflineTrigger::EnableHistograms()
59 {
60   // creates the monitoring histograms
61   
62   fHistSPD = new TH1F("fHistSPD", "SPD GFO;number of fired chips;events", 1202, -1.5, 1200.5);
63   fHistV0A = new TH1F("fHistV0A", "V0A;number of BB triggers;events", 34, -1.5, 32.5);
64   fHistV0C = new TH1F("fHistV0C", "V0C;number of BB triggers;events", 34, -1.5, 32.5);
65   fHistZDC = new TH1F("fHistZDC", "ZDC;trigger bits;events", 8, -1.5, 6.5);
66   
67   // TODO check limits
68   fHistFMDA = new TH1F("fHistFMDA", "FMDA;combinations above threshold;events", 102, -1.5, 100.5);
69   fHistFMDC = new TH1F("fHistFMDC", "FMDC;combinations above threshold;events", 102, -1.5, 100.5);
70   fHistFMDSingle = new TH1F("fHistFMDSingle", "FMD single;multiplicity value;counts", 1000, 0, 10);
71   fHistFMDSum = new TH1F("fHistFMDSum", "FMD sum;multiplicity value;counts", 1000, 0, 10);
72 }
73
74 Bool_t AliOfflineTrigger::IsEventTriggered(const AliESDEvent* aEsd, AliPWG0Helper::Trigger trigger) const
75 {
76   // checks if an event has been triggered "offline"
77
78   UInt_t triggerNoFlags = (UInt_t) trigger % (UInt_t) AliPWG0Helper::kStartOfFlags;
79   
80   switch (triggerNoFlags)
81   {
82     case AliPWG0Helper::kAcceptAll:
83     {
84       return kTRUE;
85       break;
86     }
87     case AliPWG0Helper::kMB1:
88     {
89       if (SPDGFOTrigger(aEsd) || V0Trigger(aEsd, kASide) || V0Trigger(aEsd, kCSide))
90         return kTRUE;
91       break;
92     }
93     case AliPWG0Helper::kMB2:
94     {
95       if (SPDGFOTrigger(aEsd) && (V0Trigger(aEsd, kASide) || V0Trigger(aEsd, kCSide)))
96         return kTRUE;
97       break;
98     }
99     case AliPWG0Helper::kMB3:
100     {
101       if (SPDGFOTrigger(aEsd) && V0Trigger(aEsd, kASide) && V0Trigger(aEsd, kCSide))
102         return kTRUE;
103       break;
104     }
105     case AliPWG0Helper::kSPDGFO:
106     {
107       if (SPDGFOTrigger(aEsd))
108         return kTRUE;
109       break;
110     }
111     case AliPWG0Helper::kV0A:
112     {
113       if (V0Trigger(aEsd, kASide))
114         return kTRUE;
115       break;
116     }
117     case AliPWG0Helper::kV0C:
118     {
119       if (V0Trigger(aEsd, kCSide))
120         return kTRUE;
121       break;
122     }
123     case AliPWG0Helper::kZDC:
124     {
125       if (ZDCTrigger(aEsd, kASide) || ZDCTrigger(aEsd, kCentralBarrel) || ZDCTrigger(aEsd, kCSide))
126         return kTRUE;
127       break;
128     }
129     case AliPWG0Helper::kZDCA:
130     {
131       if (ZDCTrigger(aEsd, kASide))
132         return kTRUE;
133       break;
134     }
135     case AliPWG0Helper::kZDCC:
136     {
137       if (ZDCTrigger(aEsd, kCSide))
138         return kTRUE;
139       break;
140     }
141     case AliPWG0Helper::kFMDA:
142     {
143       if (FMDTrigger(aEsd, kASide))
144         return kTRUE;
145       break;
146     }
147     case AliPWG0Helper::kFMDC:
148     {
149       if (FMDTrigger(aEsd, kCSide))
150         return kTRUE;
151       break;
152     }
153     case AliPWG0Helper::kFPANY:
154     {
155       if (SPDGFOTrigger(aEsd) || V0Trigger(aEsd, kASide) || V0Trigger(aEsd, kCSide) || ZDCTrigger(aEsd, kASide) || ZDCTrigger(aEsd, kCentralBarrel) || ZDCTrigger(aEsd, kCSide) || FMDTrigger(aEsd, kASide) || FMDTrigger(aEsd, kCSide))
156         return kTRUE;
157       break;
158     }
159     default:
160     {
161       AliFatal(Form("Trigger type %d not implemented", triggerNoFlags));
162     }
163   }
164   
165   return kFALSE;
166 }
167
168 void AliOfflineTrigger::FillHistograms(const AliESDEvent* aEsd)
169 {
170   // fills the histograms with the info from the ESD
171   
172   fHistSPD->Fill(SPDFiredChips(aEsd));
173   
174   fHistV0A->Fill(V0BBTriggers(aEsd, kASide));  
175   fHistV0C->Fill(V0BBTriggers(aEsd, kCSide));
176   
177   AliESDZDC* zdcData = aEsd->GetESDZDC();
178   if (zdcData)
179   {
180     UInt_t quality = zdcData->GetESDQuality();
181     
182     // from Nora's presentation, general first physics meeting 16.10.09
183     static UInt_t zpc  = 0x20;
184     static UInt_t znc  = 0x10;
185     static UInt_t zem1 = 0x08;
186     static UInt_t zem2 = 0x04;
187     static UInt_t zpa  = 0x02;
188     static UInt_t zna  = 0x01;
189    
190     fHistZDC->Fill(1, quality & zna);
191     fHistZDC->Fill(2, quality & zpa);
192     fHistZDC->Fill(3, quality & zem2);
193     fHistZDC->Fill(4, quality & zem1);
194     fHistZDC->Fill(5, quality & znc);
195     fHistZDC->Fill(6, quality & zpc);
196   }
197   else
198   {
199     fHistZDC->Fill(-1);
200     AliError("AliESDZDC not available");
201   }
202   
203   fHistFMDA->Fill(FMDHitCombinations(aEsd, kASide, kTRUE));
204   fHistFMDC->Fill(FMDHitCombinations(aEsd, kCSide, kTRUE));
205 }
206
207 Int_t AliOfflineTrigger::SPDFiredChips(const AliESDEvent* aEsd) const
208 {
209   // returns the number of fired chips in the SPD
210   
211   const AliMultiplicity* mult = aEsd->GetMultiplicity();
212   if (!mult)
213   {
214     AliError("AliMultiplicity not available");
215     return -1;
216   }
217   return mult->GetNumberOfFiredChips(0) + mult->GetNumberOfFiredChips(1);
218 }
219
220 Bool_t AliOfflineTrigger::SPDGFOTrigger(const AliESDEvent* aEsd) const
221 {
222   // Returns if the SPD gave a global Fast OR trigger
223   
224   Int_t firedChips = SPDFiredChips(aEsd);
225   
226   if (firedChips >= fSPDGFOThreshold)
227     return kTRUE;
228   return kFALSE;
229 }
230
231 Int_t AliOfflineTrigger::V0BBTriggers(const AliESDEvent* aEsd, AliceSide side) const
232 {
233   // returns the number of BB triggers in V0A | V0C
234   
235   AliESDVZERO* v0Data = aEsd->GetVZEROData();
236   if (!v0Data)
237   {
238     AliError("AliESDVZERO not available");
239     return -1;
240   }
241   
242   Int_t count = 0;
243   for (Int_t i=0; i<32; i++)
244   {
245     if (side == kASide && v0Data->BBTriggerV0A(i))
246       count++;
247     if (side == kCSide && v0Data->BBTriggerV0C(i))
248       count++;
249   }
250   
251   return count;
252 }
253
254 Bool_t AliOfflineTrigger::V0Trigger(const AliESDEvent* aEsd, AliceSide side) const
255 {
256   // Returns if the V0 triggered
257   
258   Int_t count = V0BBTriggers(aEsd, side);
259   
260   if (side == kASide && count >= fV0AThreshold)
261     return kTRUE;
262   if (side == kCSide && count >= fV0CThreshold)
263     return kTRUE;
264   return kFALSE;
265 }
266
267 Bool_t AliOfflineTrigger::ZDCTrigger(const AliESDEvent* aEsd, AliceSide side) const
268 {
269   // Returns if ZDC triggered
270   
271   AliESDZDC* zdcData = aEsd->GetESDZDC();
272   if (!zdcData)
273   {
274     AliError("AliESDZDC not available");
275     return kFALSE;
276   }
277   
278   UInt_t quality = zdcData->GetESDQuality();
279   
280   // from Nora's presentation, general first physics meeting 16.10.09
281   static UInt_t zpc  = 0x20;
282   static UInt_t znc  = 0x10;
283   static UInt_t zem1 = 0x08;
284   static UInt_t zem2 = 0x04;
285   static UInt_t zpa  = 0x02;
286   static UInt_t zna  = 0x01;
287   
288   if (side == kASide && ((quality & zpa) || (quality & zna)))
289     return kTRUE;
290   if (side == kCentralBarrel && ((quality & zem1) || (quality & zem2)))
291     return kTRUE;
292   if (side == kCSide && ((quality & zpc) || (quality & znc)))
293     return kTRUE;
294   
295   return kFALSE;
296 }
297
298 Int_t AliOfflineTrigger::FMDHitCombinations(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHistograms) const
299 {
300   // returns number of hit combinations agove threshold
301   //
302   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
303
304   // Workaround for AliESDEvent::GetFMDData is not const!
305   const AliESDFMD* fmdData = (const_cast<AliESDEvent*>(aEsd))->GetFMDData();
306   if (!fmdData)
307   {
308     AliError("AliESDFMD not available");
309     return -1;
310   }
311
312   Int_t detFrom = (side == kASide) ? 1 : 3;
313   Int_t detTo   = (side == kASide) ? 2 : 3;
314
315   Int_t triggers = 0;
316   Float_t totalMult = 0;
317   for (UShort_t det=detFrom;det<=detTo;det++) {
318     Int_t nRings = (det == 1 ? 1 : 2);
319     for (UShort_t ir = 0; ir < nRings; ir++) {
320       Char_t   ring = (ir == 0 ? 'I' : 'O');
321       UShort_t nsec = (ir == 0 ? 20  : 40);
322       UShort_t nstr = (ir == 0 ? 512 : 256);
323       for (UShort_t sec =0; sec < nsec;  sec++) {
324         for (UShort_t strip = 0; strip < nstr; strip++) {
325           Float_t mult = fmdData->Multiplicity(det,ring,sec,strip);
326           if (mult == AliESDFMD::kInvalidMult) continue;
327           
328           if (fillHistograms)
329             fHistFMDSingle->Fill(mult);
330           
331           if (mult > fFMDLowCut)
332             totalMult = totalMult + mult;
333           else
334           {
335             if (totalMult > fFMDHitCut)
336               triggers++;
337               
338             if (fillHistograms)
339               fHistFMDSum->Fill(totalMult);
340               
341             totalMult = 0;
342           }
343         }
344       }
345     }
346   }
347   
348   return triggers;
349 }
350
351 Bool_t AliOfflineTrigger::FMDTrigger(const AliESDEvent* aEsd, AliceSide side) const
352 {
353   // Returns if the FMD triggered
354   //
355   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
356
357   Int_t triggers = FMDHitCombinations(aEsd, side, kFALSE);
358   
359   if (triggers > 0)
360     return kTRUE;
361     
362   return kFALSE;
363 }
364
365 Long64_t AliOfflineTrigger::Merge(TCollection* list)
366 {
367   // Merge a list of AliMultiplicityCorrection objects with this (needed for
368   // PROOF).
369   // Returns the number of merged objects (including this).
370
371   if (!list)
372     return 0;
373
374   if (list->IsEmpty())
375     return 1;
376
377   TIterator* iter = list->MakeIterator();
378   TObject* obj;
379
380   // collections of all histograms
381   const Int_t nHists = 8;
382   TList collections[nHists];
383
384   Int_t count = 0;
385   while ((obj = iter->Next())) {
386
387     AliOfflineTrigger* entry = dynamic_cast<AliOfflineTrigger*> (obj);
388     if (entry == 0) 
389       continue;
390
391     collections[0].Add(entry->fHistSPD);
392     collections[1].Add(entry->fHistV0A);
393     collections[2].Add(entry->fHistV0C);
394     collections[3].Add(entry->fHistZDC);
395     collections[4].Add(entry->fHistFMDA);
396     collections[5].Add(entry->fHistFMDC);
397     collections[6].Add(entry->fHistFMDSingle);
398     collections[7].Add(entry->fHistFMDSum);
399
400     count++;
401   }
402
403   fHistSPD->Merge(&collections[0]);
404   fHistV0A->Merge(&collections[1]);
405   fHistV0C->Merge(&collections[2]);
406   fHistZDC->Merge(&collections[3]);
407   fHistFMDA->Merge(&collections[4]);
408   fHistFMDC->Merge(&collections[5]);
409   fHistFMDSingle->Merge(&collections[6]);
410   fHistFMDSum->Merge(&collections[7]);
411
412   delete iter;
413
414   return count+1;
415 }
416
417 void AliOfflineTrigger::WriteHistograms() const
418 {
419   // write histograms to current directory
420   
421   if (!fHistSPD)
422     return;
423     
424   fHistSPD->Write();
425   fHistV0A->Write();
426   fHistV0C->Write();
427   fHistZDC->Write();
428   fHistFMDA->Write();
429   fHistFMDC->Write();
430   fHistFMDSingle->Write();
431   fHistFMDSum->Write();
432 }