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