]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONQADataMakerRec.cxx
- Fix bug preventing reset of MTR histograms at SOR in DQM plots (Diego)
[u/mrichter/AliRoot.git] / MUON / AliMUONQADataMakerRec.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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
16 // $Id$
17
18 // --- MUON header files ---
19 #include "AliMUONQADataMakerRec.h"
20
21 //-----------------------------------------------------------------------------
22 /// \class AliMUONQADataMakerRec
23 ///
24 /// MUON base class for quality assurance data (histo) maker
25 ///
26 /// It is now only a steering class for the two subclasses AliMUONTrackerQADataMakerRec
27 /// and AliMUONTriggerQADataMakerRec
28 ///
29 /// \author C. Finck, D. Stocco, L. Aphecetche, A. Blanc
30
31 #include "AliDAQ.h"
32 #include "AliMUONTrackerQADataMakerRec.h"
33 #include "AliMUONTriggerQADataMakerRec.h"
34 #include "AliQAChecker.h"
35 #include "AliRawReader.h"
36 #include "AliRawEventHeaderBase.h"
37
38 /// \cond CLASSIMP
39 ClassImp(AliMUONQADataMakerRec)
40 /// \endcond
41            
42 //____________________________________________________________________________ 
43 AliMUONQADataMakerRec::AliMUONQADataMakerRec(Bool_t tracker, Bool_t trigger) : 
44 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kMUON), "MUON Quality Assurance Data Maker"),
45 fTracker(tracker ? new AliMUONTrackerQADataMakerRec(this) : 0x0),
46 fTrigger(trigger ? new AliMUONTriggerQADataMakerRec(this) : 0x0)
47 {
48   /// ctor
49 }
50
51 //__________________________________________________________________
52 AliMUONQADataMakerRec::~AliMUONQADataMakerRec()
53 {
54     /// dtor
55   delete fTracker;
56   delete fTrigger;
57 }
58
59 //____________________________________________________________________________ 
60 Int_t AliMUONQADataMakerRec::Add2List(TH1 * hist, const Int_t index, AliQAv1::TASKINDEX_t task, const Bool_t expert, const Bool_t image, const Bool_t saveForCorr)
61 {
62   TObjArray** list = GetList(task);
63   if (list)
64   {
65     return Add2List(hist,index,list,expert,image,saveForCorr);
66   }
67   return -1;
68 }
69
70 //____________________________________________________________________________ 
71 void AliMUONQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray** list)
72 {
73   /// Detector specific actions at end of cycle
74   //
75   ResetEventTrigClasses(); // RS
76   //
77   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) 
78   {
79     if (! IsValidEventSpecie(specie, list)  ) continue;
80     
81     SetEventSpecie(AliRecoParam::ConvertIndex(specie)); // needed by the GetXXXData methods
82     
83     if ( task == AliQAv1::kRAWS ) 
84     {
85       if ( fTracker ) fTracker->EndOfDetectorCycleRaws(specie,list);
86       if ( fTrigger ) fTrigger->EndOfDetectorCycleRaws(specie,list);
87     }  
88     else if ( task == AliQAv1::kRECPOINTS )
89     {
90       // normalize recpoints histograms
91       if ( fTracker ) fTracker->EndOfDetectorCycleRecPoints(specie,list);
92       if ( fTrigger ) fTrigger->EndOfDetectorCycleRecPoints(specie,list);
93     }
94     else if ( task == AliQAv1::kESDS ) 
95     {
96       // normalize esds histograms
97       if ( fTracker ) fTracker->EndOfDetectorCycleESDs(specie,list);
98       if ( fTrigger ) fTrigger->EndOfDetectorCycleESDs(specie,list);
99     }
100     else if ( task == AliQAv1::kDIGITSR ) 
101     {
102       if ( fTracker ) fTracker->EndOfDetectorCycleDigits(specie,list);        
103       if ( fTrigger ) fTrigger->EndOfDetectorCycleDigits(specie,list);
104     }
105     else
106     {
107       AliFatal(Form("Not implemented for task %s",AliQAv1::GetTaskName(task).Data()));
108     }
109   } // loop on specie
110     
111   // do the QA checking
112   AliQAChecker::Instance()->Run(AliQAv1::kMUON,task,list,const_cast<AliDetectorRecoParam*>(GetRecoParam()));
113 }
114
115 //____________________________________________________________________________ 
116 TObject* AliMUONQADataMakerRec::GetData(AliQAv1::TASKINDEX_t task, const Int_t index)
117 {
118   TObjArray** list = GetList(task);
119   if (list) return GetData(list,index);
120   return 0x0;
121 }
122
123 //____________________________________________________________________________ 
124 TObjArray** AliMUONQADataMakerRec::GetList(AliQAv1::TASKINDEX_t task)
125 {
126   //  enum TASKINDEX_t {
127   //    kNULLTASKINDEX=-1, kRAWS, kHITS, kSDIGITS, kDIGITS, kDIGITSR, kRECPOINTS, kTRACKSEGMENTS, kRECPARTICLES, kESDS, kNTASKINDEX };
128   if ( task == AliQAv1::kRAWS ) 
129   {
130       return fRawsQAList;
131   }
132   else if ( task == AliQAv1::kDIGITS || task == AliQAv1::kDIGITSR )
133   {
134     return fDigitsQAList;
135   }
136   else if ( task == AliQAv1::kRECPOINTS ) 
137   {
138     return fRecPointsQAList;
139   }
140   else
141   {
142       AliFatal(Form("task %s not supported here yet",AliQAv1::GetTaskName(task).Data()));
143   }
144   return 0x0;
145 }
146
147 //____________________________________________________________________________ 
148 void AliMUONQADataMakerRec::InitRaws()
149 {
150   /// create Raws histograms in Raws subdir
151         
152   if ( fTracker ) fTracker->InitRaws();
153   if ( fTrigger ) fTrigger->InitRaws();
154   //
155   ClonePerTrigClass(AliQAv1::kRAWS); // this should be the last line
156 }
157
158 //__________________________________________________________________
159 void AliMUONQADataMakerRec::InitDigits() 
160 {
161   /// Initialized Digits spectra 
162   if ( fTracker ) fTracker->InitDigits();
163   if ( fTrigger ) fTrigger->InitDigits();
164   //
165   ClonePerTrigClass(AliQAv1::kDIGITS); // this should be the last line
166
167
168 //____________________________________________________________________________ 
169 void AliMUONQADataMakerRec::InitRecPoints()
170 {
171         /// create Reconstructed Points histograms in RecPoints subdir
172   if ( fTracker ) fTracker->InitRecPoints();
173   if ( fTrigger ) fTrigger->InitRecPoints();
174   //
175   ClonePerTrigClass(AliQAv1::kRECPOINTS); // this should be the last line
176 }
177
178
179 //____________________________________________________________________________ 
180 void AliMUONQADataMakerRec::InitESDs()
181 {
182   ///create ESDs histograms in ESDs subdir
183   if ( fTracker ) fTracker->InitESDs();
184   if ( fTrigger ) fTrigger->InitESDs();
185   //
186   ClonePerTrigClass(AliQAv1::kESDS); // this should be the last line
187 }
188
189 //____________________________________________________________________________
190 void AliMUONQADataMakerRec::MakeRaws(AliRawReader* rawReader)
191 {
192   /// make QA for rawdata
193   /// Note that we do not call the sub-datamaker MakeRaws method
194   /// for events where the MCH or MTR is not part of the readout...
195   
196   if ( !rawReader || !rawReader->GetDetectorPattern() ) return;
197   
198   UInt_t clmask = rawReader->GetDetectorPattern()[0];
199     
200   if ( fTracker && rawReader->GetType() == AliRawEventHeaderBase::kPhysicsEvent ) 
201   {
202     UInt_t mchMask = AliDAQ::DetectorPattern(" MUONTRK ");
203     if ( clmask & mchMask ) 
204     {
205       rawReader->Reset();
206       fTracker->MakeRaws(rawReader);
207     }
208   }
209   
210   if ( fTrigger && (rawReader->GetType() == AliRawEventHeaderBase::kPhysicsEvent ||
211                     rawReader->GetType() == AliRawEventHeaderBase::kCalibrationEvent ) )
212   {
213     UInt_t mtrMask = AliDAQ::DetectorPattern(" MUONTRG ");
214     if ( clmask & mtrMask )
215     {
216       rawReader->Reset();    
217       fTrigger->MakeRaws(rawReader);
218     }
219   }
220   //
221   IncEvCountCycleRaws();
222   IncEvCountTotalRaws();
223   //
224 }
225
226 //__________________________________________________________________
227 void AliMUONQADataMakerRec::MakeDigits()         
228 {
229   /// makes data from Digits
230   
231   AliFatal("Not implemented");
232 }
233
234 //__________________________________________________________________
235 void AliMUONQADataMakerRec::MakeDigits(TTree* digitsTree)         
236 {
237   /// makes data from Digits
238
239   // Do nothing in case of calibration event
240   if ( GetEventSpecie() == AliRecoParam::kCalib ) return;
241
242   if ( fTracker ) fTracker->MakeDigits(digitsTree);
243   if ( fTrigger ) fTrigger->MakeDigits(digitsTree);  
244   //
245   IncEvCountCycleDigits();
246   IncEvCountTotalDigits();
247   //
248 }
249
250 //____________________________________________________________________________
251 void AliMUONQADataMakerRec::MakeRecPoints(TTree* clustersTree)
252 {
253         /// Fill histograms from treeR
254
255   // Do nothing in case of calibration event
256   if ( GetEventSpecie() == AliRecoParam::kCalib ) return;
257         
258   if ( fTracker ) fTracker->MakeRecPoints(clustersTree);
259   if ( fTrigger ) fTrigger->MakeRecPoints(clustersTree);  
260   //
261   IncEvCountCycleRecPoints();
262   IncEvCountTotalRecPoints();
263   //
264 }
265
266 //____________________________________________________________________________
267 void AliMUONQADataMakerRec::MakeESDs(AliESDEvent* esd)
268 {
269   /// make QA data from ESDs
270
271   // Do nothing in case of calibration event
272   if ( GetEventSpecie() == AliRecoParam::kCalib ) return;
273   
274   if ( fTracker ) fTracker->MakeESDs(esd);
275   if ( fTrigger ) fTrigger->MakeESDs(esd);  
276   //
277   IncEvCountCycleESDs();
278   IncEvCountTotalESDs();
279   //
280  }
281
282 //____________________________________________________________________________ 
283 void AliMUONQADataMakerRec::ResetDetector(AliQAv1::TASKINDEX_t task)
284 {
285   /// Reset internals
286   
287   for (int spec = 0; spec < AliRecoParam::kNSpecies; ++spec) 
288   {
289     if (!AliQAv1::Instance()->IsEventSpecieSet(AliRecoParam::ConvertIndex(spec)))
290       continue;
291     
292     if ( task == AliQAv1::kRAWS ) 
293     {
294       if (fTracker) fTracker->ResetDetectorRaws(fRawsQAList[spec]);
295       if (fTrigger) fTrigger->ResetDetectorRaws(fRawsQAList[spec]);
296     }
297     else if ( task == AliQAv1::kRECPOINTS )
298     {
299       if (fTracker) fTracker->ResetDetectorRecPoints(fRecPointsQAList[spec]);
300       if (fTrigger) fTrigger->ResetDetectorRecPoints(fRecPointsQAList[spec]);
301     }
302     else if ( task == AliQAv1::kESDS ) 
303     {
304       if (fTracker) fTracker->ResetDetectorESDs(fESDsQAList[spec]);
305       if (fTrigger) fTrigger->ResetDetectorESDs(fESDsQAList[spec]);
306     }
307     else if ( task == AliQAv1::kDIGITS ) 
308     {
309       if (fTracker) fTracker->ResetDetectorDigits(fDigitsQAList[spec]);
310       if (fTrigger) fTrigger->ResetDetectorDigits(fDigitsQAList[spec]);
311     }
312     else
313     {
314       AliFatal(Form("Not implemented for task %s",AliQAv1::GetTaskName(task).Data()));
315     }
316   }
317 }
318
319 //____________________________________________________________________________ 
320 void AliMUONQADataMakerRec::StartOfDetectorCycle()
321 {
322   /// Detector specific actions at start of cycle  
323   
324 }