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