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"
28 // --- Standard library ---
30 // --- AliRoot header files ---
33 #include "AliRawReader.h"
34 #include "AliPMDRawStream.h"
35 #include "AliPMDddldata.h"
36 #include "AliPMDCalibGain.h"
38 ClassImp(AliPMDCalibGain)
40 AliPMDCalibGain::AliPMDCalibGain(): TObject()
42 // Standard Constructor
43 for(Int_t idet = 0; idet < kDet; idet++)
45 for(Int_t ismn = 0; ismn < kMaxSMN; ismn++)
47 fSMIso[idet][ismn] = 0.;
48 fSMCount[idet][ismn] = 0.;
49 for(Int_t jrow = 0; jrow < kMaxRow; jrow++)
51 for(Int_t kcol = 0; kcol < kMaxCol; kcol++)
53 fCellIso[idet][ismn][jrow][kcol] = 0.;
54 fCellCount[idet][ismn][jrow][kcol] = 0.;
61 // ------------------------------------------------------------------------ //
62 AliPMDCalibGain::AliPMDCalibGain(const AliPMDCalibGain &pmdcalibgain):
65 for(Int_t idet = 0; idet < kDet; idet++)
67 for(Int_t ismn = 0; ismn < kMaxSMN; ismn++)
69 fSMIso[idet][ismn] = pmdcalibgain.fSMIso[idet][ismn] ;
70 fSMCount[idet][ismn] = pmdcalibgain.fSMCount[idet][ismn] ;
71 for(Int_t jrow = 0; jrow < kMaxRow; jrow++)
73 for(Int_t kcol = 0; kcol < kMaxCol; kcol++)
75 fCellIso[idet][ismn][jrow][kcol] = pmdcalibgain.fCellIso[idet][ismn][jrow][kcol];
76 fCellCount[idet][ismn][jrow][kcol] = pmdcalibgain.fCellCount[idet][ismn][jrow][kcol];
83 // ------------------------------------------------------------------------ //
84 AliPMDCalibGain &AliPMDCalibGain::operator=(const AliPMDCalibGain &pmdcalibgain)
86 if(this != &pmdcalibgain)
88 for(Int_t idet = 0; idet < kDet; idet++)
90 for(Int_t ismn = 0; ismn < kMaxSMN; ismn++)
92 fSMIso[idet][ismn] = pmdcalibgain.fSMIso[idet][ismn];
93 fSMCount[idet][ismn] = pmdcalibgain.fSMCount[idet][ismn];
94 for(Int_t jrow = 0; jrow < kMaxRow;jrow++)
96 for(Int_t kcol = 0; kcol < kMaxCol; kcol++)
98 fCellIso[idet][ismn][jrow][kcol] =
99 pmdcalibgain.fCellIso[idet][ismn][jrow][kcol];
100 fCellCount[idet][ismn][jrow][kcol] =
101 pmdcalibgain.fCellCount[idet][ismn][jrow][kcol];
109 // ------------------------------------------------------------------------ //
110 AliPMDCalibGain::~AliPMDCalibGain()
115 // ------------------------------------------------------------------------ //
116 Bool_t AliPMDCalibGain::ProcessEvent(AliRawReader *rawReader)
118 // Calculates the ADC of isolated cell
120 TObjArray pmdddlcont;
122 const Int_t kDDL = AliDAQ::NumberOfDdls("PMD");
123 const Int_t kCellNeighbour = 6;
125 Int_t neibx[6] = {1,0,-1,-1,0,1};
126 Int_t neiby[6] = {0,1,1,0,-1,-1};
128 Int_t id1,jd1; //neighbour row/col
129 Int_t isocount; //number of neighbours with 0 signal
131 Float_t d1[kDet][kMaxSMN][kMaxRow][kMaxCol];
132 Bool_t streamout = kFALSE;
134 for(Int_t idet = 0; idet < kDet; idet++)
136 for(Int_t ismn = 0; ismn < kMaxSMN; ismn++)
138 for(Int_t irow = 0; irow < kMaxRow; irow++)
140 for(Int_t icol = 0; icol < kMaxCol; icol++)
142 d1[idet][ismn][irow][icol] = 0.;
148 AliPMDRawStream rawStream(rawReader);
150 for (Int_t iddl = 0; iddl < kDDL; iddl++)
152 rawReader->Select("PMD", iddl, iddl);
153 streamout = rawStream.DdlData(iddl,&pmdddlcont);
155 Int_t ientries = pmdddlcont.GetEntries();
156 for (Int_t ient = 0; ient < ientries; ient++)
158 AliPMDddldata *pmdddl = (AliPMDddldata*)pmdddlcont.UncheckedAt(ient);
160 Int_t idet = pmdddl->GetDetector();
161 Int_t ismn = pmdddl->GetSMN();
162 //Int_t mcm = pmdddl->GetMCM();
163 //Int_t ichno = pmdddl->GetChannel();
164 Int_t irow = pmdddl->GetRow();
165 Int_t icol = pmdddl->GetColumn();
166 Int_t isig = pmdddl->GetSignal();
170 d1[idet][ismn][irow][icol] = (Float_t) isig;
176 for(Int_t idet=0; idet < kDet; idet++)
178 for(Int_t ismn = 0; ismn < kMaxSMN; ismn++)
180 for(Int_t irow = 0; irow < kMaxRow; irow++)
182 for(Int_t icol = 0; icol < kMaxCol; icol++)
184 if(d1[idet][ismn][irow][icol] > 0)
187 for(Int_t ii = 0; ii < kCellNeighbour; ii++)
189 id1 = irow + neibx[ii];
190 jd1 = icol + neiby[ii];
191 if(d1[idet][ismn][id1][jd1] == 0)
194 if(isocount == kCellNeighbour)
197 fSMIso[idet][ismn] += d1[idet][ismn][irow][icol];
198 fCellIso[idet][ismn][irow][icol] += d1[idet][ismn][irow][icol];
200 fSMCount[idet][ismn]++;
201 fCellCount[idet][ismn][irow][icol]++;
205 } // neigh cell cond.
214 // ------------------------------------------------------------------------ //
215 void AliPMDCalibGain::Analyse(TTree *gaintree)
217 // Calculates the mean
218 Int_t det, sm, row, col;
220 Float_t modmean = 0.;
221 Float_t cellmean = 0.;
223 gaintree->Branch("det",&det,"det/I");
224 gaintree->Branch("sm",&sm,"sm/I");
225 gaintree->Branch("row",&row,"row/I");
226 gaintree->Branch("col",&col,"col/I");
227 gaintree->Branch("gain",&gain,"gain/F");
229 for(Int_t idet = 0; idet < kDet; idet++)
231 for(Int_t ism = 0; ism < kMaxSMN; ism++)
234 if (fSMCount[idet][ism] > 0)
235 modmean = fSMIso[idet][ism]/fSMCount[idet][ism];
237 for(Int_t irow = 0; irow < kMaxRow; irow++)
239 for(Int_t icol = 0; icol < kMaxCol; icol++)
241 if (fCellCount[idet][ism][irow][icol] > 0)
242 cellmean = fCellIso[idet][ism][irow][icol]/fCellCount[idet][ism][irow][icol];
253 gain = cellmean/modmean;
262 // ------------------------------------------------------------------------ //