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