1 //////////////////////////////////////////////////////////////////////////////
3 // * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5 // * Author: The ALICE Off-line Project. *
6 // * Contributors are mentioned in the code where appropriate. *
8 // * Permission to use, copy, modify and distribute this software and its *
9 // * documentation strictly for non-commercial purposes is hereby granted *
10 // * without fee, provided that the above copyright notice appears in all *
11 // * copies and that both the copyright notice and this permission notice *
12 // * appear in the supporting documentation. The authors make no claims *
13 // * about the suitability of this software for any purpose. It is *
14 // * provided "as is" without express or implied warranty. *
15 // **************************************************************************/
17 //////////////////////////////////////////////////////////////////////////////
21 #include "TObjString.h"
23 #include "TClonesArray.h"
25 #include "TObjArray.h"
29 // --- Standard library ---
31 // --- AliRoot header files ---
34 #include "AliRawReader.h"
35 #include "AliPMDRawStream.h"
36 #include "AliPMDddldata.h"
37 #include "AliPMDCalibGain.h"
39 ClassImp(AliPMDCalibGain)
41 AliPMDCalibGain::AliPMDCalibGain():
45 // Standard Constructor
46 for(Int_t idet = 0; idet < kDet; idet++)
48 for(Int_t ismn = 0; ismn < kMaxSMN; ismn++)
50 fSMIso[idet][ismn] = 0.;
51 fSMCount[idet][ismn] = 0.;
52 fCountSm[idet][ismn]=0.;
53 fTempnhit[idet][ismn]=0.;
54 fTempnhitSq[idet][ismn]=0.;
55 for(Int_t jrow = 0; jrow < kMaxRow; jrow++)
57 for(Int_t kcol = 0; kcol < kMaxCol; kcol++)
59 fCellIso[idet][ismn][jrow][kcol] = 0.;
60 fCellCount[idet][ismn][jrow][kcol] = 0.;
61 fNhitCell[idet][ismn][jrow][kcol] = 0.;
62 fPedMeanRMS[idet][ismn][jrow][kcol] = 0.;
63 fHotFlag[idet][ismn][jrow][kcol] = 0.;
71 // ------------------------------------------------------------------------ //
72 AliPMDCalibGain::AliPMDCalibGain(const AliPMDCalibGain &pmdcalibgain):
73 TObject(pmdcalibgain),
76 for(Int_t idet = 0; idet < kDet; idet++)
78 for(Int_t ismn = 0; ismn < kMaxSMN; ismn++)
80 fSMIso[idet][ismn] = pmdcalibgain.fSMIso[idet][ismn] ;
81 fSMCount[idet][ismn] = pmdcalibgain.fSMCount[idet][ismn] ;
82 fCountSm[idet][ismn] = pmdcalibgain.fCountSm[idet][ismn];
83 fTempnhit[idet][ismn] = pmdcalibgain.fTempnhit[idet][ismn];
84 fTempnhitSq[idet][ismn] = pmdcalibgain.fTempnhitSq[idet][ismn];
85 for(Int_t jrow = 0; jrow < kMaxRow; jrow++)
87 for(Int_t kcol = 0; kcol < kMaxCol; kcol++)
89 fCellIso[idet][ismn][jrow][kcol] = pmdcalibgain.fCellIso[idet][ismn][jrow][kcol];
90 fCellCount[idet][ismn][jrow][kcol] = pmdcalibgain.fCellCount[idet][ismn][jrow][kcol];
91 fNhitCell[idet][ismn][jrow][kcol] = pmdcalibgain.fNhitCell[idet][ismn][jrow][kcol];
92 fPedMeanRMS[idet][ismn][jrow][kcol] = pmdcalibgain.fPedMeanRMS[idet][ismn][jrow][kcol];
93 fHotFlag[idet][ismn][jrow][kcol] = pmdcalibgain.fHotFlag[idet][ismn][jrow][kcol];
101 // ------------------------------------------------------------------------ //
102 AliPMDCalibGain &AliPMDCalibGain::operator=(const AliPMDCalibGain &pmdcalibgain)
104 if(this != &pmdcalibgain)
106 this->fpw = pmdcalibgain.fpw;
107 for(Int_t idet = 0; idet < kDet; idet++)
109 for(Int_t ismn = 0; ismn < kMaxSMN; ismn++)
111 fSMIso[idet][ismn] = pmdcalibgain.fSMIso[idet][ismn];
112 fSMCount[idet][ismn] = pmdcalibgain.fSMCount[idet][ismn];
113 fCountSm[idet][ismn] = pmdcalibgain.fCountSm[idet][ismn];
114 fTempnhit[idet][ismn] = pmdcalibgain.fTempnhit[idet][ismn];
115 fTempnhitSq[idet][ismn] = pmdcalibgain.fTempnhitSq[idet][ismn];
116 for(Int_t jrow = 0; jrow < kMaxRow;jrow++)
118 for(Int_t kcol = 0; kcol < kMaxCol; kcol++)
120 fCellIso[idet][ismn][jrow][kcol] = pmdcalibgain.fCellIso[idet][ismn][jrow][kcol];
121 fCellCount[idet][ismn][jrow][kcol] = pmdcalibgain.fCellCount[idet][ismn][jrow][kcol];
122 fNhitCell[idet][ismn][jrow][kcol] = pmdcalibgain.fNhitCell[idet][ismn][jrow][kcol]; //za
123 fPedMeanRMS[idet][ismn][jrow][kcol] = pmdcalibgain.fPedMeanRMS[idet][ismn][jrow][kcol];
124 fHotFlag[idet][ismn][jrow][kcol] = pmdcalibgain.fHotFlag[idet][ismn][jrow][kcol];
133 // ------------------------------------------------------------------------ //
134 AliPMDCalibGain::~AliPMDCalibGain()
139 // ------------------------------------------------------------------------ //
141 Int_t AliPMDCalibGain::ExtractPedestal(const Char_t *rootFile)
143 // Pedestal extraction from the PMD_PED.root file
144 // To be called once at the beginning
146 Int_t det, sm, row, col;
149 TFile *pedfile = new TFile(rootFile);
153 printf("ERROR --- NO PEDESTAL (PMD_PED1.root) FILE IS FOUND --- STOP GAIN DA\n");
158 TTree *ped =(TTree*)pedfile->Get("ped");
160 ped->SetBranchAddress("det",&det);
161 ped->SetBranchAddress("sm",&sm);
162 ped->SetBranchAddress("row",&row);
163 ped->SetBranchAddress("col",&col);
164 ped->SetBranchAddress("mean",&mean);
165 ped->SetBranchAddress("rms",&rms);
167 Int_t nentries = (Int_t)ped->GetEntries();
169 for (Int_t ient = 0; ient < nentries; ient++)
172 fPedMeanRMS[det][sm][row][col] = mean + 3.*rms;
173 //printf("Mean= %f, RMS= %f, PedMeanRMS=%f\n",mean,rms,fPedMeanRMS[det][sm][row][col]);
183 //------------------------------------------------------------------------------------------------
185 Int_t AliPMDCalibGain::ExtractHotChannel(const Char_t *rootFile)
187 // HotChannel extraction from the PMD_HOT.root file
188 // To be called once at the beginning
190 Int_t det, sm, row, col;
193 TFile *hotmapfile = new TFile(rootFile);
197 printf(" NO HOTCHANNEL MAP FOUND (PMD_HOR.root) FILE IS FOUND \n");
198 fHotFlag[kDet][kMaxSMN][kMaxRow][kMaxCol] = 0.;
202 TTree *hot =(TTree*)hotmapfile->Get("hot");
204 hot->SetBranchAddress("det",&det);
205 hot->SetBranchAddress("sm",&sm);
206 hot->SetBranchAddress("row",&row);
207 hot->SetBranchAddress("col",&col);
208 hot->SetBranchAddress("flag",&flag);
210 Int_t nentries = (Int_t)hot->GetEntries();
212 for (Int_t ient = 0; ient < nentries; ient++)
215 fHotFlag[det][sm][row][col] = flag;
217 //printf(" HotFlag=%f\n",fHotFlag[det][sm][row][col]);
228 // ------------------------------------------------------------------------ //
230 void AliPMDCalibGain::ReadTempFile(const Char_t *tempFile)
232 // Read the variables from the file
234 fpw = fopen(tempFile,"r");
236 Float_t smcount, smiso;
237 Float_t cellcount, celliso;
240 for (Int_t idet = 0; idet < kDet; idet++)
242 for (Int_t ism = 0; ism < kMaxSMN; ism++)
244 fscanf(fpw,"%d %d %f %f",&idet,&ism,&smcount,&smiso);
246 fSMCount[idet][ism] = smcount;
247 fSMIso[idet][ism] = smiso;
251 for (Int_t idet = 0; idet < kDet; idet++)
253 for (Int_t ism = 0; ism < kMaxSMN; ism++)
255 for (Int_t irow = 0; irow < kMaxRow; irow++)
257 for (Int_t icol = 0; icol < kMaxCol; icol++)
259 fscanf(fpw,"%d %d %d %d %f %f",&idet,&ism,&irow,&icol,
260 &cellcount,&celliso);
262 fCellCount[idet][ism][irow][icol] = cellcount;
263 fCellIso[idet][ism][irow][icol] = celliso;
272 // ------------------------------------------------------------------------ //
273 void AliPMDCalibGain::WriteTempFile(const Char_t *tempFile)
275 // Write the Temporary file if the required statics is not achieved
279 Following variables to be written to a file
282 fCellIso[idet][ismn][irow][icol];
285 fSMCount[idet][ismn];
286 fCellCount[idet][ismn][irow][icol];
290 fpw = fopen(tempFile,"w+");
292 for (Int_t idet = 0; idet < kDet; idet++)
294 // fprintf(fpw,"%d %f %f\n",idet,fDetCount[idet],fDetIso[idet]);
297 for (Int_t idet = 0; idet < kDet; idet++)
299 for (Int_t ism = 0; ism < kMaxSMN; ism++)
301 fprintf(fpw,"%d %d %f %f\n",idet,ism, fSMCount[idet][ism],fSMIso[idet][ism]);
305 for (Int_t idet = 0; idet < kDet; idet++)
307 for (Int_t ism = 0; ism < kMaxSMN; ism++)
309 for (Int_t irow = 0; irow < kMaxRow; irow++)
311 for (Int_t icol = 0; icol < kMaxCol; icol++)
313 fprintf(fpw,"%d %d %d %d %f %f\n",idet,ism,irow,icol,
314 fCellCount[idet][ism][irow][icol],
315 fCellIso[idet][ism][irow][icol]);
325 // ------------------------------------------------------------------------ //
327 Bool_t AliPMDCalibGain::ProcessEvent(AliRawReader *rawReader, TObjArray *pmdddlcont)
329 // Calculates the ADC of isolated cell
331 const Int_t kDDL = AliDAQ::NumberOfDdls("PMD");
332 const Int_t kCellNeighbour = 6;
333 Int_t neibx[6] = {1,0,-1,-1,0,1};
334 Int_t neiby[6] = {0,1,1,0,-1,-1};
336 Int_t id1,jd1; //neighbour row/col
337 Int_t isocount; //number of neighbours with 0 signal
339 Float_t d1[kDet][kMaxSMN][kMaxRow][kMaxCol];
341 for(Int_t idet = 0; idet < kDet; idet++)
343 for(Int_t ismn = 0; ismn < kMaxSMN; ismn++)
345 for(Int_t irow = 0; irow < kMaxRow; irow++)
347 for(Int_t icol = 0; icol < kMaxCol; icol++)
349 d1[idet][ismn][irow][icol] = 0.;
355 AliPMDRawStream rawStream(rawReader);
359 Int_t numberofDDLs = 0;
361 while ((iddl = rawStream.DdlData(pmdddlcont)) >=0) {
364 Int_t ientries = pmdddlcont->GetEntries();
366 for (Int_t ient = 0; ient < ientries; ient++)
368 AliPMDddldata *pmdddl = (AliPMDddldata*)pmdddlcont->UncheckedAt(ient);
370 Int_t idet = pmdddl->GetDetector();
371 Int_t ismn = pmdddl->GetSMN();
372 Int_t mcm = pmdddl->GetMCM();
373 //Int_t ichno = pmdddl->GetChannel();
374 Int_t irow = pmdddl->GetRow();
375 Int_t icol = pmdddl->GetColumn();
376 Int_t isig = pmdddl->GetSignal();
378 // This is the protection not to crash the code
380 if(mcm == 0) continue;
381 if (irow < 0 || icol < 0 || irow > 47 || icol > 95) continue;
383 // Pedestal subtraction
384 if(fHotFlag[idet][ismn][irow][icol] == 1.0) isig = 0;
388 d1[idet][ismn][irow][icol] =
389 (Float_t) isig - fPedMeanRMS[idet][ismn][irow][icol];
391 //printf("Signal_ped_subtracted=%f, pedestal=%f\n",d1[idet][ismn][irow][icol]),fPedMeanRMS[idet][ismn][irow][icol];
393 fNhitCell[idet][ismn][irow][icol]++; //cell hit frequency
398 pmdddlcont->Delete();
401 for(Int_t idet=0; idet < kDet; idet++)
403 for(Int_t ismn = 0; ismn < kMaxSMN; ismn++)
405 for(Int_t irow = 0; irow < kMaxRow; irow++)
407 for(Int_t icol = 0; icol < kMaxCol; icol++)
409 if(d1[idet][ismn][irow][icol] > 0)
412 for(Int_t ii = 0; ii < kCellNeighbour; ii++)
414 id1 = irow + neibx[ii];
415 jd1 = icol + neiby[ii];
416 if(d1[idet][ismn][id1][jd1] == 0)
419 if(isocount == kCellNeighbour)
421 //fDetIso[idet] += d1[idet][ismn][irow][icol];
422 fSMIso[idet][ismn] += d1[idet][ismn][irow][icol];
423 fCellIso[idet][ismn][irow][icol] += d1[idet][ismn][irow][icol];
425 fSMCount[idet][ismn]++;
426 fCellCount[idet][ismn][irow][icol]++;
430 } // neigh cell cond.
437 for(Int_t idet=0; idet < kDet; idet++)
439 for(Int_t ismn = 0; ismn < kMaxSMN; ismn++)
441 for(Int_t irow = 0; irow < kMaxRow; irow++)
443 for(Int_t icol = 0; icol < kMaxCol; icol++)
445 if(fNhitCell[idet][ismn][irow][icol]>0)
447 fCountSm[idet][ismn] += 1;
448 fTempnhit[idet][ismn] += fNhitCell[idet][ismn][irow][icol];
449 fTempnhitSq[idet][ismn] += fNhitCell[idet][ismn][irow][icol]
450 *fNhitCell[idet][ismn][irow][icol];
457 if (numberofDDLs < kDDL) return kFALSE;
461 // ------------------------------------------------------------------------ //
462 void AliPMDCalibGain::Analyse(TTree *gaintree)
464 // Calculates the mean
465 Int_t det, sm, row, col;
467 Float_t modmean = 0.;
468 Float_t cellmean = 0.;
470 gaintree->Branch("det",&det,"det/I");
471 gaintree->Branch("sm",&sm,"sm/I");
472 gaintree->Branch("row",&row,"row/I");
473 gaintree->Branch("col",&col,"col/I");
474 gaintree->Branch("gain",&gain,"gain/F");
476 for(Int_t idet = 0; idet < kDet; idet++)
478 for(Int_t ism = 0; ism < kMaxSMN; ism++)
480 if (fSMCount[idet][ism] > 0)
481 modmean = fSMIso[idet][ism]/fSMCount[idet][ism];
482 for(Int_t irow = 0; irow < kMaxRow; irow++)
484 for(Int_t icol = 0; icol < kMaxCol; icol++)
486 if (fCellCount[idet][ism][irow][icol] > 0.)
488 cellmean = fCellIso[idet][ism][irow][icol]/fCellCount[idet][ism][irow][icol];
494 if (cellmean > 0.0 && fCellCount[idet][ism][irow][icol]>0.)
496 gain = cellmean/modmean;
502 //if(fCellCount[idet][ism][irow][icol]>0.) printf("CellCount =%f, gain= %f\n",fCellCount[idet][ism][irow][icol],gain);
510 // ------------------------------------------------------------------------ //
511 void AliPMDCalibGain::AnalyseHotCell(TTree *hottree)
513 // Calculates the mean
514 Int_t det, sm, row, col;
518 Float_t sigmanhit,nhitcut;
520 hottree->Branch("det",&det,"det/I");
521 hottree->Branch("sm",&sm,"sm/I");
522 hottree->Branch("row",&row,"row/I");
523 hottree->Branch("col",&col,"col/I");
524 hottree->Branch("flag",&flag,"flag/F");
526 for(Int_t idet = 0; idet < kDet; idet++)
528 for(Int_t ism = 0; ism < kMaxSMN; ism++)
530 if (fCountSm[idet][ism]> 0)
532 meannhit = fTempnhit[idet][ism]/fCountSm[idet][ism];
533 meanSqnhit = fTempnhitSq[idet][ism]/fCountSm[idet][ism];
534 sigmanhit = sqrt(meanSqnhit-(meannhit*meannhit));
535 nhitcut = meannhit + 6.*sigmanhit;
537 for(Int_t irow = 0; irow < kMaxRow; irow++)
539 for(Int_t icol = 0; icol < kMaxCol; icol++)
546 if(fNhitCell[idet][ism][irow][icol] > nhitcut)