]> git.uio.no Git - u/mrichter/AliRoot.git/blame - RICH/AliRICHDigitizer.cxx
Digitizer class created. Code cleaning and bug fixes (J.Chudoba)
[u/mrichter/AliRoot.git] / RICH / AliRICHDigitizer.cxx
CommitLineData
b762c2f6 1/**************************************************************************
2 * Copyright(c) 1998-2000, 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/* $Log$ */
17#include <iostream>
18
19#include <TTree.h>
20#include <TObjArray.h>
21#include <TFile.h>
22#include <TDirectory.h>
23#include <TParticle.h>
24
25#include "AliRICHDigitizer.h"
26#include "AliRICHChamber.h"
27#include "AliHitMap.h"
28#include "AliRICHHitMapA1.h"
29#include "AliRICH.h"
30#include "AliRICHHit.h"
31#include "AliRICHSDigit.h"
32#include "AliRICHDigit.h"
33#include "AliRICHTransientDigit.h"
34#include "AliRun.h"
35#include "AliPDG.h"
36#include "AliRunDigitizer.h"
37
38ClassImp(AliRICHDigitizer)
39
40//___________________________________________
41 AliRICHDigitizer::AliRICHDigitizer()
42{
43// Default constructor - don't use it
44 fHits = 0;
45 fSDigits = 0;
46 fHitMap = 0;
47 fTDList = 0;
48}
49
50////////////////////////////////////////////////////////////////////////
51AliRICHDigitizer::AliRICHDigitizer(AliRunDigitizer* manager)
52 :AliDigitizer(manager)
53{
54// ctor which should be used
55 fHits = 0;
56 fSDigits = 0;
57 fHitMap = 0;
58 fTDList = 0;
59 fDebug = 0;
60 if (GetDebug()>2)
61 cerr<<"AliRICHDigitizer::AliRICHDigitizer"
62 <<"(AliRunDigitizer* manager) was processed"<<endl;
63}
64
65////////////////////////////////////////////////////////////////////////
66AliRICHDigitizer::~AliRICHDigitizer()
67{
68// Destructor
69 if (fHits) delete fHits;
70 if (fSDigits) delete fSDigits;
71 if (fHitMap) delete fHitMap;
72 if (fTDList) delete fTDList;
73}
74
75//------------------------------------------------------------------------
76Bool_t AliRICHDigitizer::Exists(const AliRICHSDigit *padhit)
77{
78 return (fHitMap[fNch]->TestHit(padhit->PadX(),padhit->PadY()));
79}
80
81//------------------------------------------------------------------------
82void AliRICHDigitizer::Update(AliRICHSDigit *padhit)
83{
84 AliRICHTransientDigit *pdigit =
85 static_cast<AliRICHTransientDigit*>(
86 fHitMap[fNch]->GetHit(padhit->PadX(),padhit->PadY()));
87
88 // update charge
89 //
90 Int_t iqpad = Int_t(padhit->QPad()); // charge per pad
91 pdigit->AddSignal(iqpad);
92 pdigit->AddPhysicsSignal(iqpad);
93
94 // update list of tracks
95 //
96 Int_t track, charge;
97 track = fTrack+fMask;
98 if (fSignal) {
99 charge = iqpad;
100 } else {
101 charge = kBgTag;
102 }
103 pdigit->UpdateTrackList(track,charge);
104}
105
106//------------------------------------------------------------------------
107void AliRICHDigitizer::CreateNew(AliRICHSDigit *padhit)
108{
109 fTDList->AddAtAndExpand(
110 new AliRICHTransientDigit(fNch,fDigits),fCounter);
111 fHitMap[fNch]->SetHit(padhit->PadX(),padhit->PadY(),fCounter);
112
113 AliRICHTransientDigit* pdigit =
114 static_cast<AliRICHTransientDigit*>(fTDList->Last());
115 // list of tracks
116 Int_t track, charge;
117 if (fSignal) {
118 track = fTrack;
119 charge = padhit->QPad();
120 } else {
121 track = kBgTag;
122 charge = kBgTag;
123 }
124 pdigit->AddToTrackList(track,charge);
125 fCounter++;
126}
127
128
129////////////////////////////////////////////////////////////////////////
130Bool_t AliRICHDigitizer::Init()
131{
132// Initialisation
133 fHits = new TClonesArray("AliRICHHit",1000);
134 fSDigits = new TClonesArray("AliRICHSDigit",1000);
135 return kTRUE;
136}
137
138////////////////////////////////////////////////////////////////////////
139//void AliRICHDigitizer::Digitise(Int_t nev, Int_t flag)
140void AliRICHDigitizer::Exec(Option_t* option)
141{
142 TString optionString = option;
143 if (optionString.Data() == "deb") {
144 cout<<"AliMUONDigitizer::Exec: called with option deb "<<endl;
145 fDebug = 3;
146 }
147
148 AliRICHChamber* iChamber;
149 AliSegmentation* segmentation;
150
151 fTDList = new TObjArray;
152
153 AliRICH *pRICH = (AliRICH *) gAlice->GetDetector("RICH");
154 pRICH->MakeBranchInTreeD(fManager->GetTreeD());
155 fHitMap= new AliHitMap* [kNCH];
156
157 for (Int_t i =0; i<kNCH; i++) {
158 iChamber= &(pRICH->Chamber(i));
159 segmentation=iChamber->GetSegmentationModel(1);
160 fHitMap[i] = new AliRICHHitMapA1(segmentation, fTDList);
161 }
162
163// Loop over files to digitize
164 fSignal = kTRUE;
165 fCounter = 0;
166 for (Int_t inputFile=0; inputFile<fManager->GetNinputs();
167 inputFile++) {
168
169// Connect RICH branches
170
171 if (inputFile > 0 ) fSignal = kFALSE;
172 TBranch *branchHits = 0;
173 TBranch *branchSDigits = 0;
174 TTree *treeH = fManager->GetInputTreeH(inputFile);
175 if (GetDebug()>2) {
176 cerr<<" inputFile "<<inputFile<<endl;
177 cerr<<" treeH, fHits "<<treeH<<" "<<fHits<<endl;
178 }
179 if (treeH && fHits) {
180 branchHits = treeH->GetBranch("RICH");
181 if (branchHits) {
182 fHits->Clear();
183 branchHits->SetAddress(&fHits);
184 }
185 else
186 Error("Exec","branch RICH was not found");
187 }
188 if (GetDebug()>2) cerr<<" branchHits = "<<branchHits<<endl;
189
190 if (treeH && fSDigits) {
191 branchSDigits = treeH->GetBranch("RICHSDigits");
192 if (branchSDigits)
193 branchSDigits->SetAddress(&fSDigits);
194 else
195 Error("exec","branch RICHSDigits was not found");
196 }
197 if (GetDebug()>2) cerr<<" branchSDigits = "<<branchSDigits<<endl;
198
199
200 //
201 // Loop over tracks
202 //
203
204 Int_t ntracks =(Int_t) treeH->GetEntries();
205 for (fTrack=0; fTrack<ntracks; fTrack++) {
206 fHits->Clear();
207 fSDigits->Clear();
208 branchHits->GetEntry(fTrack);
209 branchSDigits->GetEntry(fTrack);
210
211
212 //
213 // Loop over hits
214 for(Int_t i = 0; i < fHits->GetEntriesFast(); ++i) {
215 AliRICHHit* mHit = static_cast<AliRICHHit*>(fHits->At(i));
216 fNch = mHit->Chamber()-1; // chamber number
217 if (fNch >= kNCH) {
218 cerr<<"AliRICHDigitizer: chamber nr. fNch out of range: "<<fNch<<endl;
219 cerr<<" track: "<<fTrack<<endl;
220 continue;
221 }
222 iChamber = &(pRICH->Chamber(fNch));
223
224 //
225 // Loop over pad hits
226 for (AliRICHSDigit* mPad=
227 (AliRICHSDigit*)pRICH->FirstPad(mHit,fSDigits);
228 mPad;
229 mPad=(AliRICHSDigit*)pRICH->NextPad(fSDigits))
230 {
231 Int_t iqpad = mPad->QPad(); // charge per pad
232 fDigits[0]=mPad->PadX();
233 fDigits[1]=mPad->PadY();
234 fDigits[2]=iqpad;
235 fDigits[3]=iqpad;
236 fDigits[4]=mPad->HitNumber();
237
238 // build the list of fired pads and update the info
239 if (Exists(mPad)) {
240 Update(mPad);
241 } else {
242 CreateNew(mPad);
243 }
244 } //end loop over clusters
245 } // hit loop
246 } // track loop
247 } // end file loop
248 if (GetDebug()>2) cerr<<"END OF FILE LOOP"<<endl;
249
250 Int_t tracks[kMAXTRACKSPERRICHDIGIT];
251 Int_t charges[kMAXTRACKSPERRICHDIGIT];
252 Int_t nentries=fTDList->GetEntriesFast();
253
254 for (Int_t nent=0;nent<nentries;nent++) {
255 AliRICHTransientDigit *transDigit=(AliRICHTransientDigit*)fTDList->At(nent);
256 if (transDigit==0) continue;
257 Int_t ich=transDigit->GetChamber();
258 Int_t q=transDigit->Signal();
259 iChamber=&(pRICH->Chamber(ich));
260 AliRICHResponse * response=iChamber->GetResponseModel();
261 Int_t adcmax= (Int_t) response->MaxAdc();
262
263
264 // add white noise and do zero-suppression and signal truncation (new electronics,old electronics gaus 1.2,0.2)
265 //printf("Treshold: %d\n",iChamber->fTresh->GetHitIndex(transDigit->PadX(),transDigit->PadY()));
266
267// tmp change
268// Int_t pedestal = iChamber->fTresh->GetHitIndex(transDigit->PadX(),transDigit->PadY());
269 Int_t pedestal = 0;
270
271
272 //printf("Pedestal:%d\n",pedestal);
273 //Int_t pedestal=0;
274 Float_t treshold = (pedestal + 4*2.2);
275
276 Float_t meanNoise = gRandom->Gaus(2.2, 0.3);
277 Float_t noise = gRandom->Gaus(0, meanNoise);
278 q+=(Int_t)(noise + pedestal);
279 //q+=(Int_t)(noise);
280 // magic number to be parametrised !!!
281 if ( q <= treshold)
282 {
283 q = q - pedestal;
284 continue;
285 }
286 q = q - pedestal;
287 if ( q >= adcmax) q=adcmax;
288 fDigits[0]=transDigit->PadX();
289 fDigits[1]=transDigit->PadY();
290 fDigits[2]=q;
291 fDigits[3]=transDigit->Physics();
292 fDigits[4]=transDigit->Hit();
293
294 Int_t nptracks = transDigit->GetNTracks();
295
296 // this was changed to accomodate the real number of tracks
297 if (nptracks > kMAXTRACKSPERRICHDIGIT) {
298 printf("Attention - tracks > 10 %d\n",nptracks);
299 nptracks=kMAXTRACKSPERRICHDIGIT;
300 }
301 if (nptracks > 2) {
302 printf("Attention - tracks > 2 %d \n",nptracks);
303 }
304 for (Int_t tr=0;tr<nptracks;tr++) {
305 tracks[tr]=transDigit->GetTrack(tr);
306 charges[tr]=transDigit->GetCharge(tr);
307 } //end loop over list of tracks for one pad
308 if (nptracks < kMAXTRACKSPERRICHDIGIT ) {
309 for (Int_t t=nptracks; t<kMAXTRACKSPERRICHDIGIT; t++) {
310 tracks[t]=0;
311 charges[t]=0;
312 }
313 }
314 //write file
315 //if (ich==2)
316 //fprintf(points,"%4d, %4d, %4d\n",digits[0],digits[1],digits[2]);
317
318 // fill digits
319 pRICH->AddDigits(ich,tracks,charges,fDigits);
320 }
321 fManager->GetTreeD()->Fill();
322
323 pRICH->ResetDigits();
324 fTDList->Delete(); // or fTDList->Clear(); ???
325 for(Int_t ii=0;ii<kNCH;++ii) {
326 if (fHitMap[ii]) {
327 delete fHitMap[ii];
328 fHitMap[ii]=0;
329 }
330 }
331
332 TClonesArray *richDigits;
333 for (Int_t k=0;k<kNCH;k++) {
334 richDigits = pRICH->DigitsAddress(k);
335 Int_t ndigit=richDigits->GetEntriesFast();
336 printf ("Chamber %d digits %d \n",k,ndigit);
337 }
338 pRICH->ResetDigits(); /// ??? should it be here???
339 fManager->GetTreeD()->Write(0,TObject::kOverwrite);
340 delete [] fHitMap;
341 delete fTDList;
342
343 if (fHits) fHits->Delete();
344 if (fSDigits) fSDigits->Delete();
345
346}