1 /**************************************************************************
2 * Copyright(c) 1998-2000, 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 **************************************************************************/
18 Revision 1.5 2002/07/09 13:11:26 hristov
19 Old style include files needed on HP (aCC)
21 Revision 1.4 2002/05/28 07:53:10 morsch
22 Wrong order of arguments in for-statement corrected.
24 Revision 1.3 2001/12/05 14:53:34 hristov
25 Destructor corRevision 1.60 2002/10/22 16:28:21 alibrary
26 Introducing Riostream.h
28 Revision 1.59 2002/10/14 14:57:31 hristov
29 Merging the VirtualMC branch to the main development branch (HEAD)
31 Revision 1.58.6.1 2002/06/10 15:12:46 hristov
32 Merged with v3-08-02rected
34 Revision 1.2 2001/11/07 14:50:31 hristov
35 Minor correction of the Log part
37 Revision 1.1 2001/11/02 15:37:26 hristov
38 Digitizer class created. Code cleaning and bug fixes (J.Chudoba)
40 #include <Riostream.h>
43 #include <TObjArray.h>
45 #include <TDirectory.h>
46 #include <TParticle.h>
48 #include "AliRICHDigitizer.h"
49 #include "AliRICHChamber.h"
50 #include "AliHitMap.h"
51 #include "AliRICHHitMapA1.h"
53 #include "AliRICHHit.h"
54 #include "AliRICHSDigit.h"
55 #include "AliRICHDigit.h"
56 #include "AliRICHTransientDigit.h"
59 #include "AliRunDigitizer.h"
61 ClassImp(AliRICHDigitizer)
63 //___________________________________________
64 AliRICHDigitizer::AliRICHDigitizer()
66 // Default constructor - don't use it
73 ////////////////////////////////////////////////////////////////////////
74 AliRICHDigitizer::AliRICHDigitizer(AliRunDigitizer* manager)
75 :AliDigitizer(manager)
77 // ctor which should be used
84 cerr<<"AliRICHDigitizer::AliRICHDigitizer"
85 <<"(AliRunDigitizer* manager) was processed"<<endl;
88 ////////////////////////////////////////////////////////////////////////
89 AliRICHDigitizer::~AliRICHDigitizer()
100 for (Int_t i=0; i<kNCH; i++ )
109 //------------------------------------------------------------------------
110 Bool_t AliRICHDigitizer::Exists(const AliRICHSDigit *padhit)
112 return (fHitMap[fNch]->TestHit(padhit->PadX(),padhit->PadY()));
115 //------------------------------------------------------------------------
116 void AliRICHDigitizer::Update(AliRICHSDigit *padhit)
118 AliRICHTransientDigit *pdigit =
119 static_cast<AliRICHTransientDigit*>(
120 fHitMap[fNch]->GetHit(padhit->PadX(),padhit->PadY()));
124 Int_t iqpad = Int_t(padhit->QPad()); // charge per pad
125 pdigit->AddSignal(iqpad);
126 pdigit->AddPhysicsSignal(iqpad);
128 // update list of tracks
131 track = fTrack+fMask;
137 pdigit->UpdateTrackList(track,charge);
140 //------------------------------------------------------------------------
141 void AliRICHDigitizer::CreateNew(AliRICHSDigit *padhit)
143 fTDList->AddAtAndExpand(
144 new AliRICHTransientDigit(fNch,fDigits),fCounter);
145 fHitMap[fNch]->SetHit(padhit->PadX(),padhit->PadY(),fCounter);
147 AliRICHTransientDigit* pdigit =
148 static_cast<AliRICHTransientDigit*>(fTDList->Last());
153 charge = padhit->QPad();
158 pdigit->AddToTrackList(track,charge);
163 ////////////////////////////////////////////////////////////////////////
164 Bool_t AliRICHDigitizer::Init()
167 fHits = new TClonesArray("AliRICHHit",1000);
168 fSDigits = new TClonesArray("AliRICHSDigit",1000);
172 ////////////////////////////////////////////////////////////////////////
173 //void AliRICHDigitizer::Digitise(Int_t nev, Int_t flag)
174 void AliRICHDigitizer::Exec(Option_t* option)
176 TString optionString = option;
177 if (optionString.Data() == "deb") {
178 cout<<"AliMUONDigitizer::Exec: called with option deb "<<endl;
182 AliRICHChamber* iChamber;
183 AliSegmentation* segmentation;
185 fTDList = new TObjArray;
187 AliRICH *pRICH = (AliRICH *) gAlice->GetDetector("RICH");
188 pRICH->MakeBranchInTreeD(fManager->GetTreeD());
189 fHitMap= new AliHitMap* [kNCH];
191 for (Int_t i =0; i<kNCH; i++) {
192 iChamber= &(pRICH->Chamber(i));
193 segmentation=iChamber->GetSegmentationModel(1);
194 fHitMap[i] = new AliRICHHitMapA1(segmentation, fTDList);
197 // Loop over files to digitize
200 for (Int_t inputFile=0; inputFile<fManager->GetNinputs();
203 // Connect RICH branches
205 if (inputFile > 0 ) fSignal = kFALSE;
206 TBranch *branchHits = 0;
207 TBranch *branchSDigits = 0;
208 TTree *treeH = fManager->GetInputTreeH(inputFile);
210 cerr<<" inputFile "<<inputFile<<endl;
211 cerr<<" treeH, fHits "<<treeH<<" "<<fHits<<endl;
213 if (treeH && fHits) {
214 branchHits = treeH->GetBranch("RICH");
217 branchHits->SetAddress(&fHits);
220 Error("Exec","branch RICH was not found");
222 if (GetDebug()>2) cerr<<" branchHits = "<<branchHits<<endl;
224 if (treeH && fSDigits) {
225 branchSDigits = treeH->GetBranch("RICHSDigits");
227 branchSDigits->SetAddress(&fSDigits);
229 Error("exec","branch RICHSDigits was not found");
231 if (GetDebug()>2) cerr<<" branchSDigits = "<<branchSDigits<<endl;
238 Int_t ntracks =(Int_t) treeH->GetEntries();
239 for (fTrack=0; fTrack<ntracks; fTrack++) {
242 branchHits->GetEntry(fTrack);
243 branchSDigits->GetEntry(fTrack);
248 for(Int_t i = 0; i < fHits->GetEntriesFast(); ++i) {
249 AliRICHHit* mHit = static_cast<AliRICHHit*>(fHits->At(i));
250 fNch = mHit->Chamber()-1; // chamber number
252 cerr<<"AliRICHDigitizer: chamber nr. fNch out of range: "<<fNch<<endl;
253 cerr<<" track: "<<fTrack<<endl;
256 iChamber = &(pRICH->Chamber(fNch));
259 // Loop over pad hits
260 for (AliRICHSDigit* mPad=
261 (AliRICHSDigit*)pRICH->FirstPad(mHit,fSDigits);
263 mPad=(AliRICHSDigit*)pRICH->NextPad(fSDigits))
265 Int_t iqpad = mPad->QPad(); // charge per pad
266 fDigits[0]=mPad->PadX();
267 fDigits[1]=mPad->PadY();
270 fDigits[4]=mPad->HitNumber();
274 // build the list of fired pads and update the info
280 } //end loop over clusters
284 if (GetDebug()>2) cerr<<"END OF FILE LOOP"<<endl;
286 Int_t tracks[kMAXTRACKSPERRICHDIGIT];
287 Int_t charges[kMAXTRACKSPERRICHDIGIT];
288 Int_t nentries=fTDList->GetEntriesFast();
290 for (Int_t nent=0;nent<nentries;nent++) {
291 AliRICHTransientDigit *transDigit=(AliRICHTransientDigit*)fTDList->At(nent);
292 if (transDigit==0) continue;
293 Int_t ich=transDigit->GetChamber();
294 Int_t q=transDigit->Signal();
295 iChamber=&(pRICH->Chamber(ich));
296 AliRICHResponse * response=iChamber->GetResponseModel();
297 Int_t adcmax= (Int_t) response->MaxAdc();
300 // add white noise and do zero-suppression and signal truncation (new electronics,old electronics gaus 1.2,0.2)
301 //printf("Treshold: %d\n",iChamber->fTresh->GetHitIndex(transDigit->PadX(),transDigit->PadY()));
304 // Int_t pedestal = iChamber->fTresh->GetHitIndex(transDigit->PadX(),transDigit->PadY());
308 //printf("Pedestal:%d\n",pedestal);
310 Float_t treshold = (pedestal + 4*2.2);
312 Float_t meanNoise = gRandom->Gaus(2.2, 0.3);
313 Float_t noise = gRandom->Gaus(0, meanNoise);
314 q+=(Int_t)(noise + pedestal);
316 // magic number to be parametrised !!!
323 if ( q >= adcmax) q=adcmax;
324 fDigits[0]=transDigit->PadX();
325 fDigits[1]=transDigit->PadY();
327 fDigits[3]=transDigit->Physics();
328 fDigits[4]=transDigit->Hit();
330 Int_t nptracks = transDigit->GetNTracks();
332 // this was changed to accomodate the real number of tracks
333 if (nptracks > kMAXTRACKSPERRICHDIGIT) {
334 printf("Attention - tracks > 10 %d\n",nptracks);
335 nptracks=kMAXTRACKSPERRICHDIGIT;
338 printf("Attention - tracks > 2 %d \n",nptracks);
340 for (Int_t tr=0;tr<nptracks;tr++) {
341 tracks[tr]=transDigit->GetTrack(tr);
342 charges[tr]=transDigit->GetCharge(tr);
343 //printf("%f \n",charges[tr]);
344 } //end loop over list of tracks for one pad
345 if (nptracks < kMAXTRACKSPERRICHDIGIT ) {
346 for (Int_t t=nptracks; t<kMAXTRACKSPERRICHDIGIT; t++) {
353 //fprintf(points,"%4d, %4d, %4d\n",digits[0],digits[1],digits[2]);
356 pRICH->AddDigits(ich,tracks,charges,fDigits);
360 fManager->GetTreeD()->Fill();
362 //pRICH->ResetDigits();
363 fTDList->Delete(); // or fTDList->Clear(); ???
364 for(Int_t ii=0;ii<kNCH;++ii) {
371 TClonesArray *richDigits;
372 for (Int_t k=0;k<kNCH;k++) {
373 richDigits = pRICH->DigitsAddress(k);
374 Int_t ndigit=richDigits->GetEntriesFast();
375 printf ("Chamber %d digits %d \n",k,ndigit);
377 fManager->GetTreeD()->Write(0,TObject::kOverwrite);
378 pRICH->ResetDigits();
382 if (fHits) fHits->Delete();
383 if (fSDigits) fSDigits->Delete();