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 "AliPHOSDigit.h"
22 #include "AliPHOSCpvRawStream.h"
27 ClassImp(AliPHOSCpvRawDigiProducer);
29 //--------------------------------------------------------------------------------------
30 AliPHOSCpvRawDigiProducer::AliPHOSCpvRawDigiProducer():
36 fPedFilesRLoaded(kFALSE)
39 // create a 2d array to store the pedestals
40 for (Int_t iDDL=0;iDDL<AliPHOSCpvParam::kNDDL;iDDL++){
41 ped[0][iDDL] = new Int_t *[AliPHOSCpvParam::kPadPcX];
42 ped[1][iDDL] = new Int_t *[AliPHOSCpvParam::kPadPcX];
43 for(Int_t ix=0; ix<AliPHOSCpvParam::kPadPcX; ix++) {
44 ped[0][iDDL][ix] = new Int_t [AliPHOSCpvParam::kPadPcY];
45 ped[1][iDDL][ix] = new Int_t [AliPHOSCpvParam::kPadPcY];
49 //-------------------------------------------------------------------------------------
50 AliPHOSCpvRawDigiProducer::AliPHOSCpvRawDigiProducer(AliRawReader *& rawReader):
56 fPedFilesRLoaded(kFALSE)
58 fRawStream = new AliPHOSCpvRawStream(rawReader);
59 fRawStream->SetTurbo(fTurbo);
61 // create a 2d array to store the pedestals
62 for (Int_t iDDL=0;iDDL<AliPHOSCpvParam::kNDDL;iDDL++) {
63 ped[0][iDDL] = new Int_t *[AliPHOSCpvParam::kPadPcX];
64 ped[1][iDDL] = new Int_t *[AliPHOSCpvParam::kPadPcX];
65 for(Int_t ix=0; ix<AliPHOSCpvParam::kPadPcX; ix++) {
66 ped[0][iDDL][ix] = new Int_t [AliPHOSCpvParam::kPadPcY];
67 ped[1][iDDL][ix] = new Int_t [AliPHOSCpvParam::kPadPcY];
71 //--------------------------------------------------------------------------------------
72 AliPHOSCpvRawDigiProducer::~AliPHOSCpvRawDigiProducer()
74 if(fRawStream) delete fRawStream;
75 if(fhErrors) delete fhErrors;
76 for(Int_t iDDL = 0;iDDL<AliPHOSCpvParam::kNDDL;iDDL++) {
77 for(Int_t ix=0; ix<AliPHOSCpvParam::kPadPcX; ix++) {
78 delete [] ped[0][iDDL][ix];
79 delete [] ped[1][iDDL][ix];
81 delete [] ped[0][iDDL];
82 delete [] ped[1][iDDL];
85 //--------------------------------------------------------------------------------------
86 Bool_t AliPHOSCpvRawDigiProducer::LoadPedFiles() {
87 // read pedestals from file
88 for(Int_t iDDL = 0;iDDL<AliPHOSCpvParam::kNDDL;iDDL++)
89 for(Int_t iCC=0; iCC<AliPHOSCpvParam::kNRows; iCC++) {
91 pedFile = fopen(Form("thr%d_%02d.dat",iDDL,iCC),"r");
93 Printf("AliPHOSCpvRawDigiProducer::LoadPedFiles: Error, file thr%d_%02d.dat could not be open",iDDL,iCC);
97 Int_t i3g = 0, iPad = 0;
99 while(!feof(pedFile)) {
100 Int_t abs = AliPHOSCpvParam::Abs(iDDL,iCC,i3g,iPad);
102 if(AliPHOSCpvParam::A2DDL(abs)!=iDDL)
103 cout<<"AliPHOSCpvRawDigiProducer::LoadPedFiles(): wrong connection table! abs = "
104 <<abs<<", DDL = "<<iDDL<<", A2DDL = "<<AliPHOSCpvParam::A2DDL(abs)<<endl;
105 if(AliPHOSCpvParam::A2CC(abs)!=iCC)
106 cout<<"AliPHOSCpvRawDigiProducer::LoadPedFiles(): wrong connection table! abs = "
107 <<abs<<", CC = "<< iCC <<", A2CC = "<<AliPHOSCpvParam::A2CC(abs)<<endl;
108 if(AliPHOSCpvParam::A23G(abs)!=i3g)
109 cout<<"AliPHOSCpvRawDigiProducer::LoadPedFiles(): wrong connection table! abs = "
110 <<abs<<", 3G = "<< i3g <<", A23G = "<<AliPHOSCpvParam::A23G(abs)<<endl;
111 if(AliPHOSCpvParam::A2Pad(abs)!=iPad)
112 cout<<"AliPHOSCpvRawDigiProducer::LoadPedFiles(): wrong connection table! abs = "
113 <<abs<<", Pad = "<< iPad <<", A2Pad = "<<AliPHOSCpvParam::A2Pad(abs)<<endl;
116 fscanf(pedFile,"%x",&thr);
117 if(AliPHOSCpvParam::IsValidAbs(abs)) {
118 Int_t s = thr & 0x1ff;
120 ped[0][iDDL][AliPHOSCpvParam::A2X(abs)][AliPHOSCpvParam::A2Y(abs)] = p-s;
121 ped[1][iDDL][AliPHOSCpvParam::A2X(abs)][AliPHOSCpvParam::A2Y(abs)] = s;
122 int testAbs = AliPHOSCpvParam::XY2A(iDDL,AliPHOSCpvParam::A2X(abs),AliPHOSCpvParam::A2Y(abs));
124 cout<<"AliPHOSCpvRawDigiProducer::LoadPedFiles(): wrong connection table! abs = "
125 <<abs<<", testAbs = "<<testAbs<<endl;
126 //Printf("ped[%d][%d] = %d, pad = %d, abs = %d",AliPHOSCpvParam::A2X(abs),AliPHOSCpvParam::A2Y(abs), p + s, iPad, abs);
129 if(iPad == 64) {iPad = 0; i3g++;}
132 if(lineCnt < AliPHOSCpvParam::kN3GAdd * 64) return kFALSE;
135 fPedFilesRLoaded = kTRUE;
138 //--------------------------------------------------------------------------------------
139 void AliPHOSCpvRawDigiProducer::SetTurbo(Bool_t turbo)
142 if(fRawStream) fRawStream->SetTurbo(fTurbo);
144 //--------------------------------------------------------------------------------------
145 Bool_t AliPHOSCpvRawDigiProducer::LoadNewEvent(AliRawReader *& rawReader)
147 if(fRawStream) delete fRawStream;
148 fRawStream = new AliPHOSCpvRawStream(rawReader);
150 fRawStream->SetTurbo(fTurbo);
156 //--------------------------------------------------------------------------------------
157 void AliPHOSCpvRawDigiProducer::MakeDigits(TClonesArray *& digits) const
159 // returns histogram of error types
163 digits = new TClonesArray("AliPHOSDigit", AliPHOSCpvParam::kNDDL * AliPHOSCpvParam::kNRows * AliPHOSCpvParam::kN3GAdd * AliPHOSCpvParam::kNPadAdd);
165 while(fRawStream->Next()) {
166 for(Int_t iPad=0;iPad<fRawStream->GetNPads();iPad++) {
167 Int_t charge = fRawStream->GetChargeArray()[iPad];
168 Int_t aPad = fRawStream->GetPadArray()[iPad];
169 //cout<<"AliPHOSCpvRawDigiProducer::MakeDigits(): I've got pad "<<aPad<< "with amplitude "<<charge<<endl;
170 if(fPedFilesRLoaded) {
171 Int_t ix = AliPHOSCpvParam::A2X(aPad);
172 Int_t iy = AliPHOSCpvParam::A2Y(aPad);
173 Int_t iddl = AliPHOSCpvParam::A2DDL(aPad);
174 if (charge>ped[0][iddl][ix][iy]+ped[1][iddl][ix][iy]){
175 charge -=ped[0][iddl][ix][iy];
178 if(charge < fCpvMinE) charge = 0;
180 // if(charge) new((*digits)[iDigit++]) AliPHOSDigit(AliPHOSCpvParam::A2X(aPad),AliPHOSCpvParam::A2Y(aPad),charge);
181 // Check what is aPad! YK 31.12.2014
182 if(charge) new((*digits)[iDigit++]) AliPHOSDigit(-1,aPad,charge,0);
185 } // while(fRawStream->Next())
186 //cout<<"AliPHOSCpvRawDigiProducer::MakeDigits(): I've created "<<iDigit<<" digits."<<endl;
187 // fill histogram of errors
188 for(Int_t iDDL=0; iDDL<AliPHOSCpvParam::kNDDL; iDDL++) {
189 Int_t nErrors = AliPHOSCpvRawStream::GetNErrors();
190 for(Int_t iType=0; iType<nErrors; iType++) { // iType - type of error
191 fhErrors -> Fill(iType+1,fRawStream -> GetErrors(iDDL,iType));
195 //--------------------------------------------------------------------------------------
196 void AliPHOSCpvRawDigiProducer::CreateErrHist()
198 Int_t nErrors = AliPHOSCpvRawStream::GetNErrors();
199 const char * errNames[nErrors];
200 for(Int_t i=0; i<nErrors; i++) {
201 errNames[i] = AliPHOSCpvRawStream::GetErrName(i);
203 fhErrors = new TH1I("errorTypes","Errors occured during processing",nErrors+1,0,nErrors+1);
204 TAxis* x = fhErrors->GetXaxis();
205 x->SetBinLabel(1, "Can't get event");
206 for(Int_t i=0; i<nErrors; i++) {
207 x->SetBinLabel(i+2,errNames[i]);
211 //--------------------------------------------------------------------------------------