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