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