1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ////////////////////////////////////////////////////////////////////////////
20 // Example: fill pedestal with Gaussian noise //
22 // AliTRDCalibPadStatus ped; //
23 // ped.TestEvent(numberofevent); //
25 // // Method without histo //
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) //
39 // TCanvas* c1 = new TCanvas; //
40 // c1->Divide(2,2); //
42 // histo2dm->Draw("colz"); //
44 // histo1dm->Draw(); //
46 // histo2ds->Draw("colz"); //
48 // histo1ds->Draw(); //
50 // // Method with histo //
51 // ped.AnalyseHisto(); //
53 // // Take the histo //
54 // TH1F *histo = ped.GetHisto(31); //
55 // histo->SetEntries(1); //
58 ////////////////////////////////////////////////////////////////////////////
61 #include <TObjArray.h>
66 #include <TDirectory.h>
70 #include <AliMathBase.h>
71 #include "AliRawReader.h"
72 #include "AliRawReaderRoot.h"
73 #include "AliRawReaderDate.h"
77 #include "AliTRDCalibPadStatus.h"
78 #include "AliTRDrawStreamBase.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"
91 ClassImp(AliTRDCalibPadStatus) /*FOLD00*/
93 //_____________________________________________________________________
94 AliTRDCalibPadStatus::AliTRDCalibPadStatus() : /*FOLD00*/
100 fNumberOfTimeBins(0),
101 fCalRocArrayMean(540),
102 fCalRocArrayRMS(540),
103 fCalRocArrayMeand(540),
104 fCalRocArrayRMSd(540),
108 // default constructor
111 fGeo = new AliTRDgeometry();
115 //_____________________________________________________________________
116 AliTRDCalibPadStatus::AliTRDCalibPadStatus(const AliTRDCalibPadStatus &ped) : /*FOLD00*/
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),
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);
139 if ( calRocMean != 0x0 ) fCalRocArrayMean.AddAt(new AliTRDCalROC(*calRocMean), idet);
140 if ( calRocRMS != 0x0 ) fCalRocArrayRMS.AddAt(new AliTRDCalROC(*calRocRMS), idet);
142 if ( calRocMeand != 0x0 ) fCalRocArrayMeand.AddAt(new AliTRDCalROC(*calRocMeand), idet);
143 if ( calRocRMSd != 0x0 ) fCalRocArrayRMSd.AddAt(new AliTRDCalROC(*calRocRMSd), idet);
146 TH2F *hNew = new TH2F(*hped);
147 hNew->SetDirectory(0);
148 fHistoArray.AddAt(hNew,idet);
155 fGeo = new AliTRDgeometry();
158 //_____________________________________________________________________
159 AliTRDCalibPadStatus& AliTRDCalibPadStatus::operator = (const AliTRDCalibPadStatus &source)
162 // assignment operator
164 if (&source == this) return *this;
165 new (this) AliTRDCalibPadStatus(source);
169 //_____________________________________________________________________
170 AliTRDCalibPadStatus::~AliTRDCalibPadStatus() /*FOLD00*/
175 fCalRocArrayMean.Delete();
176 fCalRocArrayRMS.Delete();
177 fCalRocArrayMeand.Delete();
178 fCalRocArrayRMSd.Delete();
179 fHistoArray.Delete();
184 //_____________________________________________________________________
185 void AliTRDCalibPadStatus::Destroy()
190 fCalRocArrayMean.Delete();
191 fCalRocArrayRMS.Delete();
192 fCalRocArrayMeand.Delete();
193 fCalRocArrayRMSd.Delete();
194 fHistoArray.Delete();
196 //_____________________________________________________________________
197 Int_t AliTRDCalibPadStatus::UpdateHisto(const Int_t icdet, /*FOLD00*/
206 // Signal filling methode
208 Int_t nbchannel = icRow+icCol*crowMax;
210 // now the case of double read channel
212 nbchannel = (((ccold-1)*8+ icMcm)*crowMax+icRow)+144*crowMax;
213 //printf("nbchannel %d, ccold %d, icMcm %d, crowMax %d, icRow %d\n",nbchannel,ccold,icMcm,crowMax,icRow);
216 // fast filling methode.
217 // Attention: the entry counter of the histogram is not increased
218 // this means that e.g. the colz draw option gives an empty plot
220 if ( !(((Int_t)csignal>=fAdcMax ) || ((Int_t)csignal<fAdcMin)) )
221 bin = (nbchannel+1)*(fAdcMax-fAdcMin+2)+((Int_t)csignal-fAdcMin+1);
223 //GetHisto(icdet,kTRUE)->Fill(csignal,nbchannel);
225 GetHisto(icdet,kTRUE)->GetArray()[bin]++;
229 //_____________________________________________________________________
230 Int_t AliTRDCalibPadStatus::ProcessEvent(AliTRDrawStreamBase *rawStream, Bool_t nocheck)
233 // Event Processing loop - AliTRDRawStreamCosmic
234 // 0 time bin problem or zero suppression
240 // Raw version number:
241 // [3,31] non zero suppressed
242 // 2,4 and [32,63] zero suppressed
247 rawStream->SetSharedPadReadout(kTRUE);
251 // Check the raw version and if all have the same number of timebins.
253 while (rawStream->Next()) {
255 Int_t rawversion = rawStream->GetRawVersion(); // current raw version
256 //printf("Raw version is %d\n",rawversion);
258 // Could eventually change, have to check with time
259 if((rawversion < 3) || (rawversion > 31)) {
260 AliInfo(Form("this is not no-zero-suppressed data, the version is %d",rawversion));
263 Int_t idetector = rawStream->GetDet(); // current detector
264 Int_t iRow = rawStream->GetRow(); // current row
265 Int_t iRowMax = rawStream->GetMaxRow(); // current rowmax
266 Int_t iCol = rawStream->GetCol(); // current col
267 Int_t iADC = 21-rawStream->GetADC(); // current ADC
269 // It goes in the opposite direction
271 if(iADC == 1) col = 1;
273 col = TMath::Max(0,(Int_t)(iADC-19));
276 Int_t mcm = (Int_t)(iCol/18); // current group of 18 col pads
277 if(col > 1) mcm -= 1;
278 if(col ==1) mcm += 1;
281 //Bool_t shared = rawStream->IsCurrentPadShared();
282 //printf("ADC %d, iCol %d, col %d, mcm %d, shared %d\n",iADC,iCol,col,mcm,(Int_t)shared);
285 Int_t *signal = rawStream->GetSignals(); // current ADC signal
286 Int_t nbtimebin = rawStream->GetNumberOfTimeBins(); // number of time bins read from data
288 if((fDetector != -1) && (nbtimebin != fNumberOfTimeBins)) {
289 AliInfo(Form("the number of time bins is %d, is different from the previous one %d",nbtimebin,fNumberOfTimeBins));
292 fNumberOfTimeBins = nbtimebin;
293 fDetector = idetector;
295 for(Int_t k = 0; k < fNumberOfTimeBins; k++){
296 if(signal[k]>0) UpdateHisto(idetector,iRow,iCol,signal[k],iRowMax,col,mcm);
304 while (rawStream->Next()) {
306 Int_t idetector = rawStream->GetDet(); // current detector
307 Int_t iRow = rawStream->GetRow(); // current row
308 Int_t iRowMax = rawStream->GetMaxRow(); // current rowmax
309 Int_t iCol = rawStream->GetCol(); // current col
310 Int_t iADC = 21-rawStream->GetADC(); // current ADC
312 // It goes in the opposite direction
314 if(iADC == 1) col = 1;
316 col = TMath::Max(0,(Int_t)(iADC-19));
319 Int_t mcm = (Int_t)(iCol/18); // current group of 18 col pads
320 if(col > 1) mcm -= 1;
321 if(col ==1) mcm += 1;
324 Int_t *signal = rawStream->GetSignals(); // current ADC signal
325 Int_t nbtimebin = rawStream->GetNumberOfTimeBins(); // number of time bins read from data
328 //printf("det %d, row %d, signal[0] %d, signal[1] %d, signal [2] %d\n", idetector, iRow, signal[0], signal[1], signal[2]);
330 for(Int_t k = 0; k < nbtimebin; k++){
332 UpdateHisto(idetector,iRow,iCol,signal[k],iRowMax,col,mcm);
333 //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);
343 //_____________________________________________________________________
344 Int_t AliTRDCalibPadStatus::ProcessEvent(AliRawReader *rawReader, Bool_t nocheck)
347 // Event processing loop - AliRawReader
352 rawReader->Select("TRD");
354 AliTRDrawStreamBase *pstream = AliTRDrawStreamBase::GetRawStream(rawReader);
356 result = ProcessEvent(pstream, nocheck);
363 //_________________________________________________________________________
364 Int_t AliTRDCalibPadStatus::ProcessEvent(
366 eventHeaderStruct *event,
369 eventHeaderStruct* /*event*/,
376 // process date event
379 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
380 Bool_t result=ProcessEvent(rawReader, nocheck);
384 Fatal("AliTRDCalibPadStatus", "this class was compiled without DATE");
390 //_____________________________________________________________________
391 Bool_t AliTRDCalibPadStatus::TestEventHisto(Int_t nevent, Int_t sm, Int_t ch) /*FOLD00*/
395 // fill one oroc and one iroc with random gaus
400 for (Int_t ism=sm; ism<sm+1; ism++){
401 for (Int_t ich=ch; ich < ch+1; ich++){
402 for (Int_t ipl=0; ipl < 6; ipl++){
403 for(Int_t irow = 0; irow < fGeo->GetRowMax(ipl,ich,ism); irow++){
404 for(Int_t icol = 0; icol < fGeo->GetColMax(ipl); icol++){
405 for (Int_t iTimeBin=0; iTimeBin<(30*nevent); iTimeBin++){
406 Int_t signal=TMath::Nint(gRandom->Gaus(10,1.5));
407 if ( signal>0 )UpdateHisto((ipl+ich*6+ism*6*5),irow,icol,signal,fGeo->GetRowMax(ipl,ich,ism),0,0);
417 //_____________________________________________________________________
418 TH2F* AliTRDCalibPadStatus::GetHisto(Int_t det, TObjArray *arr, /*FOLD00*/
419 Int_t nbinsY, Float_t ymin, Float_t ymax,
420 const Char_t *type, Bool_t force)
423 // return pointer to histogram
424 // if force is true create a new histogram if it doesn't exist allready
426 if ( !force || arr->UncheckedAt(det) )
427 return (TH2F*)arr->UncheckedAt(det);
429 // if we are forced and histogram doesn't yes exist create it
430 Char_t name[255], title[255];
432 sprintf(name,"hCalib%s%.3d",type,det);
433 sprintf(title,"%s calibration histogram detector %.2d;ADC channel;Channel (pad)",type,det);
436 Int_t nbchannels = fGeo->GetRowMax(GetLayer(det),GetStack(det),GetSector(det))*fGeo->GetColMax(GetLayer(det));
438 // we will add 3*8*rowMax channels at the end for the double counted
439 nbchannels += 3*8*(fGeo->GetRowMax(GetLayer(det),GetStack(det),GetSector(det)));
442 // new histogram with calib information. One value for each pad!
443 TH2F* hist = new TH2F(name,title,
445 nbchannels,0,nbchannels
447 hist->SetDirectory(0);
448 arr->AddAt(hist,det);
452 //_____________________________________________________________________
453 TH2F* AliTRDCalibPadStatus::GetHisto(Int_t det, Bool_t force) /*FOLD00*/
456 // return pointer to histogram
457 // if force is true create a new histogram if it doesn't exist allready
459 TObjArray *arr = &fHistoArray;
460 return GetHisto(det, arr, fAdcMax-fAdcMin, fAdcMin-0.5, fAdcMax-0.5, "Pedestal", force);
463 //_____________________________________________________________________
464 AliTRDCalROC* AliTRDCalibPadStatus::GetCalRoc(Int_t det, TObjArray* arr, Bool_t force) /*FOLD00*/
467 // return pointer to ROC Calibration
468 // if force is true create a new AliTRDCalROC if it doesn't exist allready
470 if ( !force || arr->UncheckedAt(det) )
471 return (AliTRDCalROC*)arr->UncheckedAt(det);
473 // if we are forced and histogram doesn't yes exist create it
475 // new AliTRDCalROC. One value for each pad!
476 AliTRDCalROC *croc = new AliTRDCalROC(GetLayer(det),GetStack(det));
477 arr->AddAt(croc,det);
480 //_____________________________________________________________________
481 AliTRDCalROC* AliTRDCalibPadStatus::GetCalRocMean(Int_t det, Bool_t force) /*FOLD00*/
484 // return pointer to Carge ROC Calibration
485 // if force is true create a new histogram if it doesn't exist allready
487 TObjArray *arr = &fCalRocArrayMean;
488 return GetCalRoc(det, arr, force);
491 //_____________________________________________________________________
492 AliTRDCalROC* AliTRDCalibPadStatus::GetCalRocRMS(Int_t det, Bool_t force) /*FOLD00*/
495 // return pointer to Carge ROC Calibration
496 // if force is true create a new histogram if it doesn't exist allready
498 TObjArray *arr = &fCalRocArrayRMS;
499 return GetCalRoc(det, arr, force);
501 //_____________________________________________________________________
502 AliTRDCalROC* AliTRDCalibPadStatus::GetCalRocMeand(Int_t det, Bool_t force) /*FOLD00*/
505 // return pointer to Carge ROC Calibration
506 // if force is true create a new histogram if it doesn't exist allready
508 TObjArray *arr = &fCalRocArrayMeand;
509 return GetCalRoc(det, arr, force);
512 //_____________________________________________________________________
513 AliTRDCalROC* AliTRDCalibPadStatus::GetCalRocRMSd(Int_t det, Bool_t force) /*FOLD00*/
516 // return pointer to Carge ROC Calibration
517 // if force is true create a new histogram if it doesn't exist allready
519 TObjArray *arr = &fCalRocArrayRMSd;
520 return GetCalRoc(det, arr, force);
523 //_____________________________________________________________________
524 void AliTRDCalibPadStatus::AnalyseHisto() /*FOLD00*/
527 // Calculate calibration constants
530 Int_t nbinsAdc = fAdcMax-fAdcMin;
538 for (Int_t idet=0; idet<540; idet++){
539 TH2F *hP = GetHisto(idet);
544 //printf("Entries for %d\n",idet);
546 AliTRDCalROC *rocMean = GetCalRocMean(idet,kTRUE);
547 AliTRDCalROC *rocRMS = GetCalRocRMS(idet,kTRUE);
549 arrayHP = hP->GetArray();
550 Int_t nChannels = rocMean->GetNchannels();
552 for (Int_t iChannel=0; iChannel<nChannels; iChannel++){
553 Int_t offset = (nbinsAdc+2)*(iChannel+1)+1;
554 Double_t ret = AliMathBase::FitGaus(arrayHP+offset,nbinsAdc,fAdcMin-0.5,fAdcMax-0.5,¶m,&dummy);
555 // if the fitting failed set noise and pedestal to 0
556 if ((ret==-4) || (ret==-1) || (ret==-2)) {
560 if((param[1]/10.0) > 65534.0) param[1] = 0.0;
561 if((param[2]/10.0) > 65534.0) param[2] = 0.0;
562 rocMean->SetValue(iChannel,param[1]/10.0);
563 rocRMS->SetValue(iChannel,param[2]/10.0);
566 // here we analyse doubled read channels
568 AliTRDCalROC *rocMeand = GetCalRocMeand(idet,kTRUE);
569 AliTRDCalROC *rocRMSd = GetCalRocRMSd(idet,kTRUE);
571 Int_t nrows = rocMeand->GetNrows();
572 Int_t shift = 144*nrows;
573 Int_t total = shift+3*8*nrows;
575 for (Int_t iChannel=shift; iChannel<total; iChannel++){
576 Int_t offset = (nbinsAdc+2)*(iChannel+1)+1;
577 Double_t ret = AliMathBase::FitGaus(arrayHP+offset,nbinsAdc,fAdcMin-0.5,fAdcMax-0.5,¶m,&dummy);
578 // if the fitting failed set noise and pedestal to 0
579 if ((ret==-4) || (ret==-1) || (ret==-2)) {
583 if((param[1]/10.0) > 65534.0) param[1] = 0.0;
584 if((param[2]/10.0) > 65534.0) param[2] = 0.0;
586 // here we have to recalculate backward
587 Int_t nb = iChannel-shift;
588 Int_t row = nb%nrows;
589 Int_t j = (Int_t)(nb/nrows);
591 Int_t icol = (Int_t)(j/8);
593 Int_t finalcol = 18*imcm;
594 if(icol > 0) icol += 17;
598 Int_t channel = row+finalcol*nrows;
600 //printf("iChannel %d, nrows %d, finalcol %d, row %d, channel %d\n",iChannel,nrows,finalcol,row,channel);
601 if((finalcol < 0) || (finalcol >= 144)) continue;
603 rocMeand->SetValue(channel,param[1]/10.0);
604 rocRMSd->SetValue(channel,param[2]/10.0);
611 //_______________________________________________________________________________________
612 AliTRDCalPadStatus* AliTRDCalibPadStatus::CreateCalPadStatus()
615 // Create Pad Status out of Mean and RMS values
616 // The chamber without data are masked, this is the corrected in the preprocessor
619 AliTRDCalPadStatus* obj = new AliTRDCalPadStatus("padstatus", "padstatus");
621 for (Int_t idet=0; idet<540; ++idet)
623 AliTRDCalSingleChamberStatus *calROC = obj->GetCalROC(idet);
626 if ( !GetCalRocMean(idet)) {
627 for(Int_t k = 0; k < calROC->GetNchannels(); k++){
628 calROC->SetStatus(k,AliTRDCalPadStatus::kNotConnected);
635 AliTRDCalROC *calRocMean = new AliTRDCalROC(*( (AliTRDCalROC *) GetCalRocMean(idet)));
636 AliTRDCalROC *calRocRMS = new AliTRDCalROC(*( (AliTRDCalROC *) GetCalRocRMS(idet)));
638 //Take the stuff second chance
639 AliTRDCalROC *calRocMeand = new AliTRDCalROC(*( (AliTRDCalROC *) GetCalRocMeand(idet)));
640 AliTRDCalROC *calRocRMSd = new AliTRDCalROC(*( (AliTRDCalROC *) GetCalRocRMSd(idet)));
643 calRocRMSd->Unfold();
647 Int_t row = calROC->GetNrows();
648 Int_t col = calROC->GetNcols();
650 Double_t rmsmean = calRocMean->GetRMS()*10.0;
651 Double_t meanmean = calRocMean->GetMean()*10.0;
652 Double_t meansquares = calRocRMS->GetMean();
655 for(Int_t irow = 0; irow < row; irow++){
658 Float_t meanprevious = 0.0;
659 Float_t rmsprevious = 0.0;
663 for(Int_t icol = 0; icol < col; icol++){
665 mean = calRocMean->GetValue(icol,irow)*10.0;
666 rms = calRocRMS->GetValue(icol,irow);
669 meanprevious = calRocMean->GetValue((icol -1),irow)*10.0;
670 rmsprevious = calRocRMS->GetValue((icol - 1),irow);
674 // masked if two noisy
675 if((rms <= 0.0001) || (TMath::Abs(mean-meanmean)>(5*rmsmean)) || (TMath::Abs(rms)>(5.0*TMath::Abs(meansquares)))) {
678 // look at second chance
679 Float_t meand = calRocMeand->GetValue(icol,irow)*10.0;
680 Float_t rmsd = calRocRMSd->GetValue(icol,irow);
682 if((rmsd <= 0.0001) || (TMath::Abs(meand-meanmean)>(5*rmsmean)) || (TMath::Abs(rmsd)>(5.0*TMath::Abs(meansquares)))) {
683 if((rmsd <= 0.0001) && (rms <= 0.0001)) {
684 calROC->SetStatus(icol,irow,AliTRDCalPadStatus::kNotConnected);
687 calROC->SetStatus(icol, irow, AliTRDCalPadStatus::kMasked);
691 calROC->SetStatus(icol, irow, AliTRDCalPadStatus::kReadSecond);
696 // bridge if previous pad found something
698 if((meanprevious == mean) && (rmsprevious == rms) && (mean > 0.0001)) {
699 //printf("mean previous %f, mean %f, rms %f, rmsprevious %f, col %d\n",meanprevious,mean,rms,rmsprevious,icol);
700 calROC->SetStatus(icol -1 ,irow, AliTRDCalPadStatus::kPadBridgedRight);
701 calROC->SetStatus(icol ,irow, AliTRDCalPadStatus::kPadBridgedLeft);
719 //_______________________________________________________________________________________
720 AliTRDCalPad* AliTRDCalibPadStatus::CreateCalPad()
723 // Create Pad Noise out of RMS values
726 AliTRDCalPad* obj = new AliTRDCalPad("PadNoise", "PadNoise");
729 for (Int_t det=0; det<AliTRDgeometry::kNdet; ++det) {
731 AliTRDCalROC *calROC22 = obj->GetCalROC(det);
733 AliTRDCalROC *calRocRMS = ((AliTRDCalROC *)GetCalRocRMS(det,kTRUE));
735 for(Int_t k = 0; k < calROC22->GetNchannels(); k++){
736 calROC22->SetValue(k,calRocRMS->GetValue(k));
745 //_______________________________________________________________________________________
746 AliTRDCalDet* AliTRDCalibPadStatus::CreateCalDet() const
749 // Create Det Noise correction factor
752 AliTRDCalDet* obj = new AliTRDCalDet("DetNoise", "DetNoise (correction factor)");
754 for(Int_t l = 0; l < 540; l++){
755 obj->SetValue(l,10.0);
762 //_____________________________________________________________________
763 void AliTRDCalibPadStatus::DumpToFile(const Char_t *filename, const Char_t *dir, Bool_t append) /*FOLD00*/
766 // Write class to file
777 TDirectory *backup = gDirectory;
778 TFile f(filename,option.Data());
780 if ( !sDir.IsNull() ){
781 f.mkdir(sDir.Data());
787 if ( backup ) backup->cd();
790 //_____________________________________________________________________
791 void AliTRDCalibPadStatus::SetCalRocMean(AliTRDCalROC *mean, Int_t det) /*FOLD00*/
794 // Put the AliTRDCalROC in the array fCalRocArrayMean
798 AliTRDCalROC *rocMean = GetCalRocMean(det,kTRUE);
800 Int_t nChannels = rocMean->GetNchannels();
802 for (Int_t iChannel=0; iChannel<nChannels; iChannel++){
804 rocMean->SetValue(iChannel,mean->GetValue(iChannel));
810 //_____________________________________________________________________
811 void AliTRDCalibPadStatus::SetCalRocRMS(AliTRDCalROC *rms, Int_t det) /*FOLD00*/
814 // Put the AliTRDCalROC in the array fCalRocArrayRMS
818 AliTRDCalROC *rocRms = GetCalRocRMS(det,kTRUE);
820 Int_t nChannels = rocRms->GetNchannels();
822 for (Int_t iChannel=0; iChannel<nChannels; iChannel++){
824 rocRms->SetValue(iChannel,rms->GetValue(iChannel));
829 //_____________________________________________________________________
830 void AliTRDCalibPadStatus::SetCalRocMeand(AliTRDCalROC *mean, Int_t det) /*FOLD00*/
833 // Put the AliTRDCalROC in the array fCalRocArrayMean
837 AliTRDCalROC *rocMean = GetCalRocMeand(det,kTRUE);
839 Int_t nChannels = rocMean->GetNchannels();
841 for (Int_t iChannel=0; iChannel<nChannels; iChannel++){
843 rocMean->SetValue(iChannel,mean->GetValue(iChannel));
849 //_____________________________________________________________________
850 void AliTRDCalibPadStatus::SetCalRocRMSd(AliTRDCalROC *rms, Int_t det) /*FOLD00*/
853 // Put the AliTRDCalROC in the array fCalRocArrayRMS
857 AliTRDCalROC *rocRms = GetCalRocRMSd(det,kTRUE);
859 Int_t nChannels = rocRms->GetNchannels();
861 for (Int_t iChannel=0; iChannel<nChannels; iChannel++){
863 rocRms->SetValue(iChannel,rms->GetValue(iChannel));
868 //_____________________________________________________________________________
869 Int_t AliTRDCalibPadStatus::GetLayer(Int_t d) const
872 // Reconstruct the layer number from the detector number
875 return ((Int_t) (d % 6));
879 //_____________________________________________________________________________
880 Int_t AliTRDCalibPadStatus::GetStack(Int_t d) const
883 // Reconstruct the chamber number from the detector number
886 return ((Int_t) (d % 30) / 6);
890 //_____________________________________________________________________________
891 Int_t AliTRDCalibPadStatus::GetSector(Int_t d) const
894 // Reconstruct the sector number from the detector number
897 return ((Int_t) (d / 30));