]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/PHOSrec/AliPHOSCpvRawDigiProducer.cxx
Coverity fixes
[u/mrichter/AliRoot.git] / PHOS / PHOSrec / AliPHOSCpvRawDigiProducer.cxx
1 /**************************************************************************
2  * Copyright(c) 2007, 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 // --- ROOT system ---
17 #include "TClonesArray.h"
18
19 // --- AliRoot header files ---
20 #include "AliPHOSCpvRawDigiProducer.h"
21 #include "AliPHOSDigit.h"
22 #include "AliPHOSCpvRawStream.h"
23 #include "AliLog.h"
24 #include<iostream>
25 using namespace std;
26
27 ClassImp(AliPHOSCpvRawDigiProducer);
28
29 //--------------------------------------------------------------------------------------
30 AliPHOSCpvRawDigiProducer::AliPHOSCpvRawDigiProducer():
31   TObject(),
32   fTurbo(kFALSE),
33   fCpvMinE(10.),
34   fRawStream(0),
35   fhErrors(0),
36   fPedFilesRLoaded(kFALSE)
37 {
38   CreateErrHist();
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];
46     }
47   }
48 }
49 //-------------------------------------------------------------------------------------
50 AliPHOSCpvRawDigiProducer::AliPHOSCpvRawDigiProducer(AliRawReader *& rawReader):
51   TObject(),
52   fTurbo(kFALSE),
53   fCpvMinE(10.),
54   fRawStream(0),
55   fhErrors(0),
56   fPedFilesRLoaded(kFALSE)
57 {
58   fRawStream = new AliPHOSCpvRawStream(rawReader);
59   fRawStream->SetTurbo(fTurbo);
60   CreateErrHist();
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];
68     }
69   }
70 }
71 //--------------------------------------------------------------------------------------
72 AliPHOSCpvRawDigiProducer::~AliPHOSCpvRawDigiProducer()
73 {
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];
80     }
81     delete [] ped[0][iDDL];
82     delete [] ped[1][iDDL];
83   }
84 }
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++) {
90       FILE * pedFile;
91       pedFile = fopen(Form("thr%d_%02d.dat",iDDL,iCC),"r");
92       if(!pedFile) {
93         Printf("AliPHOSCpvRawDigiProducer::LoadPedFiles: Error, file thr%d_%02d.dat could not be open",iDDL,iCC);
94         continue;
95         //return kFALSE;
96       }
97       Int_t i3g = 0, iPad = 0;
98       Int_t lineCnt = 0;
99       while(!feof(pedFile)) {
100         Int_t abs = AliPHOSCpvParam::Abs(iDDL,iCC,i3g,iPad);
101         if(iPad<48&&i3g<10){
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;
114         }
115         Int_t thr;
116         fscanf(pedFile,"%x",&thr);
117         if(AliPHOSCpvParam::IsValidAbs(abs)) {
118           Int_t s = thr & 0x1ff;
119           Int_t p = thr >> 9;
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));
123           if(abs!=testAbs)
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);
127         }
128         iPad++;
129         if(iPad == 64) {iPad = 0; i3g++;}
130         lineCnt++;
131       }
132       if(lineCnt < AliPHOSCpvParam::kN3GAdd * 64) return kFALSE;
133       fclose(pedFile);
134     }
135   fPedFilesRLoaded = kTRUE;
136   return kTRUE;
137 }
138 //--------------------------------------------------------------------------------------
139 void AliPHOSCpvRawDigiProducer::SetTurbo(Bool_t turbo) 
140 {
141   fTurbo = turbo;
142   if(fRawStream) fRawStream->SetTurbo(fTurbo);
143 }
144 //--------------------------------------------------------------------------------------
145 Bool_t AliPHOSCpvRawDigiProducer::LoadNewEvent(AliRawReader *& rawReader)
146 {
147   if(fRawStream) delete fRawStream;
148   fRawStream = new AliPHOSCpvRawStream(rawReader);
149   if(fRawStream) {
150     fRawStream->SetTurbo(fTurbo);
151     return kTRUE;
152   }
153   fhErrors->Fill(0);
154   return kFALSE;
155 }
156 //--------------------------------------------------------------------------------------
157 void AliPHOSCpvRawDigiProducer::MakeDigits(TClonesArray *& digits) const
158 {
159   // returns histogram of error types
160
161   if(digits) 
162     digits->Clear();
163   digits = new TClonesArray("AliPHOSDigit", AliPHOSCpvParam::kNDDL * AliPHOSCpvParam::kNRows * AliPHOSCpvParam::kN3GAdd * AliPHOSCpvParam::kNPadAdd);
164   Int_t iDigit = 0;
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];
176         }
177         else charge=0;
178         if(charge < fCpvMinE) charge = 0;
179       }
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);
183       
184     }
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));
192     }
193   }
194 }
195 //--------------------------------------------------------------------------------------
196 void AliPHOSCpvRawDigiProducer::CreateErrHist()
197 {
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);
202   }
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]);
208   }
209
210 }
211 //--------------------------------------------------------------------------------------