]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSOnlineSPDfoAnalyzer.cxx
coverity fix
[u/mrichter/AliRoot.git] / ITS / AliITSOnlineSPDfoAnalyzer.cxx
1 /**************************************************************************
2  * Copyright(c) 2008-2010, 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 // Author: A. Mastroserio                                 //
20 // This class is used in the detector algorithm framework //
21 // to process the data stored in special container files  //
22 // (see AliITSOnlineSPDfo). The "good" set of DAC values  //
23 // is extracted.                                          //
24 ////////////////////////////////////////////////////////////
25
26 #include <TFile.h>
27 #include <TMath.h>
28 #include <TString.h>
29 #include <TStyle.h>
30 #include <TF1.h>
31 #include <TH1D.h>
32 #include <TH2D.h>
33 #include <TArrayI.h>
34 #include <TCanvas.h>
35 #include <THnSparse.h>
36 #include <TKey.h>
37 #include <iostream>
38 #include <fstream>
39 #include "AliITSOnlineSPDfoChipConfig.h"
40 #include "AliITSOnlineSPDfoChip.h"
41 #include "AliITSOnlineSPDfoInfo.h"
42 #include "AliITSOnlineSPDfo.h"
43 #include "AliITSOnlineSPDfoAnalyzer.h"
44 #include "AliLog.h"
45
46 using std::ifstream;
47
48 AliITSOnlineSPDfoAnalyzer::AliITSOnlineSPDfoAnalyzer(const TString fileName, Bool_t readFromGridFile) :
49   fFileName(0),
50   fNdims(0),
51   fNbins(0x0),
52   fXmin(0x0),
53   fXmax(0x0),
54   fFOHandler(0x0),
55   fHighOccupancyCheck(kTRUE)
56 {
57   // constructor
58   fFileName = fileName;
59   for(Int_t iqual =0; iqual<kNqualityFlags; iqual++) fGeneralThresholds[iqual]=0;
60   
61   for (UInt_t chipNr=0; chipNr<10; chipNr++) {
62     for (UInt_t hs=0; hs<6; hs++) {
63       for(Int_t i =0; i< kNqualityFlags; i++) fNh[i][hs][chipNr]=NULL;
64     }
65   }
66   
67   Init(readFromGridFile); 
68 }
69 //-------------------------------------------------------------------
70 AliITSOnlineSPDfoAnalyzer::AliITSOnlineSPDfoAnalyzer(const AliITSOnlineSPDfoAnalyzer& foan) :
71   fFileName(foan.fFileName),
72   fNdims(foan.fNdims),
73   fNbins(foan.fNbins),
74   fXmin(foan.fXmin),
75   fXmax(foan.fXmax),
76   fFOHandler(foan.fFOHandler),
77   fHighOccupancyCheck(foan.fHighOccupancyCheck)
78 {
79   //
80   // copy constructor, only copies the filename and params (not the processed data
81   //
82   
83   for(Int_t iqual =0; iqual<3; iqual++) fGeneralThresholds[iqual] =foan.fGeneralThresholds[iqual];
84   for (UInt_t chipNr=0; chipNr<10; chipNr++) {
85     for (UInt_t hs=0; hs<6; hs++) {
86       for(Int_t i=0; i<kNqualityFlags;i++) fNh[i][hs][chipNr]=NULL;
87     }
88   }
89   
90   Init();
91 }
92 //-------------------------------------------------------------------
93 AliITSOnlineSPDfoAnalyzer::~AliITSOnlineSPDfoAnalyzer() 
94 {
95   //
96   // destructor
97   // 
98   
99   for (UInt_t hs=0; hs<6; hs++) {
100     for (UInt_t chipNr=0; chipNr<10; chipNr++) {
101       for(Int_t i=0; i<kNqualityFlags ; i++ ) if(fNh[i][hs][chipNr]!=NULL) delete fNh[i][hs][chipNr];          
102     }
103   }
104   delete fFOHandler;
105 }
106 //-------------------------------------------------------------------
107 AliITSOnlineSPDfoAnalyzer& AliITSOnlineSPDfoAnalyzer::operator=(const AliITSOnlineSPDfoAnalyzer& foan) 
108 {
109   // assignment operator, only copies the filename and params (not the processed data)
110   if (this!=&foan) {
111     for (UInt_t hs=0; hs<6; hs++) {
112       for (UInt_t chipNr=0; chipNr<10; chipNr++) {
113         for(Int_t i=0; i<kNqualityFlags ; i++ ) if(fNh[i][hs][chipNr]!=NULL) delete fNh[i][hs][chipNr];
114       }
115     }
116     
117     fFileName=foan.fFileName;
118     fHighOccupancyCheck=foan.fHighOccupancyCheck;
119     Init();    
120   }
121   return *this;
122 }
123 //-------------------------------------------------------------------
124 void AliITSOnlineSPDfoAnalyzer::Init(Bool_t readFromGridFile) 
125 {
126   //
127   // first checks type of container and then initializes container obj
128   //
129   if (!readFromGridFile) {
130     TFile* p0 = TFile::Open(fFileName.Data());
131     if (p0 == NULL) {
132       printf("no file open!"); 
133       return;
134     }
135     else { 
136       fFOHandler = new AliITSOnlineSPDfo();
137       fFOHandler->SetFile(fFileName);
138     }
139     p0->Close();   
140   }
141 }
142 //-------------------------------------------------------------------
143 void AliITSOnlineSPDfoAnalyzer::SetGeneralThresholds(Float_t thre[3])
144 {
145   // here the settings for the 3 quality flags are defined
146   for(Int_t i=0; i< 3; i++) fGeneralThresholds[i]=thre[i];  
147 }
148 //-------------------------------------------------------------------
149 void AliITSOnlineSPDfoAnalyzer::SetNdimensions()
150 {
151   //
152   // here the axis of the N-dim histograms are setted 
153   // 
154   if(!fFOHandler) {
155     printf("no fo object. Exiting...\n"); 
156     return;
157   }
158   
159   TArrayI array = fFOHandler->GetDACscanParams();
160   fNdims = array.GetSize()/3;
161   fNbins = new Int_t[fNdims];
162   fXmin = new Double_t[fNdims];
163   fXmax = new Double_t[fNdims];
164   for(Int_t i=0; i< fNdims; i++){
165     fXmin[i] = array.At(3*i)-0.25;
166     fXmax[i] = array.At(3*i+1)+0.25;
167     fNbins[i] = (Int_t)((fXmax[i]-fXmin[i])/0.5)+1;//to avoid Int->Double conversion problems when checking results. 
168   } 
169 }
170 //-------------------------------------------------------------------
171 void AliITSOnlineSPDfoAnalyzer::BuildTHnSparse(Int_t ihs, Int_t ichip)
172 {
173   //
174   // here the N-dim histogram is booked per chip in one HS 
175   //
176   
177   if(!fNdims) SetNdimensions();
178   
179   for(Int_t iqual =0; iqual < 3; iqual++) fNh[iqual][ihs][ichip] = new THnSparseI(Form("h%i_HS%i_C%i",iqual,ihs,ichip), Form("h%i_HS%i_C%i",iqual,ihs,ichip), fNdims, fNbins, fXmin, fXmax);
180 }
181 //-------------------------------------------------------------------
182 Int_t AliITSOnlineSPDfoAnalyzer::Select(const AliITSOnlineSPDfoChip *chip) const
183 {
184   //
185   // Selects the DAC values if in the chip: the I configuration corresponds to
186   // 0% efficiency ( no FO response case). All the others shoud be within the 
187   // predefined thresholds
188   
189   if(!fFOHandler->GetFOscanInfo()) {
190     printf("no general information object in the file. Exiting...\n");
191     return -1;
192   }
193   
194   Double_t npulses = (fFOHandler->GetFOscanInfo())->GetNumTriggers();
195   
196   if(npulses == 0. ) {
197     Info("AliITSOnlineSPDfoAnalyzer::Select","no trigger pulses set. Exiting...");
198     return -999;
199   } 
200   
201   TObjArray *array = chip->GetChipConfigInfo();
202   if(!array) {
203     printf("No measurement array found in the chip!!\n");
204     return 0;
205   }
206   
207   Int_t quality = -1;
208   
209   Float_t counts = 0;
210   
211   Int_t processedconfigurations = chip->GetNumberOfChipConfigs();
212   
213  
214   
215   
216   for(Int_t isteps =0; isteps < processedconfigurations; isteps++){ 
217     
218     Int_t matrixId = ((AliITSOnlineSPDfoChipConfig*)array->At(isteps))->GetChipConfigMatrixId();
219     counts = (Float_t)(((AliITSOnlineSPDfoChipConfig*)array->At(isteps))->GetChipConfigCounter());
220     if(matrixId==0 && counts > 0) return -1;
221     if(fHighOccupancyCheck &&  matrixId ==6) continue;
222     
223     Float_t efficiency = counts/npulses;
224     
225     if(matrixId > 0){
226       Int_t response = IsSelected(efficiency);
227       if( response >=0) {
228         if(quality < response) quality = response;
229       }
230       else return -1; 
231     }
232   }
233   return quality;
234 }
235 //-----------------------------------------------------
236 Int_t AliITSOnlineSPDfoAnalyzer::IsSelected(Float_t eff) const
237 {
238   //
239   // returns the quality of the selection 
240   //
241   
242   for(Int_t i=0; i<3; i++){  
243     if(eff <= 1.+ fGeneralThresholds[i] && eff >= 1. - fGeneralThresholds[i]  ) return i;
244   }
245   return -1;
246 }
247 //----------------------------------------------------
248 void AliITSOnlineSPDfoAnalyzer::Process()
249
250   //
251   // The procedure is the following:
252   // - DAC 4-tuples are checked  
253   // - if the 4-tuple survives the selection, the chip-related histograms are filled.
254   // (- Per each histogram the mean values of each axis are taken)
255   //
256   
257   if(!fFOHandler) { 
258     Warning("AliITSOnlineSPDfoAnalyzer::Process","no fo object. Exiting.. \n");
259     return;
260   } 
261   
262   TKey *key;
263   TIter iter((fFOHandler->GetFile())->GetListOfKeys());  
264   while ((key = (TKey*)(iter.Next()))) {
265     TString classname = key->GetClassName();
266     if(classname.Contains("OnlineSPD")) break;
267     
268     TObjArray *array = (TObjArray*)(fFOHandler->GetFile())->Get(key->GetName()); // array of chips corresponding to the DACS (size 1-60)
269     if(!array){
270       printf("no array found! Exiting...\n");
271       break; 
272     }
273     
274     Double_t *dacvalues = fFOHandler->GetDACvaluesD(key->GetName(), GetFOHandler()->GetFOscanInfo()->GetNumDACindex());
275     
276     for(Int_t i=0; i< array->GetSize(); i++){
277       AliITSOnlineSPDfoChip *chip = (AliITSOnlineSPDfoChip *)array->At(i); 
278           
279       if(!chip) continue;
280       Int_t hs = chip->GetActiveHS();
281       Int_t chipid = chip->GetChipId();     
282       Int_t quality = Select(chip);
283       if(quality<0) continue; 
284       if(!fNh[quality][hs][chipid]) BuildTHnSparse(hs,chipid);
285       fNh[quality][hs][chipid]->Fill(dacvalues);       
286     } 
287
288     if(dacvalues) delete [] dacvalues;
289   } 
290 }
291 //---------------------------------------------
292 void AliITSOnlineSPDfoAnalyzer::WriteToFile(TString outputfile)
293 {
294   //
295   // The 4-dim histograms are stored into a file
296   //
297   
298   TFile * f = TFile::Open(outputfile.Data(),"recreate");
299   for(Int_t ihs =0; ihs < 6; ihs++) {
300     for(Int_t ichip =0; ichip < 10; ichip++){
301       for(Int_t i=0; i<kNqualityFlags ; i++ ) {
302       if(fNh[i][ihs][ichip]) f->WriteObjectAny(fNh[i][ihs][ichip],"THnSparse",Form("h%i_hs%i_chip%i",i,ihs,ichip));
303       }
304       
305     }
306   }
307   f->Close();
308 }
309 //---------------------------------------------
310 void AliITSOnlineSPDfoAnalyzer::CheckResults(TString filename, Int_t hs, Int_t ichip, Int_t iqual) const
311 {
312   //  
313   //The chip related 4-dim histograms are produced and stored into eps files
314   //    
315   
316   TFile *f = TFile::Open(filename.Data());
317   if(!f) {
318     Info("AliITSOnlineSPDfoAnalyzer::CehckResults","no file open!. Exiting...\n");
319     return; 
320   }
321   
322   THnSparse *hn;
323   TObject *obj;
324   
325   TIter iter(f->GetListOfKeys());
326   while((obj=iter.Next())){
327     TString name = obj->GetName();
328     hn=(THnSparse*)f->Get(name.Data()); 
329     if(name.Contains(Form("hs%i",hs)) && name.Contains(Form("chip%i",ichip)) && name.Contains(Form("h%i",iqual))) GetCanvases(hn,hs,ichip,iqual);
330   } 
331 }
332 //-----------------------------------------------------------------------------------------------
333 void AliITSOnlineSPDfoAnalyzer::GetCanvases(const THnSparse *hn,Int_t hs, Int_t chip, Int_t iqual) const
334 {
335   //
336   // 1-Dim and 2 Dim Projections of 4-dim histograms are visualized in canvases per quality selection
337   //
338   //
339   
340   if(!hn) {printf("no thnsparse...exiting!\n"); return;} 
341   
342   gStyle->SetPalette(1);
343   
344   TString dacname[4];
345   
346   if(  (fFOHandler->GetFOscanInfo())->GetDACindex(0) == 20 ) dacname[0] = "FOPOL";
347   else dacname[0] = "possible DAC-name/ DAC-register-number mismatch";  
348   if(  (fFOHandler->GetFOscanInfo())->GetDACindex(1) == 17 ) dacname[1] = "CONVPOL";
349   else dacname[1] = "possible DAC-name/ DAC-register-number mismatch";
350   if(  (fFOHandler->GetFOscanInfo())->GetDACindex(2) == 16 ) dacname[2] = "COMPREF";
351   else dacname[2] = "possible DAC-name/ DAC-register-number mismatch";
352   if(  (fFOHandler->GetFOscanInfo())->GetDACindex(3) == 39 ) dacname[3] = "pre_VTH";
353   else dacname[3] = "possible DAC-name/ DAC-register-number mismatch";
354   
355   TString titles = Form("|eff-1| <= %f   for CHIP %i  in HS  %i",fGeneralThresholds[iqual],hs,chip);
356   TCanvas *c[3];
357   
358   for(Int_t i=0; i<2; i++) {
359     c[i] = new TCanvas(Form("c%i",i+1),Form(" %i DIM plots. %s ",i+1,titles.Data()),1200,800); 
360     c[i]->Divide(2,2);
361   }
362   
363   
364   for(Int_t idim =0; idim<2; idim++){
365     for(Int_t k=1; k<5; k++){
366       c[idim]->cd(k);
367       
368       TH1D *h1 =0x0;
369       TH2D *h2=0x0;       
370       if(idim == 0) {
371         h1 = hn->Projection(k-1);
372         if(!h1) {
373          printf("no histogram!!...\n\n\n");
374         } else {
375         h1->SetXTitle(Form("DAC %i  ( %s )",k-1,dacname[k-1].Data()));
376         h1->SetYTitle("entries (eff within thresholds)");
377         h1->Draw();
378         }
379       } 
380       
381       if(idim==1) {      
382         if(k<4){
383           h2  = hn->Projection(k-1,k);
384           h2->SetXTitle(Form("DAC %i ( %s )",k,dacname[k].Data())); h2->SetYTitle(Form("DAC %i  ( %s )",k-1,dacname[k-1].Data()));
385           h2->SetTitleOffset(2,"Y"); h2->SetTitleOffset(1.5,"X");
386         } else {
387           h2 = hn->Projection(0,3);
388           h2->SetXTitle(Form("DAC %i ( %s )",3,dacname[3].Data())); h2->SetYTitle(Form("DAC %i  ( %s )",0, dacname[0].Data()));
389           h2->SetTitleOffset(2,"Y"); h2->SetTitleOffset(1.5,"X"); 
390         }    
391         
392         h2->Draw("lego2");
393       }  
394     }// k loop
395     
396     c[idim]->SaveAs(Form("c%iDIM_%i_%i_%i.eps",idim,iqual,hs,chip));   
397   }// idim loop  
398 }
399 //-----------------------------------------------------
400 TArrayI AliITSOnlineSPDfoAnalyzer::ChooseDACValues(Int_t hs, Int_t chip) const
401 {
402   //
403   // here the "best" 4 dacs are chosen. The most present are taken. 
404   // If the n-tuple does not correspond to a meaningful entry, the
405   // closest value to the mean point in the n-dimensional histogram
406   // is taken.
407   
408   TH1D *tmp[5];
409   for(Int_t jj=0;jj<5;jj++)tmp[jj]=NULL;
410   if(fNdims > 5) printf("AliITSOnlineSPDfoAnalyzer::ChooseDACValues -> N. of dimensions are more than expected! Break! \n");
411   TArrayI dacs(fNdims+1);
412   
413   for(Int_t i=0; i<fNdims+1; i++) dacs.AddAt(-1,i);
414   
415   for(Int_t iqual =0; iqual < kNqualityFlags; iqual++){
416     if(!fNh[iqual][hs][chip]) continue;    
417     if(fNh[iqual][hs][chip]->GetEntries()==0) continue;
418     for(Int_t idim =0; idim<fNdims; idim++){
419       if(dacs.At(idim)>=0) continue;
420       tmp[idim] = fNh[iqual][hs][chip]->Projection(idim);
421       dacs.AddAt((Int_t)tmp[idim]->GetBinLowEdge(tmp[idim]->GetMaximumBin()+1),idim);
422       Int_t bin=-1;
423       if(fFOHandler->GetFOscanInfo()->GetDACindex(idim)==fFOHandler->kIdPreVTH && CorrectPreVTHChioce(tmp[idim],bin)) {
424        dacs.AddAt((Int_t)tmp[idim]->GetBinLowEdge(bin+1),idim);
425       }
426       dacs.AddAt(iqual,fNdims);
427     }//idim
428   }//iqual
429   
430   if(!IsExisting(dacs,hs,chip)  && dacs.At(fNdims)>-1) {   
431    TArrayI centraldacs = GetCentralDACS(dacs.At(fNdims),hs,chip,tmp);
432     dacs = centraldacs;
433   }
434   return dacs; 
435 }
436 //_____________________________________________________
437 Bool_t AliITSOnlineSPDfoAnalyzer::IsExisting(TArrayI dacs,Int_t hs, Int_t chip) const
438 {
439   //
440   // check the N-tuple corresponds to a real one
441   //  
442   
443   if(dacs.At(fNdims)<0) return kFALSE;
444   Bool_t isOk = kFALSE;
445   
446   Int_t size = dacs.GetSize()-1;
447   Double_t *entry = new Double_t[size];
448   for(Int_t i=0; i<size; i++) entry[i] = dacs.At(i);
449   Int_t checkbin = fNh[dacs.At(dacs.GetSize()-1)][hs][chip]->GetBin(entry,kFALSE); // kFALSE does not allocate another bin
450   if(checkbin > -1) isOk = kTRUE;
451   delete [] entry;
452   return isOk; 
453 }
454 //-----------------------------------------------------------
455 TArrayI AliITSOnlineSPDfoAnalyzer::GetCentralDACS(Int_t qualityflag, Int_t hs, Int_t chip, TH1D **hd) const
456 {
457   //
458   // This method gets the DAC values which are closest to the mean point in the N-dim histogram
459   // It is called when the most probable DAC set does not correspond to a real entry in the N-dim
460   // histogram.
461   //
462  
463   TArrayI dacs(fNdims+1);
464    
465   Double_t *mean=new Double_t[fNdims];
466   Int_t *goodbins=new Int_t[fNdims];
467   Double_t distance = 9999999;
468     for(Int_t i=0; i<fNdims ;i++){ 
469     mean[i]=hd[i]->GetMean();
470     goodbins[i]=0;
471     dacs.AddAt(-1,i);
472   }
473   
474   Int_t *bins = new Int_t[fNdims];
475   Double_t *val=new Double_t[fNdims];
476   
477   for(Int_t in=0; in< fNh[qualityflag][hs][chip]->GetNbins() ; in++){
478     
479     fNh[qualityflag][hs][chip]->GetBinContent(in,bins);
480     Double_t r2 = 0;  
481     for(Int_t j=0; j<fNdims; j++) {
482      val[j] = hd[j]->GetBinCenter(bins[j]);
483       r2+=TMath::Power(val[j]-mean[j],2);
484     }
485     
486     if(r2<distance) {
487       distance = r2;
488       fNh[qualityflag][hs][chip]->GetBinContent(in,goodbins);
489     }    
490   }
491   
492   
493   for(Int_t k=0; k<fNdims; k++) {  
494    dacs.AddAt((Int_t)(hd[k]->GetBinCenter(goodbins[k]) + 0.5*hd[k]->GetBinWidth(goodbins[k])),k);
495   }
496  
497   dacs.AddAt(qualityflag,fNdims);
498   
499   
500   delete [] mean;
501   delete [] goodbins;
502   delete [] bins;
503   delete [] val;
504   return dacs;
505 }
506 //-------------------------------------------------------
507 void AliITSOnlineSPDfoAnalyzer::ReadParamsFromLocation(const Char_t *dirName) 
508 {
509   //
510   // opens file (default name) in dir dirName and reads parameters from it
511   // The file should be in database
512   //
513   
514   TString paramsFileName = Form("./%s/focalib_params.txt",dirName);
515   
516   ifstream paramsFile;
517   paramsFile.open(paramsFileName, ifstream::in);
518   if (paramsFile.fail()) {
519     printf("No config file (%s) present. Using default tuning parameters.\n",paramsFileName.Data());
520   }
521   else {
522     while(1) {
523       Char_t paramN[50];
524       Char_t paramV[50];
525       paramsFile >> paramN;
526       if (paramsFile.eof()) break;
527       paramsFile >> paramV;
528       SetParam(paramN,paramV);
529       if (paramsFile.eof()) break;
530     }
531     paramsFile.close();
532   }
533 }
534 //---------------------------------------------------------
535 void AliITSOnlineSPDfoAnalyzer::SetParam(const Char_t *pname, const Char_t *pval) 
536 {
537   //
538   // sets a single parameter when reading from a file in the database
539   //
540   
541   TString name = pname;
542   TString val = pval;
543   
544   
545   if (name.CompareTo("fGeneralThresholds0")==0) {
546     fGeneralThresholds[0] = val.Atof();
547   }
548   else if (name.CompareTo("fGeneralThresholds1")==0) {
549     fGeneralThresholds[1] = val.Atof();
550   }
551   else if (name.CompareTo("fGeneralThresholds2")==0) {
552     fGeneralThresholds[2] = val.Atof();
553   }
554   else {
555     Error("AliITSOnlineSPDscanAnalyzer::SetParam","Parameter %s in configuration file unknown.",name.Data());
556   }
557 }
558 //--------------------------------------------------------
559 Bool_t AliITSOnlineSPDfoAnalyzer::CorrectPreVTHChioce(const TH1D *h,Int_t &bin) const
560 {
561   //
562   // Checks if more maxima of the same height are present in the pre_VTH case
563   //
564   
565   
566   Int_t maxbin = (Int_t)h->GetMaximumBin();
567   Double_t maxentries = h->GetBinContent(maxbin);
568   
569   Int_t binindex = -1;
570   Bool_t check=kFALSE;
571
572   for(Int_t i=0; i< h->GetNbinsX(); i++){    
573      if(h->GetBinContent(i) == maxentries){
574       if(binindex <= i) binindex =i;
575     }
576   }
577   
578   
579   if(binindex>-1) {
580     bin=binindex;
581     check = kTRUE; 
582   }
583   
584   
585   return check; 
586 }