]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/scripts/RunQATest.C
Non-implemented private copy constructor and assignment operator
[u/mrichter/AliRoot.git] / FMD / scripts / RunQATest.C
1 #ifdef _CINT__
2 class AliFMDInput;
3 class AliQADataMaker;
4 #endif
5 #ifdef BUILD 
6 # ifndef _CINT__
7 #  include <AliCDBManager.h>
8 #  include <AliQADataMaker.h>
9 #  include "AliFMDQADataMakerRec.h"
10 #  include "AliFMDQADataMakerSim.h"
11 #  include <AliQAv1.h>
12 #  include <TObjArray.h>
13 #  include <AliRecoParam.h>
14 #  include <AliRawReader.h>
15 #  include <AliQACheckerBase.h>
16 #  include <AliFMDQAChecker.h>
17 #  include <AliQAChecker.h>
18 #  include <TCanvas.h>
19 #  include <TMath.h>
20 #  include <AliFMDInput.h>
21 # endif
22
23 /** 
24  * Class to test the QA code 
25  * 
26  * 
27  */
28 class QATest : public AliFMDInput
29 {
30 public:
31   /** 
32    * Constructor 
33    */
34   QATest()
35     : AliFMDInput("galice.root"),
36       fMaker(0), 
37       fArray(0),
38       fSpecie(AliRecoParam::kLowMult), 
39       fCycleLength(-1)
40   {
41     for (Int_t i = 0; i < AliQAv1::kNTASKINDEX; i++)  
42       fTasks[i] = AliQAv1::kNULLTASKINDEX;
43     Int_t nArray = AliQAv1::kNTASKINDEX * AliRecoParam::kNSpecies;
44     Info("QAtest", "Allocating %dx%d=%d TObjArrays",
45          AliQAv1::kNTASKINDEX, AliRecoParam::kNSpecies, nArray);
46     fArray = new TObjArray*[nArray];
47     for (Int_t i = 0; i < nArray; i++) fArray[i] = 0;
48   }
49   /** 
50    * Set the location of the QA reference storage 
51    * 
52    * @param refloc Location of QA storage (e.g., local://${ALICE_ROOT}/QAref)
53    */
54   void SetQARefStorage(const char* refloc)
55   {
56     AliQAv1::SetQARefStorage(refloc);
57   }
58   /** 
59    * Calculate index of TObjArray* in fArray
60    * 
61    * @param specie Event species
62    * @param task   Task index 
63    * 
64    * @return Index. 
65    */
66   Int_t CalcArrayIndex(AliQAv1::TASKINDEX_t task,
67                        AliRecoParam::EventSpecie_t specie) const
68   {
69     Int_t es = AliRecoParam::AConvert(specie);
70     if (es >= AliRecoParam::kNSpecies) return -1;
71     
72     Int_t base = CalcSpeciesIndex(task);
73     if (base < 0) return base;
74     return base + es;
75   }
76   /** 
77    * Get index to species array 
78    * 
79    * @param task Task 
80    * 
81    * @return 
82    */
83   Int_t CalcSpeciesIndex(AliQAv1::TASKINDEX_t task) const
84   {
85     if (task >= AliQAv1::kNTASKINDEX) return -1;
86     // Species are consequtive 
87     return task * AliRecoParam::kNSpecies;
88   }
89     
90   /** 
91    * Initialise the code
92    * 
93    * 
94    * @return true on success, false otherwise 
95    */
96   Bool_t Init()
97   {
98     // --- Create the maker ------------------------------------------
99     if (IsLoaded(kHits)       || 
100         IsLoaded(kDigits)     || 
101         IsLoaded(kKinematics) ||
102         IsLoaded(kSDigits)) {
103       AddLoad(kHeader);
104       fMaker = new AliFMDQADataMakerSim();
105     }
106     else 
107       fMaker = new AliFMDQADataMakerRec();
108
109     // --- Figure out tasks ------------------------------------------
110     Int_t j = 0;
111     if (IsLoaded(kHits))      fTasks[j++]         = AliQAv1::kHITS;
112     if (IsLoaded(kDigits))    fTasks[j++]         = AliQAv1::kDIGITS;
113     if (IsLoaded(kSDigits))   fTasks[j++]         = AliQAv1::kSDIGITS;
114     if (IsLoaded(kRecPoints)) fTasks[j++]         = AliQAv1::kRECPOINTS;
115     if (IsLoaded(kESD))       fTasks[j++]         = AliQAv1::kESDS;
116     if (IsLoaded(kRaw))       fTasks[j++]         = AliQAv1::kRAWS;
117     if (j == 0) { 
118       AliError(Form("Loaded trees (%s) cannot be processed by QA", 
119                     LoadedString(true)));
120       return false;
121     }
122     // Int_t nTasks = j;
123
124     // --- Data maker ------------------------------------------------
125     Info("TestQA", "Creating data maker");
126     fMaker = new AliFMDQADataMakerRec();
127     
128     // --- Init all species histograms -------------------------------
129     Info("TestQA", "Setup data species");
130     AliQAv1* qa = AliQAv1::Instance();
131     for (unsigned int es = 0; es < AliRecoParam::kNSpecies; es++) {
132       AliRecoParam::EventSpecie_t specie = AliRecoParam::ConvertIndex(es);
133       fMaker->SetEventSpecie(specie);
134       qa->SetEventSpecie(specie);
135       for (Int_t i = 0; i < AliQAv1::kNTASKINDEX; i++) {
136         if (fTasks[i] == AliQAv1::kNULLTASKINDEX) continue;
137         Int_t k = CalcArrayIndex(fTasks[i], specie);
138         Info("Init", "Array for task %d (%s), specie %d (%s) @ %d/%d", 
139              fTasks[i], AliQAv1::GetTaskName(fTasks[i]).Data(), 
140              specie, AliRecoParam::GetEventSpecieName(specie), k,
141              AliQAv1::kNTASKINDEX * AliRecoParam::kNSpecies);
142         fArray[k] = fMaker->Init(fTasks[i], specie);
143         fArray[k]->ls();
144       }
145     }
146
147     // --- Start of cycle --------------------------------------------
148     Int_t  run  = AliCDBManager::Instance()->GetRun();
149     Bool_t same = false;
150     for (Int_t i = 0; i < AliQAv1::kNTASKINDEX; i++) {
151       if (fTasks[i] == AliQAv1::kNULLTASKINDEX) continue;
152       fMaker->StartOfCycle(fTasks[i], run, same);
153       same = !same;
154     }
155     
156     
157     return AliFMDInput::Init();
158   }
159   /** 
160    * Process the hits 
161    * 
162    * @return true on success 
163    */
164   virtual Bool_t ProcessHits()
165   {
166     fMaker->SetEventSpecie(fSpecie);
167     fMaker->MakeHits(fTreeH);
168     return true;
169   }
170   /** 
171    * Process the digits
172    * 
173    * @return true on success 
174    */
175   virtual Bool_t ProcessDigits()
176   {
177     fMaker->SetEventSpecie(fSpecie);
178     fMaker->MakeDigits(fTreeD);
179     return true;
180   }
181   /** 
182    * Process the summable digits
183    * 
184    * @return true on success 
185    */
186   virtual Bool_t ProcessSDigits()
187   {
188     fMaker->SetEventSpecie(fSpecie);
189     fMaker->MakeSDigits(fTreeS);
190     return true;
191   }
192   /** 
193    * Process the reconstructed points 
194    * 
195    * @return true on success 
196    */
197   virtual Bool_t ProcessRecPoints()
198   {
199     fMaker->SetEventSpecie(fSpecie);
200     fMaker->MakeRecPoints(fTreeR);
201     return true;
202   }
203   /** 
204    * Process the event summary data
205    * 
206    * @return true on success 
207    */
208   virtual Bool_t ProcesssESDs()
209   {
210     fMaker->SetEventSpecie(fSpecie);
211     fMaker->MakeESDs(fESDEvent);
212     return true;
213   }
214   /** 
215    * Process the raw data 
216    * 
217    * @return true on success 
218    */
219   virtual Bool_t ProcessRawDigits()
220   {
221     fMaker->SetEventSpecie(fSpecie);
222     fMaker->MakeRaws(fReader);
223     return true;
224   }
225   virtual Bool_t End()
226   {
227     Bool_t ret = AliFMDInput::End();
228     if (fCycleLength < 0) return ret;
229
230     if (!(fEventCount != 0 && (fEventCount % fCycleLength) == 0))
231       return ret;
232
233     // --- End of cycle - this calls the ecker ---------------------
234     Info("End", "End of cycle");
235     for (Int_t i = 0; i < AliQAv1::kNTASKINDEX; i++) {
236       if (fTasks[i] == AliQAv1::kNULLTASKINDEX) continue;
237       fMaker->EndOfCycle(fTasks[i]);
238     }
239
240     // --- Get the checker -------------------------------------------
241     Info("End", "Running checker");
242     AliQACheckerBase * checker = AliQAChecker::Instance()->
243       GetDetQAChecker(AliQAv1::GetDetIndex("FMD"));
244     ((AliFMDQAChecker*)checker)->SetDoScale();
245     
246     // --- Test: Remake plots ----------------------------------------
247     for (unsigned int idx = 0; idx < AliQAv1::kNTASKINDEX; idx++) {
248       // AliRecoParam::EventSpecie_t specie = AliRecoParam::ConvertIndex(es);
249       // if (!qa->IsEventSpecieSet(specie)) continue;
250       if (fTasks[idx] == AliQAv1::kNTASKINDEX) continue;
251       AliQAv1::TASKINDEX_t task = AliQAv1::TASKINDEX_t(idx); // AliQAv1::kRAWS;
252       AliQAv1::MODE_t      mode = AliQAv1::kRECMODE;
253       Int_t k = CalcSpeciesIndex(task);
254       Info("End", "Array for task %d (%s) @ %d: %p", 
255            task, AliQAv1::GetTaskName(task).Data(), k, fArray[k]);
256       if (!fArray[k]) continue;
257       fArray[k]->ls();
258       checker->MakeImage(&(fArray[k]), task, mode);
259     }
260     return ret;
261   }
262   /** 
263    * Called at the end of the job.  Runs the checkers
264    * 
265    * @return true on success 
266    */
267   virtual Bool_t Finish()
268   {
269     // --- Finish maker ----------------------------------------------
270     fMaker->Finish();
271   
272     // --- Get the checker ------------------------------------------- 
273     AliQACheckerBase * checker = AliQAChecker::Instance()->
274       GetDetQAChecker(AliQAv1::GetDetIndex("FMD")); 
275
276     // --- Get images from checker -----------------------------------
277     AliQAv1* qa = AliQAv1::Instance();
278     TObjArray* canvases = new TObjArray();
279     for (unsigned int es = 0; es < AliRecoParam::kNSpecies; es++) {
280       AliRecoParam::EventSpecie_t specie = AliRecoParam::ConvertIndex(es);
281       if (!qa->IsEventSpecieSet(specie)) continue;
282     
283       TCanvas* c = checker->GetImage(specie);
284       if (!c) continue; 
285
286       canvases->Add(c);
287     }
288   
289     // --- Create summary image --------------------------------------
290     TCanvas* out = new TCanvas("summary", "Summary", 1024, 1024);
291     out->SetFillColor(kWhite);
292     out->SetBorderSize(0);
293     out->SetBorderMode(0);
294   
295     Int_t nImgs = canvases->GetEntriesFast();
296     Int_t nX    = Int_t(TMath::Sqrt(nImgs) + .5);
297     Int_t nY    = nX;                           
298     out->Divide(nX, nY);
299     for (Int_t i = 0; i < nImgs; i++) { 
300       TVirtualPad* p = out->cd(i + 1);
301       p->SetRightMargin(0.001);
302       p->SetTopMargin(0.001);
303       if (!p) { 
304         Warning("TestQA", "No pad at index %d / %d", i+1, nImgs);
305         continue;
306       }
307
308       TCanvas* c = static_cast<TCanvas*>(canvases->At(i));
309       c->DrawClonePad();
310     }
311     out->Print("summary.png");    
312
313     return true;
314   }
315   void SetSpecie(AliRecoParam::EventSpecie_t s) { fSpecie = s; }
316 protected:
317   AliQADataMaker*             fMaker; // Data maker 
318   AliQAv1::TASKINDEX_t        fTasks[AliQAv1::kNTASKINDEX]; // Tasks to do
319   TObjArray**                 fArray;
320   AliRecoParam::EventSpecie_t fSpecie;
321   Int_t                       fCycleLength;
322   ClassDef(QATest, 0); 
323 };
324
325 #else
326 void
327 RunQATest(const char* src, const char* specie="low", Int_t runno=0)
328 {
329   gROOT->LoadMacro("$ALICE_ROOT/FMD/scripts/Compile.C");
330   gSystem->AddIncludePath("-DBUILD=1");
331   Compile("$ALICE_ROOT/../trunk/FMD/scripts/RunQATest.C", "+g");
332
333   AliCDBManager* cdb = AliCDBManager::Instance();
334   cdb->SetRun(runno);
335   
336   QATest*                     qaTest = new QATest;
337   TString                     what(src);
338   TString                     spec(specie);
339   Int_t                       colon = what.Index(":");
340   TString                     type = "";
341   if (colon != TString::kNPOS) { 
342     type = what(0, colon);
343     what = what(colon+1, what.Length()-colon-1);
344   }
345   else { 
346     if      (what.Contains("AliESD")) type = "esd";
347     else if (what.Contains("galice")) type = "sim";
348     else if (what.Contains(".raw"))   type = "raw";
349     else if (what.Contains(".root"))  type = "raw";
350   }
351   Info("RunQATest", "type=%s, what=%s specie=%s", 
352        type.Data(), what.Data(), spec.Data());
353
354   spec.ToLower();
355   if      (spec.Contains("low"))    qaTest->SetSpecie(AliRecoParam::kLowMult);
356   else if (spec.Contains("high"))   qaTest->SetSpecie(AliRecoParam::kHighMult);
357   else if (spec.Contains("cosmic")) qaTest->SetSpecie(AliRecoParam::kCosmic);
358   else if (spec.Contains("calib"))  qaTest->SetSpecie(AliRecoParam::kCalib);
359
360   type.ToLower();
361   if (type.CompareTo("esd") == 0) { 
362     qaTest->AddLoad(AliFMDInput::kESD);
363     qaTest->SetInputDir(what);
364   }
365   else if (type.CompareTo("raw") == 0) {
366     qaTest->AddLoad(AliFMDInput::kRaw);
367     qaTest->SetRawFile(what);
368   }
369   else if (type.CompareTo("sim") == 0) { 
370     qaTest->AddLoad(AliFMDInput::kHits);
371     qaTest->AddLoad(AliFMDInput::kSDigits);
372     qaTest->AddLoad(AliFMDInput::kDigits);
373   }
374   else if (type.CompareTo("rec") == 0) { 
375     qaTest->AddLoad(AliFMDInput::kDigits);
376     qaTest->AddLoad(AliFMDInput::kRecPoints);
377     qaTest->AddLoad(AliFMDInput::kESD);
378   }
379   else { 
380     Error("RunQATest", "Unknown type='%s' in '%s'", type.Data(), src);
381     return;
382   }
383   qaTest->Run(1000);
384 }
385 #endif
386   
387   
388
389 //
390 // EOF
391 //