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