]> git.uio.no Git - u/mrichter/AliRoot.git/blob - VZERO/AliVZERODigitizer.cxx
Modified plots and made jet finder use SDigits
[u/mrichter/AliRoot.git] / VZERO / AliVZERODigitizer.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 //_________________________________________________________________________
17 //
18 // This  constructs Digits out of Hits
19 //
20
21 // --- Standard library ---
22 #include <Riostream.h>
23 #include <stdlib.h>
24
25 // --- ROOT system ---
26 #include <TFile.h>
27 #include <TFolder.h>
28 #include <TROOT.h>
29 #include <TSystem.h>
30 #include <TTree.h>
31 #include <TVirtualMC.h>
32
33 // --- AliRoot header files ---
34 #include "AliVZEROConst.h"
35 #include "AliRun.h"
36 #include "AliVZERO.h"
37 #include "AliVZEROhit.h"
38 #include "AliHit.h"
39 #include "AliDetector.h"
40 #include "AliRunLoader.h"
41 #include "AliLoader.h"
42 #include "AliConfig.h"
43 #include "AliRunDigitizer.h"
44 #include "AliDigitizer.h"
45 #include "AliHeader.h"
46 #include "AliStack.h"
47 #include "AliVZERODigitizer.h"
48 #include "AliVZEROdigit.h"
49 #include "AliMC.h"
50
51 ClassImp(AliVZERODigitizer)
52
53  AliVZERODigitizer::AliVZERODigitizer()
54 {
55    if (!fDigits) fDigits = new TClonesArray("AliVZEROdigit", 1000);
56    
57    fNevents = 0 ;     
58    fDigits = 0 ;
59    fNdigits = 0;
60    fHits = 0 ;
61    fRunLoader = 0;
62   
63    fPhotoCathodeEfficiency =   0.18;
64    fPMVoltage              =  768.0;
65    fPMGain = TMath::Power((fPMVoltage / 112.5) ,7.04277);     
66 }
67
68 //____________________________________________________________________________ 
69   AliVZERODigitizer::AliVZERODigitizer(AliRunDigitizer* manager)
70                     :AliDigitizer(manager)
71                     
72 {
73   // constructor
74   
75   fNevents = 0;     
76   fDigits = 0;
77   fNdigits = 0;
78   fHits = 0;
79   fRunLoader = 0;
80   
81   fPhotoCathodeEfficiency =   0.18;
82   fPMVoltage              =  768.0;
83   fPMGain = TMath::Power( (fPMVoltage / 112.5) ,7.04277 );
84 }
85            
86 //____________________________________________________________________________ 
87   AliVZERODigitizer::~AliVZERODigitizer()
88 {
89   // destructor
90   
91   if (fDigits) {
92     fDigits->Delete();
93     delete fDigits;
94     fDigits=0; }
95 }
96
97 //____________________________________________________________________________
98 void AliVZERODigitizer::OpengAliceFile(const char *file)
99 {
100   // Loads galice.root file and corresponding header, kinematics
101   // hits and  digits 
102   
103   fRunLoader = AliRunLoader::Open(file,AliConfig::GetDefaultEventFolderName(),
104                                   "UPDATE");
105   
106   if (!fRunLoader)
107    {
108      Error("Open","Can not open session for file %s.",file);
109    }
110   
111   fRunLoader->LoadgAlice();
112   fRunLoader->LoadHeader();
113   fRunLoader->LoadKinematics();
114
115   gAlice = fRunLoader->GetAliRun();
116   
117   if (gAlice)
118     { printf("<AliVZEROdigitizer::Open> ");
119       printf("AliRun object found on file.\n");}
120   else
121     { printf("<AliVZEROdigitizer::Open> ");
122       printf("Could not find AliRun object.\n");}
123     
124   // Initialise Hit and Digit arrays
125   fHits   = new TClonesArray ("AliVZEROhit", 1000);
126   fDigits = new TClonesArray ("AliVZEROdigit", 1000);
127
128   fVZERO  = (AliVZERO*) gAlice->GetDetector("VZERO");
129   fVZEROLoader = fRunLoader->GetLoader("VZEROLoader");
130   
131   if (fVZEROLoader == 0x0){
132       cerr<<"Hits2Digits : Can not find VZERO or VZEROLoader\n";}
133     
134   Int_t retval = fVZEROLoader->LoadHits("read");
135   if (retval){
136      Error("Open","Error occured while loading hits... Exiting.");
137      return;}
138       
139   fVZEROLoader->LoadDigits("recreate");  
140 }
141
142 //____________________________________________________________________________
143 void AliVZERODigitizer::Exec() 
144  { 
145
146   fNdigits = 0;
147   Int_t N; 
148   Int_t map[96];
149   Int_t cell = 0;
150   Float_t cPM = fPhotoCathodeEfficiency * fPMGain;
151              
152   for(Int_t i=0; i<96; i++) map[i] = 0; 
153           
154   fNevents = (Int_t) fRunLoader->TreeE()->GetEntries ();  
155   printf(" Number of events in file =  %d \n", fNevents);
156   
157   for (Int_t ievent = 0; ievent < fNevents; ievent++){
158       
159       for(Int_t i=0; i<96; i++) map[i] = 0;     
160       
161       fRunLoader->GetEvent(ievent);
162       
163       fTreeH = fVZEROLoader->TreeH();
164       if (fTreeH == 0x0)
165        { Error("Exec","Cannot get TreeH");
166          return; }
167              
168       fTreeD = fVZEROLoader->TreeD();
169       if (fTreeD == 0x0)
170       { fVZEROLoader->MakeTree("D");
171         fVZEROLoader->MakeDigitsContainer();
172         fTreeD = fVZEROLoader->TreeD(); }
173         
174       Int_t bufsize = 16000;
175       fTreeD->Branch("VZERODigit", &fDigits, bufsize); 
176       
177 //    Now make Digits from hits
178
179       if (fVZERO)
180        {
181          fHits = fVZERO->Hits();
182          
183          Int_t ntracks = (Int_t) fTreeH->GetEntries ();
184 //       printf(" Number of Tracks in the TreeH = %d \n", ntracks);
185          for (Int_t track = 0; track < ntracks; track++)
186            {
187              gAlice->ResetHits ();
188              fTreeH->GetEvent(track);
189              fParticle = fRunLoader->Stack()->Particle(track);
190              Int_t nhits = fHits->GetEntriesFast();
191              for (Int_t hit = 0; hit < nhits; hit++)
192                  {
193                    fVZEROHit = (AliVZEROhit *)fHits->UncheckedAt(hit);
194                    N    = fVZEROHit->Nphot();
195                    cell = fVZEROHit->Cell();                                    
196                    map[cell] = map[cell] + N;
197                  }           // hit   loop
198            }                 // track loop
199             
200            Int_t icount = 0; 
201            
202            for(Int_t i=0; i<96; i++) {
203               Float_t q1 = Float_t ( map[i] )* cPM * kQe;
204               Float_t noise = gRandom->Gaus(10.5,3.22);
205               Float_t PMresponse  =  q1/kC*TMath::Power(ktheta/kthau,1/(1-ktheta/kthau)) 
206                                   + noise*1e-3;
207               map[i] = Int_t( PMresponse * 200.0);
208               if(map[i] > 3) {
209                  icount++;
210 //               printf(" Event, cell, adc = %d %d %d\n", ievent, i, map[i]);
211                  AddDigit(ievent, i, map[i]);} 
212             }
213
214            fTreeD->Reset();
215            fTreeD->Fill();
216            ResetDigit();
217            fVZEROLoader->WriteDigits("OVERWRITE");  
218            fVZEROLoader->UnloadDigits();     
219       }                     // VZERO loop
220     }  //event loop
221 }
222
223 //____________________________________________________________________________
224 void AliVZERODigitizer::AddDigit(Int_t eventnumber, Int_t cellnumber, Int_t adc) 
225  { 
226  
227 // Adds Digit 
228  
229   TClonesArray &ldigits = *fDigits;  
230   new(ldigits[fNdigits++]) AliVZEROdigit(eventnumber,cellnumber,adc);
231 }
232 //____________________________________________________________________________
233 void AliVZERODigitizer::ResetDigit()
234 {
235 // Clears Digits
236   fNdigits = 0;
237   if (fDigits) fDigits->Clear();
238 }