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