1 /**************************************************************************
2 * Copyright(c) 2007, 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 **************************************************************************/
16 // --- ROOT system ---
17 #include "TClonesArray.h"
19 // --- AliRoot header files ---
20 #include "AliPHOSCpvRawDigiProducer.h"
21 #include "AliPHOSCpvRawStream.h"
22 #include "AliPHOSDigit.h"
23 #include "AliPHOSGeometry.h"
28 ClassImp(AliPHOSCpvRawDigiProducer);
30 //--------------------------------------------------------------------------------------
31 AliPHOSCpvRawDigiProducer::AliPHOSCpvRawDigiProducer():
38 fPedFilesRLoaded(kFALSE)
40 fGeom=AliPHOSGeometry::GetInstance() ;
41 if(!fGeom) fGeom = AliPHOSGeometry::GetInstance("IHEP");
44 // create a 2d array to store the pedestals
45 for (Int_t iDDL=0;iDDL<2*AliPHOSCpvParam::kNDDL;iDDL++){
46 fPed[0][iDDL] = new Int_t *[AliPHOSCpvParam::kPadPcX];
47 fPed[1][iDDL] = new Int_t *[AliPHOSCpvParam::kPadPcX];
48 for(Int_t ix=0; ix<AliPHOSCpvParam::kPadPcX; ix++) {
49 fPed[0][iDDL][ix] = new Int_t [AliPHOSCpvParam::kPadPcY];
50 fPed[1][iDDL][ix] = new Int_t [AliPHOSCpvParam::kPadPcY];
54 //-------------------------------------------------------------------------------------
55 AliPHOSCpvRawDigiProducer::AliPHOSCpvRawDigiProducer(AliRawReader * rawReader):
62 fPedFilesRLoaded(kFALSE)
64 fGeom=AliPHOSGeometry::GetInstance() ;
65 if(!fGeom) fGeom = AliPHOSGeometry::GetInstance("IHEP");
67 fRawStream = new AliPHOSCpvRawStream(rawReader);
68 fRawStream->SetTurbo(fTurbo);
70 // create a 2d array to store the pedestals
71 for (Int_t iDDL=0;iDDL<2*AliPHOSCpvParam::kNDDL;iDDL++) {
72 fPed[0][iDDL] = new Int_t *[AliPHOSCpvParam::kPadPcX];
73 fPed[1][iDDL] = new Int_t *[AliPHOSCpvParam::kPadPcX];
74 for(Int_t ix=0; ix<AliPHOSCpvParam::kPadPcX; ix++) {
75 fPed[0][iDDL][ix] = new Int_t [AliPHOSCpvParam::kPadPcY];
76 fPed[1][iDDL][ix] = new Int_t [AliPHOSCpvParam::kPadPcY];
80 //--------------------------------------------------------------------------------------
81 AliPHOSCpvRawDigiProducer::~AliPHOSCpvRawDigiProducer()
83 if(fRawStream) delete fRawStream;
84 if(fhErrors) delete fhErrors;
85 for(Int_t iDDL = 0;iDDL<2*AliPHOSCpvParam::kNDDL;iDDL++) {
86 for(Int_t ix=0; ix<AliPHOSCpvParam::kPadPcX; ix++) {
87 delete [] fPed[0][iDDL][ix];
88 delete [] fPed[1][iDDL][ix];
90 delete [] fPed[0][iDDL];
91 delete [] fPed[1][iDDL];
94 //--------------------------------------------------------------------------------------
95 Bool_t AliPHOSCpvRawDigiProducer::LoadPedFiles() {
96 // read pedestals from file
97 for(Int_t iDDL = 0;iDDL<2*AliPHOSCpvParam::kNDDL;iDDL+=2)
98 for(Int_t iCC=0; iCC<AliPHOSCpvParam::kNRows; iCC++) {
100 pedFile = fopen(Form("thr%d_%02d.dat",iDDL,iCC),"r");
102 Printf("AliPHOSCpvRawDigiProducer::LoadPedFiles: Error, file thr%d_%02d.dat could not be open",iDDL,iCC);
105 Int_t i3g = 0, iPad = 0;
107 while(!feof(pedFile)) {
108 Int_t abs = AliPHOSCpvParam::Abs(iDDL,iCC,i3g,iPad);
110 if(AliPHOSCpvParam::A2DDL(abs)!=iDDL)
111 cout<<"AliPHOSCpvRawDigiProducer::LoadPedFiles(): wrong connection table! abs = "
112 <<abs<<", DDL = "<<iDDL<<", A2DDL = "<<AliPHOSCpvParam::A2DDL(abs)<<endl;
113 if(AliPHOSCpvParam::A2CC(abs)!=iCC)
114 cout<<"AliPHOSCpvRawDigiProducer::LoadPedFiles(): wrong connection table! abs = "
115 <<abs<<", CC = "<< iCC <<", A2CC = "<<AliPHOSCpvParam::A2CC(abs)<<endl;
116 if(AliPHOSCpvParam::A23G(abs)!=i3g)
117 cout<<"AliPHOSCpvRawDigiProducer::LoadPedFiles(): wrong connection table! abs = "
118 <<abs<<", 3G = "<< i3g <<", A23G = "<<AliPHOSCpvParam::A23G(abs)<<endl;
119 if(AliPHOSCpvParam::A2Pad(abs)!=iPad)
120 cout<<"AliPHOSCpvRawDigiProducer::LoadPedFiles(): wrong connection table! abs = "
121 <<abs<<", Pad = "<< iPad <<", A2Pad = "<<AliPHOSCpvParam::A2Pad(abs)<<endl;
124 fscanf(pedFile,"%x",&thr);
125 if(AliPHOSCpvParam::IsValidAbs(abs)) {
126 Int_t s = thr & 0x1ff;
128 fPed[0][iDDL][AliPHOSCpvParam::A2X(abs)][AliPHOSCpvParam::A2Y(abs)] = p-s;
129 fPed[1][iDDL][AliPHOSCpvParam::A2X(abs)][AliPHOSCpvParam::A2Y(abs)] = s;
130 int testAbs = AliPHOSCpvParam::XY2A(iDDL,AliPHOSCpvParam::A2X(abs),AliPHOSCpvParam::A2Y(abs));
132 cout<<"AliPHOSCpvRawDigiProducer::LoadPedFiles(): wrong connection table! abs = "
133 <<abs<<", testAbs = "<<testAbs<<endl;
136 if(iPad == 64) {iPad = 0; i3g++;}
140 if(lineCnt < AliPHOSCpvParam::kN3GAdd * 64) return kFALSE;
142 fPedFilesRLoaded = kTRUE;
145 //--------------------------------------------------------------------------------------
146 void AliPHOSCpvRawDigiProducer::SetTurbo(Bool_t turbo)
149 if(fRawStream) fRawStream->SetTurbo(fTurbo);
151 //--------------------------------------------------------------------------------------
152 Bool_t AliPHOSCpvRawDigiProducer::LoadNewEvent(AliRawReader * rawReader)
154 if(fRawStream) delete fRawStream;
155 fRawStream = new AliPHOSCpvRawStream(rawReader);
157 fRawStream->SetTurbo(fTurbo);
163 //--------------------------------------------------------------------------------------
164 void AliPHOSCpvRawDigiProducer::MakeDigits(TClonesArray * digits) const
166 // Fills TClonesArray of AliPHOSDigits and
167 // returns histogram of error types
169 Int_t relId[4], absId=-1;
171 while (fRawStream->Next()) {
172 for (Int_t iPad=0; iPad<fRawStream->GetNPads(); iPad++) {
173 Float_t charge = fRawStream->GetChargeArray()[iPad];
174 Int_t aPad = fRawStream->GetPadArray()[iPad];
175 Int_t ix = AliPHOSCpvParam::A2X(aPad);
176 Int_t iy = AliPHOSCpvParam::A2Y(aPad);
177 Int_t iddl = AliPHOSCpvParam::A2DDL(aPad);
179 relId[0] = AliPHOSCpvParam::DDL2Mod(iddl) ; // counts from 1 to 5
180 relId[1] = -1; // -1=CPV
181 relId[2] = ix + 1; // counts from 1 to 128
182 relId[3] = iy + 1; // counts from 1 to 60
183 fGeom->RelToAbsNumbering(relId, absId);
185 AliDebug(2,Form("CPV digit: pad=%d, (x,z)=(%3d,%2d), DDL=%d, charge=%.0f",
186 aPad,ix,iy,iddl,charge));
188 if(fPedFilesRLoaded) {
189 if (charge > fPed[0][iddl][ix][iy] + fPed[1][iddl][ix][iy])
190 charge -= fPed[0][iddl][ix][iy];
194 if(charge < fCpvMinE) charge = 0;
195 if (charge>0) new((*digits)[iDigit++]) AliPHOSDigit(-1,absId,charge,0.);
198 } // while(fRawStream->Next())
200 for (Int_t i = 0 ; i < digits->GetEntriesFast() ; i++) {
201 AliPHOSDigit *digit = static_cast<AliPHOSDigit*>( digits->At(i) ) ;
202 digit->SetIndexInList(i) ;
205 AliDebug(1,Form("Array of %d CPV digits is created",digits->GetEntriesFast()));
207 // fill histogram of errors
208 for(Int_t iDDL=0; iDDL<2*AliPHOSCpvParam::kNDDL; iDDL++) {
209 Int_t nErrors = AliPHOSCpvRawStream::GetNErrors();
210 for(Int_t iType=0; iType<nErrors; iType++) { // iType - type of error
211 fhErrors -> Fill(iType+1,fRawStream -> GetErrors(iDDL,iType));
215 //--------------------------------------------------------------------------------------
216 void AliPHOSCpvRawDigiProducer::CreateErrHist()
218 Int_t nErrors = AliPHOSCpvRawStream::GetNErrors();
219 const char * errNames[nErrors];
220 for(Int_t i=0; i<nErrors; i++) {
221 errNames[i] = AliPHOSCpvRawStream::GetErrName(i);
223 fhErrors = new TH1I("errorTypes","Errors occured during processing",nErrors+1,0,nErrors+1);
224 TAxis* x = fhErrors->GetXaxis();
225 x->SetBinLabel(1, "Can't get event");
226 for(Int_t i=0; i<nErrors; i++) {
227 x->SetBinLabel(i+2,errNames[i]);
231 //--------------------------------------------------------------------------------------