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); //
59 // R. Bailhache (R.Bailhache@gsi.de, rbailhache@ikf.uni-frankfurt.de)
60 // J. Book (jbook@ikf.uni-frankfurt.de)
62 ////////////////////////////////////////////////////////////////////////////
66 #include <TObjArray.h>
72 //#include <TRandom.h>
73 #include <TDirectory.h>
77 #include <AliMathBase.h>
78 #include "AliRawReader.h"
79 #include "AliRawReaderRoot.h"
80 #include "AliRawReaderDate.h"
84 #include "AliTRDCalibPadStatus.h"
85 #include "AliTRDgeometry.h"
86 #include "AliTRDCommonParam.h"
87 #include "./Cal/AliTRDCalROC.h"
88 #include "./Cal/AliTRDCalPadStatus.h"
89 #include "./Cal/AliTRDCalDet.h"
90 #include "./Cal/AliTRDCalPad.h"
91 #include "./Cal/AliTRDCalSingleChamberStatus.h"
93 #include "AliTRDdigitsManager.h"
94 #include "AliTRDdigitsParam.h"
95 #include "AliTRDSignalIndex.h"
96 #include "AliTRDarraySignal.h"
97 #include "AliTRDarrayADC.h"
98 #include "AliTRDfeeParam.h"
100 #include "AliTRDrawStream.h"
106 ClassImp(AliTRDCalibPadStatus) /*FOLD00*/
108 //_____________________________________________________________________
109 AliTRDCalibPadStatus::AliTRDCalibPadStatus() : /*FOLD00*/
115 fNumberOfTimeBins(0),
116 fCalRocArrayMean(540),
117 fCalRocArrayRMS(540),
118 fCalRocArrayMeand(540),
119 fCalRocArrayRMSd(540),
123 // default constructor
126 fGeo = new AliTRDgeometry();
130 //_____________________________________________________________________
131 AliTRDCalibPadStatus::AliTRDCalibPadStatus(const AliTRDCalibPadStatus &ped) : /*FOLD00*/
134 fAdcMin(ped.GetAdcMin()),
135 fAdcMax(ped.GetAdcMax()),
136 fDetector(ped.fDetector),
137 fNumberOfTimeBins(ped.fNumberOfTimeBins),
138 fCalRocArrayMean(540),
139 fCalRocArrayRMS(540),
140 fCalRocArrayMeand(540),
141 fCalRocArrayRMSd(540),
147 for (Int_t idet = 0; idet < 540; idet++){
148 const AliTRDCalROC *calRocMean = (AliTRDCalROC*)ped.fCalRocArrayMean.UncheckedAt(idet);
149 const AliTRDCalROC *calRocRMS = (AliTRDCalROC*)ped.fCalRocArrayRMS.UncheckedAt(idet);
150 const AliTRDCalROC *calRocMeand = (AliTRDCalROC*)ped.fCalRocArrayMeand.UncheckedAt(idet);
151 const AliTRDCalROC *calRocRMSd = (AliTRDCalROC*)ped.fCalRocArrayRMSd.UncheckedAt(idet);
152 const TH2F *hped = (TH2F*)ped.fHistoArray.UncheckedAt(idet);
154 if ( calRocMean != 0x0 ) fCalRocArrayMean.AddAt(new AliTRDCalROC(*calRocMean), idet);
155 if ( calRocRMS != 0x0 ) fCalRocArrayRMS.AddAt(new AliTRDCalROC(*calRocRMS), idet);
157 if ( calRocMeand != 0x0 ) fCalRocArrayMeand.AddAt(new AliTRDCalROC(*calRocMeand), idet);
158 if ( calRocRMSd != 0x0 ) fCalRocArrayRMSd.AddAt(new AliTRDCalROC(*calRocRMSd), idet);
161 TH2F *hNew = new TH2F(*hped);
162 hNew->SetDirectory(0);
163 fHistoArray.AddAt(hNew,idet);
170 fGeo = new AliTRDgeometry();
173 //_____________________________________________________________________
174 AliTRDCalibPadStatus& AliTRDCalibPadStatus::operator = (const AliTRDCalibPadStatus &source)
177 // assignment operator
179 if (&source == this) return *this;
180 new (this) AliTRDCalibPadStatus(source);
184 //_____________________________________________________________________
185 AliTRDCalibPadStatus::~AliTRDCalibPadStatus() /*FOLD00*/
190 fCalRocArrayMean.Delete();
191 fCalRocArrayRMS.Delete();
192 fCalRocArrayMeand.Delete();
193 fCalRocArrayRMSd.Delete();
194 fHistoArray.Delete();
199 //_____________________________________________________________________
200 void AliTRDCalibPadStatus::Destroy()
205 fCalRocArrayMean.Delete();
206 fCalRocArrayRMS.Delete();
207 fCalRocArrayMeand.Delete();
208 fCalRocArrayRMSd.Delete();
209 fHistoArray.Delete();
211 //_____________________________________________________________________
212 Int_t AliTRDCalibPadStatus::UpdateHisto(const Int_t icdet, /*FOLD00*/
223 // Signal filling methode
225 Int_t nbchannel = icRow+icCol*crowMax;
229 // now the case of double read channel
231 nbchannel = (((ccold-1)*8+ (mCm+rOb*4))*crowMax+icRow)+144*crowMax;
232 //printf("nbchannel %d, ccold %d, icMcm %d, crowMax %d, icRow %d\n",nbchannel,ccold,icMcm,crowMax,icRow);
235 // fast filling methode.
236 // Attention: the entry counter of the histogram is not increased
237 // this means that e.g. the colz draw option gives an empty plot
241 if ( !(((Int_t)csignal>=fAdcMax ) || ((Int_t)csignal<fAdcMin)) )
242 bin = (nbchannel+1)*(fAdcMax-fAdcMin+2)+((Int_t)csignal-fAdcMin+1);
244 //GetHisto(icdet,kTRUE)->Fill(csignal,nbchannel);
246 GetHisto(icdet,kTRUE)->GetArray()[bin]++;
250 //_____________________________________________________________________
251 Int_t AliTRDCalibPadStatus::ProcessEvent(AliRawReader *rawReader)
254 // RawReader = AliTRDrawStream (Jochen Klein)
256 // Event Processing loop - AliTRDRawStreamCosmic
257 // 0 time bin problem or zero suppression
260 // Raw version number:
261 // [3,31] non zero suppressed
262 // 2,4 and [32,63] zero suppressed
269 AliTRDdigitsManager *digitsManager = new AliTRDdigitsManager(kTRUE);
270 digitsManager->CreateArrays();
272 AliTRDrawStream *rawStream = new AliTRDrawStream(rawReader);
273 rawStream->SetDigitsManager(digitsManager);
274 //rawStream->SetNoErrorWarning();
275 //rawStream->SetSharedPadReadout(kTRUE);
277 AliTRDfeeParam *feeParam = AliTRDfeeParam::Instance();
280 while ((det = rawStream->NextChamber(digitsManager, NULL, NULL)) >= 0) { //idetector
281 if (digitsManager->GetIndexes(det)->HasEntry()) {//QA
282 // printf("there is ADC data on this chamber!\n");
284 AliTRDarrayADC *digits = (AliTRDarrayADC *) digitsManager->GetDigits(det); //mod
285 if (digits->HasData()) { //array
287 AliTRDSignalIndex *indexes = digitsManager->GetIndexes(det);
288 if (indexes->IsAllocated() == kFALSE) {
289 AliError("Indexes do not exist!");
294 indexes->ResetCounters();
296 while (indexes->NextRCIndex(iRow, iCol)) { //column,row
298 AliTRDdigitsParam *digitParam = (AliTRDdigitsParam *)digitsManager->GetDigitsParam();
300 Int_t mcm = 0; // MCM from AliTRDfeeParam
301 Int_t rob = 0; // ROB from AliTRDfeeParam
302 Int_t extCol = 0; // extended column from AliTRDfeeParam
303 mcm = feeParam->GetMCMfromPad(iRow,iCol);
304 rob = feeParam->GetROBfromPad(iRow,iCol);
306 Int_t idetector = det; // current detector
307 Int_t iRowMax = 16; // current rowmax
308 if(GetStack(det) == 2) iRowMax = 12;
310 Int_t adc = 20 - (iCol%18) -1; // current adc
311 Int_t col = 0; // col!=0 ->Shared Pad
312 extCol = feeParam->GetExtendedPadColFromADC(rob,mcm,adc);
313 //printf(" iCol %d iRow %d iRowMax %d rob %d mcm %d adc %d extCol %d\n",iCol,iRow,iRowMax,rob,mcm,adc,extCol);
315 // Signal for regular pads
316 Int_t nbtimebin = digitParam->GetNTimeBins(idetector); // number of time bins read from data
317 for(Int_t k = 0; k < nbtimebin; k++){
319 signal = digits->GetData(iRow,iCol,k);
322 UpdateHisto(idetector,iRow,iCol,signal,iRowMax,col,mcm,rob);
328 if((adc==3-1 || adc==20-1 || adc==19-1) && (iCol > 1 && iCol <142) ) { //SHARED PADS
332 adc = 20; //shared Pad adc
333 mcm = feeParam->GetMCMfromSharedPad(iRow,iCol); //shared Pad mcm
337 adc = 1; //shared Pad adc
338 mcm = feeParam->GetMCMfromSharedPad(iRow,iCol); //shared Pad mcm
342 adc = 0; //shared Pad adc
343 mcm = feeParam->GetMCMfromSharedPad(iRow,iCol); //shared Pad mcm
347 rob = feeParam->GetROBfromSharedPad(iRow,iCol); //shared Pad rob
350 extCol = feeParam->GetExtendedPadColFromADC(rob,mcm,adc); //extended pad col via the shared pad rob,mcm and adc
352 //printf("SHARED PAD --- iCol %d iRow %d rob %d mcm %d adc %d extCol %d col %d\n",iCol,iRow,rob,mcm,adc,extCol,col);
353 for(Int_t k = 0; k < nbtimebin; k++){
355 signal = digits->GetDataByAdcCol(iRow,extCol,k);
358 UpdateHisto(idetector,iRow,iCol,signal,iRowMax,col,mcm,rob);
369 digitsManager->ClearArrays(det);
371 delete digitsManager;
377 //_____________________________________________________________________
378 TH2F* AliTRDCalibPadStatus::GetHisto(Int_t det, TObjArray *arr, /*FOLD00*/
379 Int_t nbinsY, Float_t ymin, Float_t ymax,
380 const Char_t *type, Bool_t force)
383 // return pointer to histogram
384 // if force is true create a new histogram if it doesn't exist allready
386 if ( !force || arr->UncheckedAt(det) )
387 return (TH2F*)arr->UncheckedAt(det);
389 // if we are forced and histogram doesn't yes exist create it
390 Char_t name[255], title[255];
392 snprintf(name,255,"hCalib%s%.3d",type,det);
393 snprintf(title,255,"%s calibration histogram detector %.2d;ADC channel;Channel (pad)",type,det);
396 Int_t nbchannels = fGeo->GetRowMax(GetLayer(det),GetStack(det),GetSector(det))*fGeo->GetColMax(GetLayer(det));
398 // we will add 3*8*rowMax channels at the end for the double counted
399 nbchannels += 3*8*(fGeo->GetRowMax(GetLayer(det),GetStack(det),GetSector(det)));
402 // new histogram with calib information. One value for each pad!
403 TH2F* hist = new TH2F(name,title,
405 nbchannels,0,nbchannels
407 hist->SetDirectory(0);
408 arr->AddAt(hist,det);
412 //_____________________________________________________________________
413 TH2F* AliTRDCalibPadStatus::GetHisto(Int_t det, Bool_t force) /*FOLD00*/
416 // return pointer to histogram
417 // if force is true create a new histogram if it doesn't exist allready
419 TObjArray *arr = &fHistoArray;
420 return GetHisto(det, arr, fAdcMax-fAdcMin, fAdcMin-0.5, fAdcMax-0.5, "Pedestal", force);
423 //_____________________________________________________________________
424 AliTRDCalROC* AliTRDCalibPadStatus::GetCalRoc(Int_t det, TObjArray* arr, Bool_t force) /*FOLD00*/
427 // return pointer to ROC Calibration
428 // if force is true create a new AliTRDCalROC if it doesn't exist allready
430 if ( !force || arr->UncheckedAt(det) )
431 return (AliTRDCalROC*)arr->UncheckedAt(det);
433 // if we are forced and histogram doesn't yes exist create it
435 // new AliTRDCalROC. One value for each pad!
436 AliTRDCalROC *croc = new AliTRDCalROC(GetLayer(det),GetStack(det));
437 arr->AddAt(croc,det);
440 //_____________________________________________________________________
441 AliTRDCalROC* AliTRDCalibPadStatus::GetCalRocMean(Int_t det, Bool_t force) /*FOLD00*/
444 // return pointer to Carge ROC Calibration
445 // if force is true create a new histogram if it doesn't exist allready
447 TObjArray *arr = &fCalRocArrayMean;
448 return GetCalRoc(det, arr, force);
451 //_____________________________________________________________________
452 AliTRDCalROC* AliTRDCalibPadStatus::GetCalRocRMS(Int_t det, Bool_t force) /*FOLD00*/
455 // return pointer to Carge ROC Calibration
456 // if force is true create a new histogram if it doesn't exist allready
458 TObjArray *arr = &fCalRocArrayRMS;
459 return GetCalRoc(det, arr, force);
461 //_____________________________________________________________________
462 AliTRDCalROC* AliTRDCalibPadStatus::GetCalRocMeand(Int_t det, Bool_t force) /*FOLD00*/
465 // return pointer to Carge ROC Calibration
466 // if force is true create a new histogram if it doesn't exist allready
468 TObjArray *arr = &fCalRocArrayMeand;
469 return GetCalRoc(det, arr, force);
472 //_____________________________________________________________________
473 AliTRDCalROC* AliTRDCalibPadStatus::GetCalRocRMSd(Int_t det, Bool_t force) /*FOLD00*/
476 // return pointer to Carge ROC Calibration
477 // if force is true create a new histogram if it doesn't exist allready
479 TObjArray *arr = &fCalRocArrayRMSd;
480 return GetCalRoc(det, arr, force);
483 //_____________________________________________________________________
484 void AliTRDCalibPadStatus::AnalyseHisto() /*FOLD00*/
487 // Calculate calibration constants
490 Int_t nbinsAdc = fAdcMax-fAdcMin;
498 for (Int_t idet=0; idet<540; idet++){
499 TH2F *hP = GetHisto(idet);
504 //printf("Entries for %d\n",idet);
506 AliTRDCalROC *rocMean = GetCalRocMean(idet,kTRUE);
507 AliTRDCalROC *rocRMS = GetCalRocRMS(idet,kTRUE);
509 arrayHP = hP->GetArray();
510 Int_t nChannels = rocMean->GetNchannels();
512 for (Int_t iChannel=0; iChannel<nChannels; iChannel++){
513 Int_t offset = (nbinsAdc+2)*(iChannel+1)+1;
514 Double_t ret = AliMathBase::FitGaus(arrayHP+offset,nbinsAdc,fAdcMin-0.5,fAdcMax-0.5,¶m,&dummy);
515 // if the fitting failed set noise and pedestal to 0
516 if ((ret==-4) || (ret==-1) || (ret==-2)) {
520 if((param[1]/10.0) > 65534.0) param[1] = 0.0;
521 if((param[2]/10.0) > 65534.0) param[2] = 0.0;
522 rocMean->SetValue(iChannel,param[1]/10.0);
523 rocRMS->SetValue(iChannel,param[2]/10.0);
526 // here we analyse doubled read channels
528 AliTRDCalROC *rocMeand = GetCalRocMeand(idet,kTRUE);
529 AliTRDCalROC *rocRMSd = GetCalRocRMSd(idet,kTRUE);
531 Int_t nrows = rocMeand->GetNrows();
532 Int_t shift = 144*nrows;
533 Int_t total = shift+3*8*nrows;
535 for (Int_t iChannel=shift; iChannel<total; iChannel++){
536 Int_t offset = (nbinsAdc+2)*(iChannel+1)+1;
537 Double_t ret = AliMathBase::FitGaus(arrayHP+offset,nbinsAdc,fAdcMin-0.5,fAdcMax-0.5,¶m,&dummy);
538 // if the fitting failed set noise and pedestal to 0
539 if ((ret==-4) || (ret==-1) || (ret==-2)) {
543 if((param[1]/10.0) > 65534.0) param[1] = 0.0;
544 if((param[2]/10.0) > 65534.0) param[2] = 0.0;
546 // here we have to recalculate backward
547 Int_t nb = iChannel-shift;
548 Int_t row = nb%nrows;
549 Int_t j = (Int_t)(nb/nrows);
551 Int_t icol = (Int_t)(j/8);
553 Int_t finalcol = 18*imcm;
554 if(icol > 0) icol += 17;
558 Int_t channel = row+finalcol*nrows;
560 //printf("iChannel %d, nrows %d, finalcol %d, row %d, channel %d\n",iChannel,nrows,finalcol,row,channel);
561 if((finalcol < 0) || (finalcol >= 144)) continue;
563 rocMeand->SetValue(channel,param[1]/10.0);
564 rocRMSd->SetValue(channel,param[2]/10.0);
571 //_______________________________________________________________________________________
572 AliTRDCalPadStatus* AliTRDCalibPadStatus::CreateCalPadStatus()
575 // Create Pad Status out of Mean and RMS values
576 // The chamber without data are masked, this is the corrected in the preprocessor
579 AliTRDCalPadStatus* obj = new AliTRDCalPadStatus("padstatus", "padstatus");
581 for (Int_t idet=0; idet<540; ++idet)
583 AliTRDCalSingleChamberStatus *calROC = obj->GetCalROC(idet);
586 if ( !GetCalRocMean(idet)) {
587 for(Int_t k = 0; k < calROC->GetNchannels(); k++){
588 calROC->SetStatus(k,AliTRDCalPadStatus::kNotConnected);
595 AliTRDCalROC *calRocMean = new AliTRDCalROC(*( (AliTRDCalROC *) GetCalRocMean(idet)));
596 AliTRDCalROC *calRocRMS = new AliTRDCalROC(*( (AliTRDCalROC *) GetCalRocRMS(idet)));
598 //Take the stuff second chance
599 AliTRDCalROC *calRocMeand = new AliTRDCalROC(*( (AliTRDCalROC *) GetCalRocMeand(idet)));
600 AliTRDCalROC *calRocRMSd = new AliTRDCalROC(*( (AliTRDCalROC *) GetCalRocRMSd(idet)));
603 calRocRMSd->Unfold();
607 Int_t row = calROC->GetNrows();
608 Int_t col = calROC->GetNcols();
610 Double_t rmsmean = calRocMean->GetRMS()*10.0;
611 Double_t meanmean = calRocMean->GetMean()*10.0;
612 Double_t meansquares = calRocRMS->GetMean();
615 for(Int_t irow = 0; irow < row; irow++){
618 Float_t meanprevious = 0.0;
619 Float_t rmsprevious = 0.0;
623 for(Int_t icol = 0; icol < col; icol++){
625 mean = calRocMean->GetValue(icol,irow)*10.0;
626 rms = calRocRMS->GetValue(icol,irow);
629 meanprevious = calRocMean->GetValue((icol -1),irow)*10.0;
630 rmsprevious = calRocRMS->GetValue((icol - 1),irow);
634 // masked if two noisy
635 if((rms <= 0.0001) || (TMath::Abs(mean-meanmean)>(5*rmsmean)) || (TMath::Abs(rms)>(5.0*TMath::Abs(meansquares)))) {
638 // look at second chance
639 Float_t meand = calRocMeand->GetValue(icol,irow)*10.0;
640 Float_t rmsd = calRocRMSd->GetValue(icol,irow);
642 if((rmsd <= 0.0001) || (TMath::Abs(meand-meanmean)>(5*rmsmean)) || (TMath::Abs(rmsd)>(5.0*TMath::Abs(meansquares)))) {
643 if((rmsd <= 0.0001) && (rms <= 0.0001)) {
644 calROC->SetStatus(icol,irow,AliTRDCalPadStatus::kNotConnected);
647 calROC->SetStatus(icol, irow, AliTRDCalPadStatus::kMasked);
651 calROC->SetStatus(icol, irow, AliTRDCalPadStatus::kReadSecond);
656 // bridge if previous pad found something
658 if((TMath::Abs(meanprevious-mean)<0.000001) && (TMath::Abs(rmsprevious-rms)<0.000001) && (mean > 0.0001)) {
659 //printf("mean previous %f, mean %f, rms %f, rmsprevious %f, col %d\n",meanprevious,mean,rms,rmsprevious,icol);
660 calROC->SetStatus(icol -1 ,irow, AliTRDCalPadStatus::kPadBridgedRight);
661 calROC->SetStatus(icol ,irow, AliTRDCalPadStatus::kPadBridgedLeft);
679 //_______________________________________________________________________________________
680 AliTRDCalPad* AliTRDCalibPadStatus::CreateCalPad()
683 // Create Pad Noise out of RMS values
686 AliTRDCalPad* obj = new AliTRDCalPad("PadNoise", "PadNoise");
689 for (Int_t det=0; det<AliTRDgeometry::kNdet; ++det) {
691 AliTRDCalROC *calROC22 = obj->GetCalROC(det);
693 AliTRDCalROC *calRocRMS = ((AliTRDCalROC *)GetCalRocRMS(det,kTRUE));
695 for(Int_t k = 0; k < calROC22->GetNchannels(); k++){
696 calROC22->SetValue(k,calRocRMS->GetValue(k));
705 //_______________________________________________________________________________________
706 AliTRDCalDet* AliTRDCalibPadStatus::CreateCalDet() const
709 // Create Det Noise correction factor
712 AliTRDCalDet* obj = new AliTRDCalDet("DetNoise", "DetNoise (correction factor)");
714 for(Int_t l = 0; l < 540; l++){
715 obj->SetValue(l,10.0);
722 //_____________________________________________________________________
723 void AliTRDCalibPadStatus::DumpToFile(const Char_t *filename, const Char_t *dir, Bool_t append) /*FOLD00*/
726 // Write class to file
737 TDirectory *backup = gDirectory;
738 TFile f(filename,option.Data());
740 if ( !sDir.IsNull() ){
741 f.mkdir(sDir.Data());
747 if ( backup ) backup->cd();
750 //_____________________________________________________________________
751 void AliTRDCalibPadStatus::SetCalRocMean(AliTRDCalROC *mean, Int_t det) /*FOLD00*/
754 // Put the AliTRDCalROC in the array fCalRocArrayMean
758 AliTRDCalROC *rocMean = GetCalRocMean(det,kTRUE);
760 Int_t nChannels = rocMean->GetNchannels();
762 for (Int_t iChannel=0; iChannel<nChannels; iChannel++){
764 rocMean->SetValue(iChannel,mean->GetValue(iChannel));
770 //_____________________________________________________________________
771 void AliTRDCalibPadStatus::SetCalRocRMS(AliTRDCalROC *rms, Int_t det) /*FOLD00*/
774 // Put the AliTRDCalROC in the array fCalRocArrayRMS
778 AliTRDCalROC *rocRms = GetCalRocRMS(det,kTRUE);
780 Int_t nChannels = rocRms->GetNchannels();
782 for (Int_t iChannel=0; iChannel<nChannels; iChannel++){
784 rocRms->SetValue(iChannel,rms->GetValue(iChannel));
789 //_____________________________________________________________________
790 void AliTRDCalibPadStatus::SetCalRocMeand(AliTRDCalROC *mean, Int_t det) /*FOLD00*/
793 // Put the AliTRDCalROC in the array fCalRocArrayMean
797 AliTRDCalROC *rocMean = GetCalRocMeand(det,kTRUE);
799 Int_t nChannels = rocMean->GetNchannels();
801 for (Int_t iChannel=0; iChannel<nChannels; iChannel++){
803 rocMean->SetValue(iChannel,mean->GetValue(iChannel));
809 //_____________________________________________________________________
810 void AliTRDCalibPadStatus::SetCalRocRMSd(AliTRDCalROC *rms, Int_t det) /*FOLD00*/
813 // Put the AliTRDCalROC in the array fCalRocArrayRMS
817 AliTRDCalROC *rocRms = GetCalRocRMSd(det,kTRUE);
819 Int_t nChannels = rocRms->GetNchannels();
821 for (Int_t iChannel=0; iChannel<nChannels; iChannel++){
823 rocRms->SetValue(iChannel,rms->GetValue(iChannel));
828 //_____________________________________________________________________________
829 Int_t AliTRDCalibPadStatus::GetLayer(Int_t d) const
832 // Reconstruct the layer number from the detector number
835 return ((Int_t) (d % 6));
839 //_____________________________________________________________________________
840 Int_t AliTRDCalibPadStatus::GetStack(Int_t d) const
843 // Reconstruct the chamber number from the detector number
846 return ((Int_t) (d % 30) / 6);
850 //_____________________________________________________________________________
851 Int_t AliTRDCalibPadStatus::GetSector(Int_t d) const
854 // Reconstruct the sector number from the detector number
857 return ((Int_t) (d / 30));