]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/PHOSbase/AliPHOSCpvRawDigiProducer.cxx
593b9f8d37f1cac487c4e8be5ba4048449ad0a81
[u/mrichter/AliRoot.git] / PHOS / PHOSbase / 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 "AliPHOSCpvRawStream.h"
22 #include "AliPHOSDigit.h"
23 #include "AliPHOSGeometry.h"
24 #include "AliLog.h"
25 #include<iostream>
26 using namespace std;
27
28 ClassImp(AliPHOSCpvRawDigiProducer);
29
30 //--------------------------------------------------------------------------------------
31 AliPHOSCpvRawDigiProducer::AliPHOSCpvRawDigiProducer():
32   TObject(),
33   fGeom(0),
34   fTurbo(kFALSE),
35   fCpvMinE(10.),
36   fRawStream(0),
37   fhErrors(0),
38   fPedFilesRLoaded(kFALSE)
39 {
40   fGeom=AliPHOSGeometry::GetInstance() ;
41   if(!fGeom) fGeom = AliPHOSGeometry::GetInstance("IHEP");
42
43   CreateErrHist();
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];
51     }
52   }
53 }
54 //-------------------------------------------------------------------------------------
55 AliPHOSCpvRawDigiProducer::AliPHOSCpvRawDigiProducer(AliRawReader * rawReader):
56   TObject(),
57   fGeom(0),
58   fTurbo(kFALSE),
59   fCpvMinE(10.),
60   fRawStream(0),
61   fhErrors(0),
62   fPedFilesRLoaded(kFALSE)
63 {
64   fGeom=AliPHOSGeometry::GetInstance() ;
65   if(!fGeom) fGeom = AliPHOSGeometry::GetInstance("IHEP");
66
67   fRawStream = new AliPHOSCpvRawStream(rawReader);
68   fRawStream->SetTurbo(fTurbo);
69   CreateErrHist();
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];
77     }
78   }
79 }
80 //--------------------------------------------------------------------------------------
81 AliPHOSCpvRawDigiProducer::~AliPHOSCpvRawDigiProducer()
82 {
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];
89     }
90     delete [] fPed[0][iDDL];
91     delete [] fPed[1][iDDL];
92   }
93 }
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++) {
99       FILE * pedFile;
100       pedFile = fopen(Form("thr%d_%02d.dat",iDDL,iCC),"r");
101       if(!pedFile) {
102         Printf("AliPHOSCpvRawDigiProducer::LoadPedFiles: Error, file thr%d_%02d.dat could not be open",iDDL,iCC);
103         continue;
104       }
105       Int_t i3g = 0, iPad = 0;
106       Int_t lineCnt = 0;
107       while(!feof(pedFile)) {
108         Int_t abs = AliPHOSCpvParam::Abs(iDDL,iCC,i3g,iPad);
109         if(iPad<48&&i3g<10){
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;
122         }
123         Int_t thr;
124         fscanf(pedFile,"%x",&thr);
125         if(AliPHOSCpvParam::IsValidAbs(abs)) {
126           Int_t s = thr & 0x1ff;
127           Int_t p = thr >> 9;
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));
131           if(abs!=testAbs)
132             cout<<"AliPHOSCpvRawDigiProducer::LoadPedFiles(): wrong connection table! abs = "
133                 <<abs<<", testAbs = "<<testAbs<<endl;
134         }
135         iPad++;
136         if(iPad == 64) {iPad = 0; i3g++;}
137         lineCnt++;
138       }
139       fclose(pedFile);
140       if(lineCnt < AliPHOSCpvParam::kN3GAdd * 64) return kFALSE;
141     }
142   fPedFilesRLoaded = kTRUE;
143   return kTRUE;
144 }
145 //--------------------------------------------------------------------------------------
146 void AliPHOSCpvRawDigiProducer::SetTurbo(Bool_t turbo) 
147 {
148   fTurbo = turbo;
149   if(fRawStream) fRawStream->SetTurbo(fTurbo);
150 }
151 //--------------------------------------------------------------------------------------
152 Bool_t AliPHOSCpvRawDigiProducer::LoadNewEvent(AliRawReader * rawReader)
153 {
154   if(fRawStream) delete fRawStream;
155   fRawStream = new AliPHOSCpvRawStream(rawReader);
156   if(fRawStream) {
157     fRawStream->SetTurbo(fTurbo);
158     return kTRUE;
159   }
160   fhErrors->Fill(0);
161   return kFALSE;
162 }
163 //--------------------------------------------------------------------------------------
164 void AliPHOSCpvRawDigiProducer::MakeDigits(TClonesArray * digits) const
165 {
166   // Fills TClonesArray of AliPHOSDigits and 
167   // returns histogram of error types
168
169   Int_t relId[4], absId=-1;
170   Int_t iDigit = 0;
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);
178
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);
184
185       AliDebug(2,Form("CPV digit: pad=%d, (x,z)=(%3d,%2d), DDL=%d, charge=%.0f",
186                       aPad,ix,iy,iddl,charge)); 
187
188       if(fPedFilesRLoaded) {
189         if (charge > fPed[0][iddl][ix][iy] + fPed[1][iddl][ix][iy])
190           charge -= fPed[0][iddl][ix][iy];
191         else 
192           charge  = 0;
193       }
194       if(charge < fCpvMinE) charge = 0;
195       if (charge>0) new((*digits)[iDigit++]) AliPHOSDigit(-1,absId,charge,0.);
196       
197     }
198   } // while(fRawStream->Next())
199   digits->Sort() ;
200   for (Int_t i = 0 ; i < digits->GetEntriesFast() ; i++) { 
201     AliPHOSDigit *digit = static_cast<AliPHOSDigit*>( digits->At(i) ) ; 
202     digit->SetIndexInList(i) ;     
203   }
204
205   AliDebug(1,Form("Array of %d CPV digits is created",digits->GetEntriesFast())); 
206
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));
212     }
213   }
214 }
215 //--------------------------------------------------------------------------------------
216 void AliPHOSCpvRawDigiProducer::CreateErrHist()
217 {
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);
222   }
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]);
228   }
229
230 }
231 //--------------------------------------------------------------------------------------