]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TRD/AliTRDCalibPadStatus.cxx
New DAQ calibration DAs by Raphaelle
[u/mrichter/AliRoot.git] / TRD / AliTRDCalibPadStatus.cxx
CommitLineData
170c35f1 1/**************************************************************************
2 * Copyright(c) 1998-1999, 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 example: fill pedestal with gausschen noise
20 AliTRDCalibPadStatus ped;
21 ped.TestEvent(numberofevent);
22 // Method without histo
23 //////////////////////////
24 ped.Analyse();
25 //Create the histo of the AliTRDCalROC
26 TH2F * histo2dm = ped.GetCalRocMean(0,kFALSE)->MakeHisto2D();
27 histo2dm->Scale(10.0);
28 TH1F * histo1dm = ped.GetCalRocMean(0,kFALSE)->MakeHisto1D();
29 histo1dm->Scale(10.0);
30 TH2F * histo2ds = ped.GetCalRocSquares(0,kFALSE)->MakeHisto2D();
31 histo2ds->Scale(10.0);
32 TH1F * histo1ds = ped.GetCalRocSquares(0,kFALSE)->MakeHisto1D();
33 histo1ds->Scale(10.0)
34 //Draw output;
35 TCanvas* c1 = new TCanvas;
36 c1->Divide(2,2);
37 c1->cd(1);
38 histo2dm->Draw("colz");
39 c1->cd(2);
40 histo1dm->Draw();
41 c1->cd(3);
42 histo2ds->Draw("colz");
43 c1->cd(4);
44 histo1ds->Draw();
45// Method with histo
46/////////////////////////
47ped.AnalyseHisto();
48//Take the histo
49TH1F * histo = ped.GetHisto(31);
50histo->SetEntries(1);
51histo->Draw();
52
53*/
54
55//Root includes
56#include <TObjArray.h>
57#include <TH1F.h>
58#include <TH2F.h>
59#include <TString.h>
60#include <TMath.h>
61#include <TF1.h>
62#include <TRandom.h>
63#include <TDirectory.h>
64#include <TFile.h>
65#include <AliMathBase.h>
66#include "TTreeStream.h"
67//AliRoot includes
68#include "AliRawReader.h"
69#include "AliRawReaderRoot.h"
70#include "AliRawReaderDate.h"
71#include "AliTRDRawStream.h"
72#include "./Cal/AliTRDCalROC.h"
73#include "./Cal/AliTRDCalPadStatus.h"
74#include "./Cal/AliTRDCalSingleChamberStatus.h"
75#include "AliTRDarrayF.h"
76#include "AliTRDCommonParam.h"
77
78
79#ifdef ALI_DATE
80#include "event.h"
81#endif
82
83//header file
84#include "AliTRDCalibPadStatus.h"
85
86
87
88ClassImp(AliTRDCalibPadStatus) /*FOLD00*/
89
90//_____________________________________________________________________
91AliTRDCalibPadStatus::AliTRDCalibPadStatus() : /*FOLD00*/
92 TObject(),
93 fAdcMin(0),
94 fAdcMax(20),
95 fDetector(-1),
96 fCalArrayEntries(540),
97 fCalArrayMean(540),
98 fCalArraySquares(540),
99 fCalRocArrayMean(540),
100 fCalRocArrayRMS(540),
101 fHistoArray(540),
102 fCalEntries(0x0),
103 fCalMean(0x0),
104 fCalSquares(0x0)
105{
106 //
107 // default constructor
108 //
109}
110
111//_____________________________________________________________________
112AliTRDCalibPadStatus::AliTRDCalibPadStatus(const AliTRDCalibPadStatus &ped) : /*FOLD00*/
113 TObject(ped),
114 fAdcMin(ped.GetAdcMin()),
115 fAdcMax(ped.GetAdcMax()),
116 fDetector(ped.fDetector),
117 fCalArrayEntries(540),
118 fCalArrayMean(540),
119 fCalArraySquares(540),
120 fCalRocArrayMean(540),
121 fCalRocArrayRMS(540),
122 fHistoArray(540),
123 fCalEntries(0x0),
124 fCalMean(0x0),
125 fCalSquares(0x0)
126{
127 //
128 // copy constructor
129 //
130 for (Int_t idet = 0; idet < 540; idet++){
131 const AliTRDarrayF *calEntries = (AliTRDarrayF*)ped.fCalArrayEntries.UncheckedAt(idet);
132 const AliTRDarrayF *calMean = (AliTRDarrayF*)ped.fCalArrayMean.UncheckedAt(idet);
133 const AliTRDarrayF *calSquares = (AliTRDarrayF*)ped.fCalArraySquares.UncheckedAt(idet);
134 const AliTRDCalROC *calRocMean = (AliTRDCalROC*)ped.fCalRocArrayMean.UncheckedAt(idet);
135 const AliTRDCalROC *calRocRMS = (AliTRDCalROC*)ped.fCalRocArrayRMS.UncheckedAt(idet);
136 const TH2F *hped = (TH2F*)ped.fHistoArray.UncheckedAt(idet);
137
138 if ( calEntries != 0x0 ) fCalArrayEntries.AddAt(new AliTRDarrayF(*calEntries), idet);
139 if ( calMean != 0x0 ) fCalArrayMean.AddAt(new AliTRDarrayF(*calMean), idet);
140 if ( calSquares != 0x0 ) fCalArraySquares.AddAt(new AliTRDarrayF(*calSquares), idet);
141 if ( calRocMean != 0x0 ) fCalRocArrayMean.AddAt(new AliTRDCalROC(*calRocMean), idet);
142 if ( calRocRMS != 0x0 ) fCalRocArrayRMS.AddAt(new AliTRDCalROC(*calRocRMS), idet);
143
144 if ( hped != 0x0 ){
145 TH2F *hNew = new TH2F(*hped);
146 hNew->SetDirectory(0);
147 fHistoArray.AddAt(hNew,idet);
148 }
149
150 }
151}
152//_____________________________________________________________________
153AliTRDCalibPadStatus& AliTRDCalibPadStatus::operator = (const AliTRDCalibPadStatus &source)
154{
155 //
156 // assignment operator
157 //
158 if (&source == this) return *this;
159 new (this) AliTRDCalibPadStatus(source);
160
161 return *this;
162}
163//_____________________________________________________________________
164AliTRDCalibPadStatus::~AliTRDCalibPadStatus() /*FOLD00*/
165{
166 //
167 // destructor
168 //
169}
170//_____________________________________________________________________
171Int_t AliTRDCalibPadStatus::Update(const Int_t icdet, /*FOLD00*/
172 const Int_t icRow,
173 const Int_t icCol,
174 const Int_t csignal,
175 const Int_t crowMax)
176{
177 //
178 // Signal filling methode
179 //
180 if ( (csignal>fAdcMax) || (csignal<fAdcMin) ) return 0;
181
182 if(fDetector != icdet){
183 fCalEntries = ((AliTRDarrayF *)GetCalEntries(icdet,kTRUE));
184 fCalMean = ((AliTRDarrayF *)GetCalMean(icdet,kTRUE));
185 fCalSquares = ((AliTRDarrayF *)GetCalSquares(icdet,kTRUE));
186 }
187
188 Float_t entries = fCalEntries->At(icRow+icCol*crowMax);
189 Float_t mean = fCalMean->At(icRow+icCol*crowMax);
190 Float_t squares = fCalSquares->At(icRow+icCol*crowMax);
191
192 Float_t entriesn = entries+1.0;
193 fCalEntries->AddAt(entriesn,(icRow+icCol*crowMax));
194 Float_t meann = (mean*entries+((Float_t)(csignal+0.5)))/entriesn;
195 fCalMean->AddAt(meann,icRow+icCol*crowMax);
196 Float_t squaresn = ((squares*entries)+(((Float_t)(csignal+0.5))*((Float_t)(csignal+0.5))))/entriesn;
197 fCalSquares->AddAt(squaresn,icRow+icCol*crowMax);
198
199 //printf("icdet %d, icRow %d, icCol %d, csignal %d, crowMax %d\n",icdet,icRow,icCol,csignal,crowMax);
200 //printf("entries %f, mean %f, squares %f\n",entriesn,meann,squaresn);
201
202 fDetector = icdet;
203
204 return 0;
205}
206//_____________________________________________________________________
207Int_t AliTRDCalibPadStatus::UpdateHisto(const Int_t icdet, /*FOLD00*/
208 const Int_t icRow,
209 const Int_t icCol,
210 const Int_t csignal,
211 const Int_t crowMax)
212{
213 //
214 // Signal filling methode
215 //
216 Int_t nbchannel = icRow+icCol*crowMax;
217
218 // fast filling methode.
219 // Attention: the entry counter of the histogram is not increased
220 // this means that e.g. the colz draw option gives an empty plot
221 Int_t bin = 0;
222 if ( !(((Int_t)csignal>fAdcMax ) || ((Int_t)csignal<fAdcMin)) )
223 bin = (nbchannel+1)*(fAdcMax-fAdcMin+2)+((Int_t)csignal-fAdcMin+1);
224
225 GetHisto(icdet,kTRUE)->GetArray()[bin]++;
226
227 return 0;
228}
229//_____________________________________________________________________
230Bool_t AliTRDCalibPadStatus::ProcessEvent(AliTRDRawStream *rawStream)
231{
232 //
233 // Event Processing loop - AliTRDRawStream
234 //
235
236
237 Bool_t withInput = kFALSE;
238
239 while (rawStream->Next()) {
240
241 Int_t rawversion = rawStream->GetRawVersion(); // current raw version
242 if(rawversion != 2) return kFALSE;
243
244 Int_t idetector = rawStream->GetDet(); // current detector
245 Int_t iRow = rawStream->GetRow(); // current row
246 Int_t iRowMax = rawStream->GetMaxRow(); // current rowmax
247 Int_t iCol = rawStream->GetCol(); // current col
248 Int_t iTimeBin = rawStream->GetTimeBin(); // current time bin
249 Int_t *signal = rawStream->GetSignals(); // current ADC signal
250 Int_t nbtimebin = rawStream->GetNumberOfTimeBins(); // number of time bins read from data
251
252 Int_t fin = TMath::Min(nbtimebin,(iTimeBin+3));
253 Int_t n = 0;
254
255 for(Int_t k = iTimeBin; k < fin; k++){
256 if(signal[n]>0) UpdateHisto(idetector,iRow,iCol,signal[n],iRowMax);
257 n++;
258 }
259
260 withInput = kTRUE;
261 }
262
263 return withInput;
264}
265//_____________________________________________________________________
266Bool_t AliTRDCalibPadStatus::ProcessEvent(AliRawReader *rawReader)
267{
268 //
269 // Event processing loop - AliRawReader
270 //
271
272
273 AliTRDRawStream rawStream(rawReader);
274
275 rawReader->Select("TRD");
276
277 return ProcessEvent(&rawStream);
278}
279//_________________________________________________________________________
280Bool_t AliTRDCalibPadStatus::ProcessEvent(
281#ifdef ALI_DATE
282 eventHeaderStruct *event
283#else
284 eventHeaderStruct* /*event*/
285
286#endif
287 )
288{
289 //
290 // process date event
291 //
292#ifdef ALI_DATE
293 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
294 Bool_t result=ProcessEvent(rawReader);
295 delete rawReader;
296 return result;
297#else
298 Fatal("AliTRDCalibPadStatus", "this class was compiled without DATE");
299 return 0;
300#endif
301
302}
303//_____________________________________________________________________
304Bool_t AliTRDCalibPadStatus::TestEvent(Int_t nevent) /*FOLD00*/
305{
306 //
307 // Test event loop
308 // fill one oroc and one iroc with random gaus
309 //
310
311 AliTRDCommonParam *comParam = AliTRDCommonParam::Instance();
312 if (!comParam) {
313 return 0;
314 }
315 gRandom->SetSeed(0);
316
317 for (Int_t ism=0; ism<18; ism++){
318 for (Int_t ich=0; ich < 5; ich++){
319 for (Int_t ipl=0; ipl < 6; ipl++){
320 for(Int_t irow = 0; irow < comParam->GetRowMax(ipl,ich,ism); irow++){
321 for(Int_t icol = 0; icol < comParam->GetColMax(ipl); icol++){
322 for (Int_t iTimeBin=0; iTimeBin<(30*nevent); iTimeBin++){
323 Int_t signal=(Int_t)(ich+8+gRandom->Gaus(0,1.2));
324 if ( signal>0 )Update((ipl+ich*6+ism*6*5),irow,icol,signal,comParam->GetRowMax(ipl,ich,ism));
325 }
326 }
327 }
328 }
329 }
330 }
331 return kTRUE;
332}
333//_____________________________________________________________________
334Bool_t AliTRDCalibPadStatus::TestEventHisto(Int_t nevent) /*FOLD00*/
335{
336 //
337 // Test event loop
338 // fill one oroc and one iroc with random gaus
339 //
340
341 AliTRDCommonParam *comParam = AliTRDCommonParam::Instance();
342 if (!comParam) {
343 return 0;
344 }
345 gRandom->SetSeed(0);
346
347 for (Int_t ism=0; ism<18; ism++){
348 for (Int_t ich=0; ich < 5; ich++){
349 for (Int_t ipl=0; ipl < 6; ipl++){
350 for(Int_t irow = 0; irow < comParam->GetRowMax(ipl,ich,ism); irow++){
351 for(Int_t icol = 0; icol < comParam->GetColMax(ipl); icol++){
352 for (Int_t iTimeBin=0; iTimeBin<(30*nevent); iTimeBin++){
353 Int_t signal=(Int_t)(ich+8+gRandom->Gaus(0,1.2));
354 if ( signal>0 )UpdateHisto((ipl+ich*6+ism*6*5),irow,icol,signal,comParam->GetRowMax(ipl,ich,ism));
355 }
356 }
357 }
358 }
359 }
360 }
361 return kTRUE;
362}
363//_____________________________________________________________________
364TH2F* AliTRDCalibPadStatus::GetHisto(Int_t det, TObjArray *arr, /*FOLD00*/
365 Int_t nbinsY, Float_t ymin, Float_t ymax,
366 Char_t *type, Bool_t force)
367{
368 //
369 // return pointer to histogram
370 // if force is true create a new histogram if it doesn't exist allready
371 //
372 if ( !force || arr->UncheckedAt(det) )
373 return (TH2F*)arr->UncheckedAt(det);
374
375 // if we are forced and histogram doesn't yes exist create it
376 Char_t name[255], title[255];
377
378 sprintf(name,"hCalib%s%.2d",type,det);
379 sprintf(title,"%s calibration histogram detector %.2d;ADC channel;Channel (pad)",type,det);
380
381 AliTRDCommonParam *comParam = AliTRDCommonParam::Instance();
382 if (!comParam) {
383 return 0x0;
384 }
385 Int_t nbchannels = comParam->GetRowMax(GetPlane(det),GetChamber(det),GetSector(det))*comParam->GetColMax(GetPlane(det));
386
387 // new histogram with calib information. One value for each pad!
388 TH2F* hist = new TH2F(name,title,
389 nbinsY, ymin, ymax,
390 nbchannels,0,nbchannels
391 );
392 hist->SetDirectory(0);
393 arr->AddAt(hist,det);
394 return hist;
395}
396//_____________________________________________________________________
397TH2F* AliTRDCalibPadStatus::GetHisto(Int_t det, Bool_t force) /*FOLD00*/
398{
399 //
400 // return pointer to histogram
401 // if force is true create a new histogram if it doesn't exist allready
402 //
403 TObjArray *arr = &fHistoArray;
404 return GetHisto(det, arr, fAdcMax-fAdcMin, fAdcMin, fAdcMax, "Pedestal", force);
405}
406//_____________________________________________________________________
407AliTRDarrayF* AliTRDCalibPadStatus::GetCalEntries(Int_t det, TObjArray* arr, Bool_t force) /*FOLD00*/
408{
409 //
410 // return pointer to ROC Calibration
411 // if force is true create a new AliTRDarrayF if it doesn't exist allready
412 //
413 if ( !force || arr->UncheckedAt(det) )
414 return (AliTRDarrayF*)arr->UncheckedAt(det);
415
416 // if we are forced and histogram doesn't yes exist create it
417 AliTRDarrayF *croc = new AliTRDarrayF();
418 AliTRDCommonParam *comParam = AliTRDCommonParam::Instance();
419 if (!comParam) {
420 return croc;
421 }
422 Int_t nbpad = comParam->GetRowMax(GetPlane(det),GetChamber(det),GetSector(det))*comParam->GetColMax(GetPlane(det));
423
424 // new AliTRDCalROC. One value for each pad!
425 croc->Expand(nbpad);
426 for(Int_t k = 0; k < nbpad; k++){
427 croc->AddAt(0.0,k);
428 }
429 arr->AddAt(croc,det);
430 return croc;
431}
432//_____________________________________________________________________
433AliTRDarrayF* AliTRDCalibPadStatus::GetCalMean(Int_t det, TObjArray* arr, Bool_t force) /*FOLD00*/
434{
435 //
436 // return pointer to ROC Calibration
437 // if force is true create a new AliTRDarrayF if it doesn't exist allready
438 //
439 if ( !force || arr->UncheckedAt(det) )
440 return (AliTRDarrayF*)arr->UncheckedAt(det);
441
442 // if we are forced and histogram doesn't yes exist create it
443 AliTRDarrayF *croc = new AliTRDarrayF();
444 AliTRDCommonParam *comParam = AliTRDCommonParam::Instance();
445 if (!comParam) {
446 return croc;
447 }
448 Int_t nbpad = comParam->GetRowMax(GetPlane(det),GetChamber(det),GetSector(det))*comParam->GetColMax(GetPlane(det));
449
450 // new AliTRDCalROC. One value for each pad!
451 croc->Expand(nbpad);
452 for(Int_t k = 0; k < nbpad; k++){
453 croc->AddAt(0.0,k);
454 }
455 arr->AddAt(croc,det);
456 return croc;
457}
458//_____________________________________________________________________
459AliTRDarrayF* AliTRDCalibPadStatus::GetCalSquares(Int_t det, TObjArray* arr, Bool_t force) /*FOLD00*/
460{
461 //
462 // return pointer to ROC Calibration
463 // if force is true create a new AliTRDarrayF if it doesn't exist allready
464 //
465 if ( !force || arr->UncheckedAt(det) )
466 return (AliTRDarrayF*)arr->UncheckedAt(det);
467
468 // if we are forced and histogram doesn't yes exist create it
469 AliTRDarrayF *croc = new AliTRDarrayF();
470 AliTRDCommonParam *comParam = AliTRDCommonParam::Instance();
471 if (!comParam) {
472 return croc;
473 }
474 Int_t nbpad = comParam->GetRowMax(GetPlane(det),GetChamber(det),GetSector(det))*comParam->GetColMax(GetPlane(det));
475
476 // new AliTRDCalROC. One value for each pad!
477 croc->Expand(nbpad);
478 for(Int_t k = 0; k < nbpad; k++){
479 croc->AddAt(0.0,k);
480 }
481 arr->AddAt(croc,det);
482 return croc;
483}
484//_____________________________________________________________________
485AliTRDCalROC* AliTRDCalibPadStatus::GetCalRocMean(Int_t det, TObjArray* arr, Bool_t force) /*FOLD00*/
486{
487 //
488 // return pointer to ROC Calibration
489 // if force is true create a new AliTRDCalROC if it doesn't exist allready
490 //
491 if ( !force || arr->UncheckedAt(det) )
492 return (AliTRDCalROC*)arr->UncheckedAt(det);
493
494 // if we are forced and histogram doesn't yes exist create it
495
496 // new AliTRDCalROC. One value for each pad!
497 AliTRDCalROC *croc = new AliTRDCalROC(GetPlane(det),GetChamber(det));
498 arr->AddAt(croc,det);
499 return croc;
500}
501//_____________________________________________________________________
502AliTRDCalROC* AliTRDCalibPadStatus::GetCalRocRMS(Int_t det, TObjArray* arr, Bool_t force) /*FOLD00*/
503{
504 //
505 // return pointer to ROC Calibration
506 // if force is true create a new AliTRDCalROC if it doesn't exist allready
507 //
508 if ( !force || arr->UncheckedAt(det) )
509 return (AliTRDCalROC*)arr->UncheckedAt(det);
510
511 // if we are forced and histogram doesn't yes exist create it
512
513 // new AliTRDCalROC. One value for each pad!
514 AliTRDCalROC *croc = new AliTRDCalROC(GetPlane(det),GetChamber(det));
515 arr->AddAt(croc,det);
516 return croc;
517}
518//_____________________________________________________________________
519AliTRDarrayF* AliTRDCalibPadStatus::GetCalEntries(Int_t sector, Bool_t force) /*FOLD00*/
520{
521 //
522 // return pointer to Carge ROC Calibration
523 // if force is true create a new histogram if it doesn't exist allready
524 //
525 TObjArray *arr = &fCalArrayEntries;
526 return GetCalEntries(sector, arr, force);
527}
528//_____________________________________________________________________
529AliTRDarrayF* AliTRDCalibPadStatus::GetCalMean(Int_t sector, Bool_t force) /*FOLD00*/
530{
531 //
532 // return pointer to Carge ROC Calibration
533 // if force is true create a new histogram if it doesn't exist allready
534 //
535 TObjArray *arr = &fCalArrayMean;
536 return GetCalMean(sector, arr, force);
537}
538//_____________________________________________________________________
539AliTRDarrayF* AliTRDCalibPadStatus::GetCalSquares(Int_t sector, Bool_t force) /*FOLD00*/
540{
541 //
542 // return pointer to Carge ROC Calibration
543 // if force is true create a new histogram if it doesn't exist allready
544 //
545 TObjArray *arr = &fCalArraySquares;
546 return GetCalSquares(sector, arr, force);
547}
548//_____________________________________________________________________
549AliTRDCalROC* AliTRDCalibPadStatus::GetCalRocMean(Int_t sector, Bool_t force) /*FOLD00*/
550{
551 //
552 // return pointer to Carge ROC Calibration
553 // if force is true create a new histogram if it doesn't exist allready
554 //
555 TObjArray *arr = &fCalRocArrayMean;
556 return GetCalRocMean(sector, arr, force);
557}
558//_____________________________________________________________________
559AliTRDCalROC* AliTRDCalibPadStatus::GetCalRocRMS(Int_t sector, Bool_t force) /*FOLD00*/
560{
561 //
562 // return pointer to Carge ROC Calibration
563 // if force is true create a new histogram if it doesn't exist allready
564 //
565 TObjArray *arr = &fCalRocArrayRMS;
566 return GetCalRocRMS(sector, arr, force);
567}
568//_________________________________________________________________________
569void AliTRDCalibPadStatus::Analyse() /*FOLD00*/
570{
571 //
572 // Calcul the rms properly
573 //
574
575 for(Int_t idet = 0; idet < 540; idet++){
576
577 // Take the stuff
578 fCalEntries = ((AliTRDarrayF *)GetCalEntries(idet));
579 fCalMean = ((AliTRDarrayF *)GetCalMean(idet));
580 fCalSquares = ((AliTRDarrayF *)GetCalSquares(idet));
581
582 if(!fCalEntries) continue;
583
584 AliTRDCalROC *calRocMean = ((AliTRDCalROC *)GetCalRocMean(idet,kTRUE));
585 AliTRDCalROC *calRocRMS = ((AliTRDCalROC *)GetCalRocRMS(idet,kTRUE));
586
587 // range channels
588 Int_t channels = calRocMean->GetNchannels();
589
590 for(Int_t ichannels = 0 ; ichannels < channels; ichannels++){
591
592 Float_t entries = fCalEntries->At(ichannels);
593 Float_t mean = fCalMean->At(ichannels);
594 Float_t squares = fCalSquares->At(ichannels);
595
596 Float_t rms = 0.0;
597 if(entries > 0){
598 Double_t rm = TMath::Abs(squares-(mean*mean));
599 rms = TMath::Sqrt(rm);
600 calRocRMS->SetValue(ichannels,rms/10.0);
601 calRocMean->SetValue(ichannels,mean/10.0);
602 }
603 }
604 }
605}
606//_____________________________________________________________________
607void AliTRDCalibPadStatus::AnalyseHisto() /*FOLD00*/
608{
609 //
610 // Calculate calibration constants
611 //
612
613 Int_t nbinsAdc = fAdcMax-fAdcMin;
614
615 TVectorD param(3);
616 TMatrixD dummy(3,3);
617
618 Float_t *array_hP=0;
619
620
621 for (Int_t idet=0; idet<540; idet++){
622 TH2F *hP = GetHisto(idet);
623 if ( !hP ) {
624 continue;
625 }
626
627 AliTRDCalROC *rocMean = GetCalRocMean(idet,kTRUE);
628 AliTRDCalROC *rocRMS = GetCalRocRMS(idet,kTRUE);
629
630 array_hP = hP->GetArray();
631 Int_t nChannels = rocMean->GetNchannels();
632
633 for (Int_t iChannel=0; iChannel<nChannels; iChannel++){
634 Int_t offset = (nbinsAdc+2)*(iChannel+1)+1;
635 Double_t ret = AliMathBase::FitGaus(array_hP+offset,nbinsAdc,fAdcMin,fAdcMax,&param,&dummy);
636 // if the fitting failed set noise and pedestal to 0
637 if ((ret == -4) || (ret == -1)) {
638 param[1]=0.0;
639 param[2]=0.0;
640 }
641 if((param[1]/10.0) > 65534.0) param[1] = 0.0;
642 if((param[2]/10.0) > 65534.0) param[2] = 0.0;
643 rocMean->SetValue(iChannel,param[1]/10.0);
644 rocRMS->SetValue(iChannel,param[2]/10.0);
645 }
646 }
647
648}
649//_______________________________________________________________________________________
650AliTRDCalPadStatus* AliTRDCalibPadStatus::CreateCalPadStatus()
651{
652 //
653 //
654 //
655
656 AliTRDCalPadStatus* obj = new AliTRDCalPadStatus("padstatus", "padstatus");
657
658 for (Int_t idet=0; idet<540; ++idet)
659 {
660 AliTRDCalSingleChamberStatus *calROC = obj->GetCalROC(idet);
661
662 //Take the stuff
663 fCalEntries = ((AliTRDarrayF *)GetCalEntries(idet));
664 AliTRDCalROC *calRocMean = ((AliTRDCalROC *)GetCalRocMean(idet));
665 AliTRDCalROC *calRocRMS = ((AliTRDCalROC *)GetCalRocRMS(idet));
666
667 if ( !calRocMean ) {
668 for(Int_t k = 0; k < calROC->GetNchannels(); k++){
669 calROC->SetStatus(k,AliTRDCalPadStatus::kMasked);
670 }
671 //printf("no fCalRocMean for %d\n",idet);
672 continue;
673 }
674
675 //Range
676 Int_t channels = calROC->GetNchannels();
677
678
679 //Mean
680 Float_t meanentries = 0.0;
681 if(!fCalEntries){
682 if(GetHisto(idet)){
683 meanentries = GetHisto(idet)->GetEntries()/(channels);
684 }
685 }
686 else meanentries = TMath::Mean(channels,((TArrayF *)fCalEntries)->GetArray());
687 //Double_t meanmean = calRocMean->GetMean()*10.0;
688 //Double_t meansquares = calRocRMS->GetMean()*10.0;
689
690
691 for(Int_t ich = 0; ich < channels; ich++){
692
693 Float_t entries = 0.0;
694 if(!fCalEntries){
695 if(GetHisto(idet)){
696 for(Int_t bin = 0; bin < (fAdcMax-fAdcMin); bin++){
697 entries += GetHisto(idet)->GetArray()[(ich+1)*(fAdcMax-fAdcMin+2)+(bin+1)];
698 }
699 }
700 }
701 else entries = fCalEntries->At(ich);
702 //Float_t mean = calRocMean->GetValue(ich)*10.0;
703 Float_t rms = calRocRMS->GetValue(ich)*10.0;
704 if(ich > 1720) printf("rms %f\n",rms);
705
706 //if((entries < 0.3*meanentries) || (TMath::Abs(rms-meansquares) > ((Float_t)(fAdcMax-fAdcMin)/2.0)) || (rms == 0.0)) calROC->SetStatus(ich, AliTRDCalPadStatus::kMasked);
707 if(rms <= 0.01) calROC->SetStatus(ich, AliTRDCalPadStatus::kMasked);
708 }
709 }
710
711 return obj;
712
713}
714//_____________________________________________________________________
715void AliTRDCalibPadStatus::DumpToFile(const Char_t *filename, const Char_t *dir, Bool_t append) /*FOLD00*/
716{
717 //
718 // Write class to file
719 //
720
721 TString sDir(dir);
722 TString option;
723
724 if ( append )
725 option = "update";
726 else
727 option = "recreate";
728
729 TDirectory *backup = gDirectory;
730 TFile f(filename,option.Data());
731 f.cd();
732 if ( !sDir.IsNull() ){
733 f.mkdir(sDir.Data());
734 f.cd(sDir);
735 }
736 this->Write();
737 f.Close();
738
739 if ( backup ) backup->cd();
740}//_____________________________________________________________________________
741Int_t AliTRDCalibPadStatus::GetPlane(Int_t d) const
742{
743 //
744 // Reconstruct the plane number from the detector number
745 //
746
747 return ((Int_t) (d % 6));
748
749}
750
751//_____________________________________________________________________________
752Int_t AliTRDCalibPadStatus::GetChamber(Int_t d) const
753{
754 //
755 // Reconstruct the chamber number from the detector number
756 //
757
758 return ((Int_t) (d % 30) / 6);
759
760}
761//_____________________________________________________________________________
762Int_t AliTRDCalibPadStatus::GetSector(Int_t d) const
763{
764 //
765 // Reconstruct the sector number from the detector number
766 //
767
768 return ((Int_t) (d / 30));
769
770}
771
772