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