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