]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSQASSDChecker.cxx
Changes from S. Rossegger :
[u/mrichter/AliRoot.git] / ITS / AliITSQASSDChecker.cxx
1 /**************************************************************************
2  * Copyright(c) 2007-2009, 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 // *****************************************
19 //  Checks the quality assurance 
20 //  by comparing with reference data
21 //  P. Cerello Apr 2008
22 //  INFN Torino
23
24 // --- ROOT system ---
25 #include "TH1.h"
26 #include "TString.h"
27 //#include "Riostream.h"
28
29 // --- AliRoot header files ---
30 #include "AliITSQASSDChecker.h"
31 #include "AliITSQADataMakerRec.h"
32 #include "AliLog.h"
33
34 ClassImp(AliITSQASSDChecker)
35 //__________________________________________________________________
36 AliITSQASSDChecker& AliITSQASSDChecker::operator = (const AliITSQASSDChecker& qac ) 
37 {
38   // Equal operator.
39   this->~AliITSQASSDChecker();
40   new(this) AliITSQASSDChecker(qac);
41   return *this;
42 }
43
44 void AliITSQASSDChecker::CheckRaws(TH1* histo) {  
45   // checker for RAWS
46   Double_t minSSDDataSize = 0;
47   Double_t maxSSDDataSize = 200;
48   Double_t minDDLDataSize = 0;
49   Double_t maxDDLDataSize = 50;
50   Double_t minLDCDataSize = 0;
51   Double_t maxLDCDataSize = 100;
52   Double_t minMeanDDLDataSize = 0;
53   Double_t maxMeanDDLDataSize = 50;
54   Double_t minMeanLDCDataSize = 0;
55   Double_t maxMeanLDCDataSize = 100;
56   //  Double_t maxOccupancy = 5;
57
58   TString histname = histo->GetName();
59
60   if (histname.EndsWith("SSDEventType")) {
61     if (histo->GetEntries()==0) {
62       AliWarning("Event type histogram is empty");
63     }
64     else if (histo->GetBinContent(histo->FindBin(7))==0) AliWarning("No type 7 (physics) events in EventType");
65   }
66
67   if (histname.EndsWith("SSDDataSize")) {
68     if (histo->GetEntries()==0) AliWarning("SSD data size histogram is empty");
69     if (histo->GetMean()>maxSSDDataSize||histo->GetMean()<minSSDDataSize) AliWarning(Form("SSD mean data size is %-.2g kB", histo->GetMean()));
70   }
71
72   if (histname.EndsWith("SSDDataSizePerDDL")) {
73     if (histo->GetEntries()==0) {
74       AliWarning("Data size per DDL histogram is empty");
75     }
76     else {
77       for(Int_t i = 512; i < 528; i++) {
78         if(histo->GetBinContent(histo->FindBin(i))==0) {
79            AliWarning(Form("Data size / DDL histogram: bin for DDL %i is empty",i));
80         }
81         else if(histo->GetBinContent(histo->FindBin(i))<minDDLDataSize||histo->GetBinContent(histo->FindBin(i))>maxDDLDataSize) AliWarning(Form("Data size DDL %i is %-.2g kB",i,histo->GetBinContent(histo->FindBin(i))));
82      }
83     }
84   }
85
86   if (histname.EndsWith("SSDDataSizePerLDC")) {
87     if (histo->GetEntries()==0) {
88       AliWarning("Data size per LDC histogram is empty");
89     }    
90     else {
91       AliInfo(Form("Data size per LDC histogram has %f entries",histo->GetEntries()));
92       for(Int_t i = 170; i < 178; i++) {
93         if(histo->GetBinContent(histo->FindBin(i))==0) {
94           AliWarning(Form("Data size / LDC histogram: bin for LDC %i is empty",i));
95         }
96         else if(AliITSQADataMakerRec::AreEqual(histo->GetBinContent(histo->FindBin(i)),minLDCDataSize) ||histo->GetBinContent(histo->FindBin(i))>maxLDCDataSize) AliWarning(Form("Data size LDC %i is %-.2g kB",i,histo->GetBinContent(i)));
97       }
98     }
99   }
100
101   if (histname.EndsWith("SSDLDCId")) {
102     if (histo->GetEntries()==0) {
103       AliWarning("LDC ID histogram is empty");
104     }    
105     else {
106       for(Int_t i = 170; i < 177; i++) {
107         if(histo->GetBinContent(histo->FindBin(i))==0) {
108           AliWarning(Form("LDC ID histogram: No entries for LDC %i",i));
109         }
110         else if(histo->GetBinContent(histo->FindBin(i))!=histo->GetBinContent(histo->FindBin(i+1))) {
111           AliWarning("LDC Id distribution is not uniform");
112           i=176;
113         }
114       }
115     }
116   }
117
118   if (histname.EndsWith("SSDDDLId")) {
119     if (histo->GetEntries()==0) {
120       AliWarning("DDL ID histogram is empty");
121     }
122     else {
123       for(Int_t i = 512; i < 527; i++) {
124         if(histo->GetBinContent(histo->FindBin(i))==0) {
125           AliWarning(Form("DDL ID histogram: No entries for DDL %i",i));
126         }
127         else if(histo->GetBinContent(histo->FindBin(i))!=histo->GetBinContent(histo->FindBin(i+1))) {
128           AliWarning("DDL Id distribution is not uniform");
129           i=526;
130         }
131       }
132     }
133   }
134
135   if (histname.Contains("SSDDataSizeLDC")) {
136     if (histo->GetEntries()==0) {
137       AliWarning(Form("LDC %s data size distribution is empty", histname(histname.Length()-3,3).Data()));
138     }
139     else if (histo->GetMean()<minMeanLDCDataSize||histo->GetMean()>maxMeanLDCDataSize) AliWarning(Form("Mean data size of LDC %s is %-.2g kB",histname(histname.Length()-3,3).Data(), histo->GetMean()));
140   }
141
142   if (histname.Contains("SSDDataSizeDDL")) {
143     if (histo->GetEntries()==0) {
144       AliWarning(Form("DDL %s data size distribution is empty", histname(histname.Length()-3,3).Data()));
145     } 
146     else if (histo->GetMean()<minMeanDDLDataSize||histo->GetMean()>maxMeanDDLDataSize) AliWarning(Form("Mean data size of DDL %s is %-.2g kB",histname(histname.Length()-3,3).Data(), histo->GetMean()));
147   }
148
149   if (histname.Contains("SSDAverageOccupancy")) {
150  
151     const char* side = "";
152     int ladder = 0;
153     int layernr = 0;
154
155     if (histname.EndsWith("5")) layernr = 499;
156     if (histname.EndsWith("6")) layernr = 599;
157
158     for (Int_t i = 1; i < histo->GetNbinsY() + 1; i++) { //ladder/side loop
159       if(i==3.*int(i/3.)){
160         ladder=int(i/3.)+layernr;
161         side="P side";
162       }
163       else if(i==3.*int(i+1/3.)){
164         ladder=int((i+1)/3.)+layernr;
165         side="N side";
166       }
167
168       for (Int_t j = 1; j < histo->GetNbinsX() + 1; j++) { //module loop
169         //if(histo->GetBinContent(j,i)>maxOccupancy)
170           // AliWarning(Form("Occupancy ladder %i, module %i, %s is %-.2f %%",ladder,j,side, histo->GetBinContent(j,i)));
171       }//module loop
172     }//ladder loop
173   }
174
175 }
176
177
178
179 //__________________________________________________________________
180 Double_t AliITSQASSDChecker::Check(AliQAv1::ALITASK_t /*index*/, const TObjArray * list, const AliDetectorRecoParam * /*recoParam*/) { 
181   // main checker method 
182   AliDebug(AliQAv1::GetQADebugLevel(),Form("AliITSQASSDChecker called with offset: %d\n", fSubDetOffset));
183
184   AliInfo(Form("AliITSQASSDChecker called with offset: %d\n", fSubDetOffset) );
185   //cout<<"(AliITSQASSDChecker::Check): List name "<<list->GetName()<<endl;
186   Double_t test = 0.0  ;
187   Int_t count = 0 ;
188   TString listname = list->GetName();
189
190   if (list->GetEntries() == 0){
191     test = 1. ; // nothing to check
192   }
193   else {
194
195     TIter next(list) ;
196     TH1 * hdata ;
197     count = 0 ;
198     while ( (hdata = dynamic_cast<TH1 *>(next())) ) {
199       if (hdata) {
200         TString histname = hdata->GetName();
201         if(!histname.Contains("fHistSSD")) continue;
202         Double_t rv = 0.;
203         if(hdata->GetEntries()>0) {
204            rv = 1;
205
206            //if(histname.Contains("PerDDL")) cout << "(AliITSQASSDChecker::Check) " << histname << " has " << hdata->GetEntries() << " entries. Mean: " << hdata->GetMean() << endl;
207        
208        //    if(hdata->GetMean()>0&&!histname.Contains("_Ladder")) cout << "(AliITSQASSDChecker::Check) " << histname << " not empty! " << hdata->GetEntries() << " entries. Mean: " << hdata->GetMean() << endl;
209         }
210
211     //    if (listname.Contains("Raws")) CheckRaws(hdata);
212    //     if (listname.Contains("RecPoints")) CheckRecPoints(hdata);
213
214         //AliDebug(AliQAv1::GetQADebugLevel(), Form("%s -> %f", hdata->GetName(), rv)) ;
215         //cout<<hdata->GetName()<<" - "<<rv<<endl;
216         count++ ;
217         test += rv ;
218       }
219       else{
220         AliError("Data type cannot be processed") ;
221       }
222     }
223     if (count != 0) {
224       if (AliITSQADataMakerRec::AreEqual(test,0.)) {
225         AliWarning("Histograms are there, but they are all empty: setting flag to kWARNING");
226         test = 0.5;  //upper limit value to set kWARNING flag for a task
227       }
228       else {
229         test /= count ;
230       }
231     }
232   }
233   
234   //AliDebug(AliQAv1::GetQADebugLevel(), Form("Test Result = %f", test)) ;
235   //cout<<"Test result: "<<test<<endl;
236
237   return test ;
238
239   //return 0.;
240
241
242 }
243
244 //__________________________________________________________________
245 void AliITSQASSDChecker::SetTaskOffset(Int_t TaskOffset){
246   // defines offset for SSD
247   fSubDetOffset = TaskOffset;
248 }
249
250 //__________________________________________________________________
251 void AliITSQASSDChecker::SetStepBit(const Double_t *steprange) {
252   // defines step range
253   fStepBitSSD = new Double_t[AliQAv1::kNBIT];
254   for(Int_t bit=0;bit<AliQAv1::kNBIT;bit++)
255     {
256       fStepBitSSD[bit]=steprange[bit];
257     }
258 }
259
260 //__________________________________________________________________
261 void  AliITSQASSDChecker::SetSSDLimits(const Float_t *lowvalue, const Float_t * highvalue){
262   // defines 
263   fLowSSDValue = new Float_t[AliQAv1::kNBIT];
264   fHighSSDValue= new Float_t[AliQAv1::kNBIT];
265
266   for(Int_t bit=0;bit<AliQAv1::kNBIT;bit++)
267     {
268       fLowSSDValue[bit]=lowvalue[bit];
269       fHighSSDValue[bit]= highvalue[bit];
270     }
271
272 }