1 /**************************************************************************
2 * Copyright(c) 2008-2010, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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 //
24 ////////////////////////////////////////////////////////////
35 #include <THnSparse.h>
39 #include "AliITSOnlineSPDfoChipConfig.h"
40 #include "AliITSOnlineSPDfoChip.h"
41 #include "AliITSOnlineSPDfoInfo.h"
42 #include "AliITSOnlineSPDfo.h"
43 #include "AliITSOnlineSPDfoAnalyzer.h"
46 AliITSOnlineSPDfoAnalyzer::AliITSOnlineSPDfoAnalyzer(const TString fileName, Bool_t readFromGridFile) :
53 fHighOccupancyCheck(kTRUE)
57 for(Int_t iqual =0; iqual<kNqualityFlags; iqual++) fGeneralThresholds[iqual]=0;
59 for (UInt_t chipNr=0; chipNr<10; chipNr++) {
60 for (UInt_t hs=0; hs<6; hs++) {
61 for(Int_t i =0; i< kNqualityFlags; i++) fNh[i][hs][chipNr]=NULL;
65 Init(readFromGridFile);
67 //-------------------------------------------------------------------
68 AliITSOnlineSPDfoAnalyzer::AliITSOnlineSPDfoAnalyzer(const AliITSOnlineSPDfoAnalyzer& foan) :
69 fFileName(foan.fFileName),
74 fFOHandler(foan.fFOHandler),
75 fHighOccupancyCheck(foan.fHighOccupancyCheck)
78 // copy constructor, only copies the filename and params (not the processed data
81 for(Int_t iqual =0; iqual<3; iqual++) fGeneralThresholds[iqual] =foan.fGeneralThresholds[iqual];
82 for (UInt_t chipNr=0; chipNr<11; chipNr++) {
83 for (UInt_t hs=0; hs<6; hs++) {
84 for(Int_t i=0; i<kNqualityFlags;i++) fNh[i][hs][chipNr]=NULL;
90 //-------------------------------------------------------------------
91 AliITSOnlineSPDfoAnalyzer::~AliITSOnlineSPDfoAnalyzer()
97 for (UInt_t hs=0; hs<6; hs++) {
98 for (UInt_t chipNr=0; chipNr<11; chipNr++) {
99 for(Int_t i=0; i<kNqualityFlags ; i++ ) if(fNh[i][hs][chipNr]!=NULL) delete fNh[i][hs][chipNr];
104 //-------------------------------------------------------------------
105 AliITSOnlineSPDfoAnalyzer& AliITSOnlineSPDfoAnalyzer::operator=(const AliITSOnlineSPDfoAnalyzer& foan)
107 // assignment operator, only copies the filename and params (not the processed data)
109 for (UInt_t hs=0; hs<6; hs++) {
110 for (UInt_t chipNr=0; chipNr<11; chipNr++) {
111 for(Int_t i=0; i<kNqualityFlags ; i++ ) if(fNh[i][hs][chipNr]!=NULL) delete fNh[i][hs][chipNr];
115 fFileName=foan.fFileName;
116 fHighOccupancyCheck=foan.fHighOccupancyCheck;
121 //-------------------------------------------------------------------
122 void AliITSOnlineSPDfoAnalyzer::Init(Bool_t readFromGridFile)
125 // first checks type of container and then initializes container obj
127 if (!readFromGridFile) {
128 TFile* p0 = TFile::Open(fFileName.Data());
130 printf("no file open!");
134 fFOHandler = new AliITSOnlineSPDfo();
135 fFOHandler->SetFile(fFileName);
140 //-------------------------------------------------------------------
141 void AliITSOnlineSPDfoAnalyzer::SetGeneralThresholds(Float_t thre[3])
143 // here the settings for the 3 quality flags are defined
144 for(Int_t i=0; i< 3; i++) fGeneralThresholds[i]=thre[i];
146 //-------------------------------------------------------------------
147 void AliITSOnlineSPDfoAnalyzer::SetNdimensions()
150 // here the axis of the N-dim histograms are setted
153 printf("no fo object. Exiting...\n");
157 TArrayI array = fFOHandler->GetDACscanParams();
158 fNdims = array.GetSize()/3;
159 fNbins = new Int_t[fNdims];
160 fXmin = new Double_t[fNdims];
161 fXmax = new Double_t[fNdims];
162 for(Int_t i=0; i< fNdims; i++){
163 fXmin[i] = array.At(3*i)-0.25;
164 fXmax[i] = array.At(3*i+1)+0.25;
165 fNbins[i] = (Int_t)((fXmax[i]-fXmin[i])/0.5)+1;//to avoid Int->Double conversion problems when checking results.
168 //-------------------------------------------------------------------
169 void AliITSOnlineSPDfoAnalyzer::BuildTHnSparse(Int_t ihs, Int_t ichip)
172 // here the N-dim histogram is booked per chip in one HS
175 if(!fNdims) SetNdimensions();
177 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);
179 //-------------------------------------------------------------------
180 Int_t AliITSOnlineSPDfoAnalyzer::Select(const AliITSOnlineSPDfoChip *chip) const
183 // Selects the DAC values if in the chip: the I configuration corresponds to
184 // 0% efficiency ( no FO response case). All the others shoud be within the
185 // predefined thresholds
187 if(!fFOHandler->GetFOscanInfo()) {
188 printf("no general information object in the file. Exiting...\n");
192 Double_t npulses = (fFOHandler->GetFOscanInfo())->GetNumTriggers();
195 Info("AliITSOnlineSPDfoAnalyzer::Select","no trigger pulses set. Exiting...");
199 TObjArray *array = chip->GetChipConfigInfo();
201 printf("No measurement array found in the chip!!\n");
209 Int_t processedconfigurations = chip->GetNumberOfChipConfigs();
214 for(Int_t isteps =0; isteps < processedconfigurations; isteps++){
216 Int_t matrixId = ((AliITSOnlineSPDfoChipConfig*)array->At(isteps))->GetChipConfigMatrixId();
217 counts = (Float_t)(((AliITSOnlineSPDfoChipConfig*)array->At(isteps))->GetChipConfigCounter());
218 if(matrixId==0 && counts > 0) return -1;
219 if(fHighOccupancyCheck && matrixId ==6) continue;
221 Float_t efficiency = counts/npulses;
224 Int_t response = IsSelected(efficiency);
226 if(quality < response) quality = response;
233 //-----------------------------------------------------
234 Int_t AliITSOnlineSPDfoAnalyzer::IsSelected(Float_t eff) const
237 // returns the quality of the selection
240 for(Int_t i=0; i<3; i++){
241 if(eff <= 1.+ fGeneralThresholds[i] && eff >= 1. - fGeneralThresholds[i] ) return i;
245 //----------------------------------------------------
246 void AliITSOnlineSPDfoAnalyzer::Process()
249 // The procedure is the following:
250 // - DAC 4-tuples are checked
251 // - if the 4-tuple survives the selection, the chip-related histograms are filled.
252 // (- Per each histogram the mean values of each axis are taken)
256 Warning("AliITSOnlineSPDfoAnalyzer::Process","no fo object. Exiting.. \n");
262 TIter iter((fFOHandler->GetFile())->GetListOfKeys());
263 while ((key = (TKey*)(iter.Next()))) {
264 TString classname = key->GetClassName();
265 if(classname.Contains("OnlineSPD")) break;
267 TObjArray *array = (TObjArray*)(fFOHandler->GetFile())->Get(key->GetName()); // array of chips corresponding to the DACS (size 1-60)
269 printf("no array found! Exiting...\n");
273 dacvalues = fFOHandler->GetDACvaluesD(key->GetName(), GetFOHandler()->GetFOscanInfo()->GetNumDACindex());
275 for(Int_t i=0; i< array->GetSize(); i++){
276 AliITSOnlineSPDfoChip *chip = (AliITSOnlineSPDfoChip *)array->At(i);
279 Int_t hs = chip->GetActiveHS();
280 Int_t chipid = chip->GetChipId();
281 Int_t quality = Select(chip);
282 if(quality<0) continue;
283 if(!fNh[quality][hs][chipid]) BuildTHnSparse(hs,chipid);
284 fNh[quality][hs][chipid]->Fill(dacvalues);
288 //---------------------------------------------
289 void AliITSOnlineSPDfoAnalyzer::WriteToFile(TString outputfile)
292 // The 4-dim histograms are stored into a file
295 TFile * f = TFile::Open(outputfile.Data(),"recreate");
296 for(Int_t ihs =0; ihs < 6; ihs++) {
297 for(Int_t ichip =0; ichip < 10; ichip++){
298 for(Int_t i=0; i<kNqualityFlags ; i++ ) {
299 if(fNh[i][ihs][ichip]) f->WriteObjectAny(fNh[i][ihs][ichip],"THnSparse",Form("h%i_hs%i_chip%i",i,ihs,ichip));
306 //---------------------------------------------
307 void AliITSOnlineSPDfoAnalyzer::CheckResults(TString filename, Int_t hs, Int_t ichip, Int_t iqual) const
310 //The chip related 4-dim histograms are produced and stored into eps files
313 TFile *f = TFile::Open(filename.Data());
315 Info("AliITSOnlineSPDfoAnalyzer::CehckResults","no file open!. Exiting...\n");
322 TIter iter(f->GetListOfKeys());
323 while((obj=iter.Next())){
324 TString name = obj->GetName();
325 hn=(THnSparse*)f->Get(name.Data());
326 if(name.Contains(Form("hs%i",hs)) && name.Contains(Form("chip%i",ichip)) && name.Contains(Form("h%i",iqual))) GetCanvases(hn,hs,ichip,iqual);
329 //-----------------------------------------------------------------------------------------------
330 void AliITSOnlineSPDfoAnalyzer::GetCanvases(const THnSparse *hn,Int_t hs, Int_t chip, Int_t iqual) const
333 // 1-Dim and 2 Dim Projections of 4-dim histograms are visualized in canvases per quality selection
337 if(!hn) {printf("no thnsparse...exiting!\n"); return;}
339 gStyle->SetPalette(1);
343 if( (fFOHandler->GetFOscanInfo())->GetDACindex(0) == 20 ) dacname[0] = "FOPOL";
344 else dacname[0] = "possible DAC-name/ DAC-register-number mismatch";
345 if( (fFOHandler->GetFOscanInfo())->GetDACindex(1) == 17 ) dacname[1] = "CONVPOL";
346 else dacname[1] = "possible DAC-name/ DAC-register-number mismatch";
347 if( (fFOHandler->GetFOscanInfo())->GetDACindex(2) == 16 ) dacname[2] = "COMPREF";
348 else dacname[2] = "possible DAC-name/ DAC-register-number mismatch";
349 if( (fFOHandler->GetFOscanInfo())->GetDACindex(3) == 39 ) dacname[3] = "pre_VTH";
350 else dacname[3] = "possible DAC-name/ DAC-register-number mismatch";
352 TString titles = Form("|eff-1| <= %f for CHIP %i in HS %i",fGeneralThresholds[iqual],hs,chip);
355 for(Int_t i=0; i<2; i++) {
356 c[i] = new TCanvas(Form("c%i",i+1),Form(" %i DIM plots. %s ",i+1,titles.Data()),1200,800);
361 for(Int_t idim =0; idim<2; idim++){
362 for(Int_t k=1; k<5; k++){
368 h1 = hn->Projection(k-1);
369 if(!h1) printf("no histogram!!...\n\n\n");
370 h1->SetXTitle(Form("DAC %i ( %s )",k-1,dacname[k-1].Data()));
371 h1->SetYTitle("entries (eff within thresholds)");
377 h2 = hn->Projection(k-1,k);
378 h2->SetXTitle(Form("DAC %i ( %s )",k,dacname[k].Data())); h2->SetYTitle(Form("DAC %i ( %s )",k-1,dacname[k-1].Data()));
379 h2->SetTitleOffset(2,"Y"); h2->SetTitleOffset(1.5,"X");
381 h2 = hn->Projection(0,3);
382 h2->SetXTitle(Form("DAC %i ( %s )",3,dacname[3].Data())); h2->SetYTitle(Form("DAC %i ( %s )",0, dacname[0].Data()));
383 h2->SetTitleOffset(2,"Y"); h2->SetTitleOffset(1.5,"X");
390 c[idim]->SaveAs(Form("c%iDIM_%i_%i_%i.eps",idim,iqual,hs,chip));
393 //-----------------------------------------------------
394 TArrayI AliITSOnlineSPDfoAnalyzer::ChooseDACValues(Int_t hs, Int_t chip) const
397 // here the "best" 4 dacs are chosen. The most present are taken.
398 // If the n-tuple does not correspond to a meaningful entry, the
399 // closest value to the mean point in the n-dimensional histogram
403 if(fNdims > 5) printf("AliITSOnlineSPDfoAnalyzer::ChooseDACValues -> N. of dimensions are more than expected! Break! \n");
404 TArrayI dacs(fNdims+1);
406 for(Int_t i=0; i<fNdims+1; i++) dacs.AddAt(-1,i);
408 for(Int_t iqual =0; iqual < kNqualityFlags; iqual++){
409 if(!fNh[iqual][hs][chip]) continue;
410 if(fNh[iqual][hs][chip]->GetEntries()==0) continue;
411 for(Int_t idim =0; idim<fNdims; idim++){
412 if(dacs.At(idim)>=0) continue;
413 tmp[idim] = fNh[iqual][hs][chip]->Projection(idim);
414 dacs.AddAt((Int_t)tmp[idim]->GetBinLowEdge(tmp[idim]->GetMaximumBin()+1),idim);
416 if(fFOHandler->GetFOscanInfo()->GetDACindex(idim)==fFOHandler->kIdPreVTH && CorrectPreVTHChioce(tmp[idim],bin)) {
417 dacs.AddAt((Int_t)tmp[idim]->GetBinLowEdge(bin+1),idim);
419 dacs.AddAt(iqual,fNdims);
423 if(!IsExisting(dacs,hs,chip) && dacs.At(fNdims)>-1) {
424 TArrayI centraldacs = GetCentralDACS(dacs.At(fNdims),hs,chip,tmp);
429 //_____________________________________________________
430 Bool_t AliITSOnlineSPDfoAnalyzer::IsExisting(TArrayI dacs,Int_t hs, Int_t chip) const
433 // check the N-tuple corresponds to a real one
436 if(dacs.At(fNdims)<0) return kFALSE;
437 Bool_t isOk = kFALSE;
439 Int_t size = dacs.GetSize()-1;
440 Double_t *entry = new Double_t[size];
441 for(Int_t i=0; i<size; i++) entry[i] = dacs.At(i);
442 Int_t checkbin = fNh[dacs.At(dacs.GetSize()-1)][hs][chip]->GetBin(entry,kFALSE); // kFALSE does not allocate another bin
443 if(checkbin > -1) isOk = kTRUE;
447 //-----------------------------------------------------------
448 TArrayI AliITSOnlineSPDfoAnalyzer::GetCentralDACS(Int_t qualityflag, Int_t hs, Int_t chip, TH1D **hd) const
451 // This method gets the DAC values which are closest to the mean point in the N-dim histogram
452 // It is called when the most probable DAC set does not correspond to a real entry in the N-dim
456 TArrayI dacs(fNdims+1);
458 Double_t *mean=new Double_t[fNdims];
459 Int_t *goodbins=new Int_t[fNdims];
460 Double_t distance = 9999999;
461 for(Int_t i=0; i<fNdims ;i++){
462 mean[i]=hd[i]->GetMean();
467 Int_t *bins = new Int_t[fNdims];
468 Double_t *val=new Double_t[fNdims];
470 for(Int_t in=0; in< fNh[qualityflag][hs][chip]->GetNbins() ; in++){
472 fNh[qualityflag][hs][chip]->GetBinContent(in,bins);
474 for(Int_t j=0; j<fNdims; j++) {
475 val[j] = hd[j]->GetBinCenter(bins[j]);
476 r2+=TMath::Power(val[j]-mean[j],2);
481 fNh[qualityflag][hs][chip]->GetBinContent(in,goodbins);
486 for(Int_t k=0; k<fNdims; k++) {
487 dacs.AddAt((Int_t)(hd[k]->GetBinCenter(goodbins[k]) + 0.5*hd[k]->GetBinWidth(goodbins[k])),k);
490 dacs.AddAt(qualityflag,fNdims);
499 //-------------------------------------------------------
500 void AliITSOnlineSPDfoAnalyzer::ReadParamsFromLocation(const Char_t *dirName)
503 // opens file (default name) in dir dirName and reads parameters from it
504 // The file should be in database
507 TString paramsFileName = Form("./%s/focalib_params.txt",dirName);
510 paramsFile.open(paramsFileName, ifstream::in);
511 if (paramsFile.fail()) {
512 printf("No config file (%s) present. Using default tuning parameters.\n",paramsFileName.Data());
518 paramsFile >> paramN;
519 if (paramsFile.eof()) break;
520 paramsFile >> paramV;
521 SetParam(paramN,paramV);
522 if (paramsFile.eof()) break;
527 //---------------------------------------------------------
528 void AliITSOnlineSPDfoAnalyzer::SetParam(const Char_t *pname, const Char_t *pval)
531 // sets a single parameter when reading from a file in the database
534 TString name = pname;
538 if (name.CompareTo("fGeneralThresholds0")==0) {
539 fGeneralThresholds[0] = val.Atof();
541 else if (name.CompareTo("fGeneralThresholds1")==0) {
542 fGeneralThresholds[1] = val.Atof();
544 else if (name.CompareTo("fGeneralThresholds2")==0) {
545 fGeneralThresholds[2] = val.Atof();
548 Error("AliITSOnlineSPDscanAnalyzer::SetParam","Parameter %s in configuration file unknown.",name.Data());
551 //--------------------------------------------------------
552 Bool_t AliITSOnlineSPDfoAnalyzer::CorrectPreVTHChioce(const TH1D *h,Int_t &bin) const
555 // Checks if more maxima of the same height are present in the pre_VTH case
559 Int_t maxbin = (Int_t)h->GetMaximumBin();
560 Double_t maxentries = h->GetBinContent(maxbin);
565 for(Int_t i=0; i< h->GetNbinsX(); i++){
566 if(h->GetBinContent(i) == maxentries){
567 if(binindex <= i) binindex =i;