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