]> git.uio.no Git - u/mrichter/AliRoot.git/blob - VZERO/AliVZERO.cxx
VZERO reco moved to LocalEventReconstruction method of AliReconstruction, i.e. the...
[u/mrichter/AliRoot.git] / VZERO / AliVZERO.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 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////
19 //                                                                       //
20 //                          V-Zero   Detector                            //
21 //  This class contains the base procedures for the VZERO  detector      //
22 //  Default geometry of November 2003 :   V0R box is 4.4 cm thick        //
23 //                                  scintillators are 2 cm thick         //
24 //  All comments should be sent to Brigitte CHEYNIS :                    //
25 //                                 b.cheynis@ipnl.in2p3.fr               //
26 //                                                                       //
27 //                                                                       //
28 ///////////////////////////////////////////////////////////////////////////
29
30
31 // --- Standard libraries ---
32 #include <Riostream.h>
33 #include <stdlib.h>
34
35 // --- ROOT libraries ---
36 #include <TNamed.h>
37 #include "TROOT.h"
38 #include "TFile.h"
39 #include "TNetFile.h"
40 #include "TRandom.h"
41 #include "TTree.h"
42 #include "TBranch.h"
43 #include "TClonesArray.h"
44 #include "TStopwatch.h"
45
46 // --- AliRoot header files ---
47 #include "AliRun.h"
48 #include "AliMC.h"
49 #include "AliVZERO.h"
50 #include "AliVZEROLoader.h"
51 #include "AliVZERODigitizer.h"
52 #include "AliVZEROBuffer.h"
53 #include "AliRunDigitizer.h"
54 #include "AliVZEROdigit.h"
55 #include "AliDAQ.h"
56
57 ClassImp(AliVZERO)
58  //__________________________________________________________________
59 AliVZERO::AliVZERO(): AliDetector(),
60           fIdSens1(0),
61           fThickness(0.),
62           fThickness1(0.),
63           fMaxStepQua(0.),
64           fMaxStepAlu(0.),
65           fMaxDestepQua(0.),
66           fMaxDestepAlu(0.)
67 {
68 /// Default Constructor
69     
70     AliDebug(1,Form("default (empty) ctor this=%p",this));
71     fIshunt          = 0;         
72 }
73 //_____________________________________________________________________________
74 AliVZERO::AliVZERO(const char *name, const char *title)
75        : AliDetector(name,title),
76          fIdSens1(0),
77          fThickness(4.4),
78          fThickness1(2.0),
79          fMaxStepQua(0.05),
80          fMaxStepAlu(0.01),
81          fMaxDestepQua(-1.0),
82          fMaxDestepAlu(-1.0)
83 {
84   
85   // Standard constructor for VZERO Detector
86   
87   AliDebug(1,Form("ctor this=%p",this));
88   
89   //  fIshunt       =  1;  // All hits are associated with primary particles  
90    
91   fHits         =  new TClonesArray("AliVZEROhit", 400);
92   fDigits       =  new TClonesArray("AliVZEROdigit",400); 
93    
94   gAlice->GetMCApp()->AddHitList(fHits);
95
96 //   fThickness    =  4.4;   // total thickness of the V0R box in cm
97 //   fThickness1   =  2.0;   // thickness of scintillating cells in cm
98 //   
99 //   fMaxStepQua   =  0.05; 
100 //   fMaxStepAlu   =  0.01; 
101 //   
102 //   fMaxDestepQua =  -1.0;
103 //   fMaxDestepAlu =  -1.0;
104   
105   
106 }
107
108 //_____________________________________________________________________________
109 AliVZERO::~AliVZERO()
110 {
111   //
112   // Default destructor for VZERO Detector
113   //
114   
115     if (fHits) {
116         fHits->Delete();
117         delete fHits;
118         fHits=0; }
119     
120     if (fDigits) {
121         fDigits->Delete();
122         delete fDigits;
123         fDigits=0; }
124 }
125
126 //_____________________________________________________________________________
127 void AliVZERO::BuildGeometry()
128 {
129   //
130   // Builds simple ROOT TNode geometry for event display
131   //
132 }
133  
134 //_____________________________________________________________________________
135 void AliVZERO::CreateGeometry()
136 {
137   //
138   // Builds simple Geant3 geometry 
139   //
140 }
141 //_____________________________________________________________________________
142 void AliVZERO::CreateMaterials()
143 {
144   //
145   // Creates materials used for Geant3 geometry 
146   //
147 }
148
149 //_____________________________________________________________________________
150 Int_t AliVZERO::DistanceToPrimitive(Int_t /*px*/, Int_t /*py*/)
151 {
152   //
153   // Calculates the distance from the mouse to the VZERO on the screen
154   // Dummy routine
155   //
156   
157   return 9999;
158 }
159  
160 //_____________________________________________________________________________
161 void AliVZERO::Init()
162 {
163   //
164   // Initialises the VZERO  class after it has been built
165   //
166 }
167
168
169 //_____________________________________________________________________________
170 void AliVZERO::SetMaxStepQua(Float_t p1)
171 {
172   //
173   // Possible parametrisation of steps in active materials
174   //
175      fMaxStepQua = p1;
176 }
177
178 //_____________________________________________________________________________
179 void AliVZERO::SetMaxStepAlu(Float_t p1)
180 {
181   //
182   // Possible parametrisation of steps in Aluminum foils (not used in 
183   // version v2)
184   //
185     fMaxStepAlu = p1;
186 }
187
188 //_____________________________________________________________________________
189 void AliVZERO::SetMaxDestepQua(Float_t p1)
190 {
191   //
192   // Possible parametrisation of steps in active materials (quartz)
193   //
194     fMaxDestepQua = p1;
195 }
196
197 //_____________________________________________________________________________
198 void AliVZERO::SetMaxDestepAlu(Float_t p1)
199 {
200   //
201   // Possible parametrisation of steps in Aluminum (not used in 
202   // version v2)
203   //
204     fMaxDestepAlu = p1;
205 }
206
207 //_____________________________________________________________________________
208 AliLoader* AliVZERO::MakeLoader(const char* topfoldername)
209
210   //
211   // Builds VZEROgetter (AliLoader type)
212   // if detector wants to use customized getter, it must overload this method
213   //
214 //  Info("MakeLoader","Creating AliVZEROLoader. Top folder is %s.",topfoldername); 
215  
216   AliDebug(1,Form("Creating AliVZEROLoader, Top folder is %s ",topfoldername));
217   fLoader = new AliVZEROLoader(GetName(),topfoldername);
218   return fLoader;
219 }
220
221 //_____________________________________________________________________________
222 void AliVZERO::SetTreeAddress()
223 {
224   //
225   // Sets tree address for hits.
226   //
227   if (fLoader->TreeH() && (fHits == 0x0))
228     fHits = new  TClonesArray("AliVZEROhit", 400);
229
230   AliDetector::SetTreeAddress();
231 }
232
233 //_____________________________________________________________________________
234 AliDigitizer* AliVZERO::CreateDigitizer(AliRunDigitizer* manager) const
235 {
236   //
237   // Creates a digitizer for VZERO
238   //
239   return new AliVZERODigitizer(manager);
240 }
241
242 //_____________________________________________________________________________
243 void AliVZERO::Hits2Digits(){
244   //
245   // Converts hits to digits of the current event
246   //
247   // Inputs file name
248   const char *alifile = "galice.root";   
249
250   // Create the run digitizer 
251   AliRunDigitizer* manager = new AliRunDigitizer(1, 1);
252   manager->SetInputStream(0, alifile);
253   manager->SetOutputFile("H2Dfile");
254
255   // Creates the VZERO digitizer 
256   AliVZERODigitizer* dig = new AliVZERODigitizer(manager);
257
258   // Creates the digits
259   dig->Exec("");
260
261 }
262 //_____________________________________________________________________________
263 void AliVZERO::Digits2Raw()
264 {
265   //
266   // Converts digits of the current event to raw data
267   //
268   AliVZERO *fVZERO = (AliVZERO*)gAlice->GetDetector("VZERO");
269   fLoader->LoadDigits();
270   TTree* digits = fLoader->TreeD();
271   if (!digits) {
272     Error("Digits2Raw", "no digits tree");
273     return;
274   }
275   TClonesArray * VZEROdigits = new TClonesArray("AliVZEROdigit",1000);
276   fVZERO->SetTreeAddress();             
277   digits->GetBranch("VZERODigit")->SetAddress(&VZEROdigits); 
278   
279   const char *fileName    = AliDAQ::DdlFileName("VZERO",0);
280   AliVZEROBuffer* buffer  = new AliVZEROBuffer(fileName);
281   
282   //  Verbose level
283   //  0: Silent
284   //  1: cout messages
285   //  2: txt files with digits 
286   //  BE CAREFUL, verbose level 2 MUST be used only for debugging and
287   //  it is highly suggested to use this mode only for debugging digits files
288   //  reasonably small, because otherwise the size of the txt files can reach
289   //  quickly several MB wasting time and disk space.
290   
291   ofstream ftxt;
292   buffer->SetVerbose(0);
293   Int_t fVerbose = buffer->GetVerbose();
294
295   Int_t nEntries = Int_t(digits->GetEntries());
296   
297   for (Int_t i = 0; i < nEntries; i++) {
298   
299     fVZERO->ResetDigits();
300     digits->GetEvent(i);
301     Int_t ndig = VZEROdigits->GetEntriesFast(); 
302    
303     if(ndig == 0) continue;
304     if(fVerbose == 2) {ftxt.open("VZEROdigits.txt",ios::app);}
305     for(Int_t k=0; k<ndig; k++){
306         AliVZEROdigit* fVZERODigit = (AliVZEROdigit*) VZEROdigits->At(k);                       
307         Int_t ADC       = fVZERODigit->ADC();
308         Int_t PMNumber  = fVZERODigit->PMNumber();
309         Int_t Time      = fVZERODigit->Time();
310         if(fVerbose == 1) { cout <<"DDL: "<<fileName<< "\tdigit number: "<< k<<"\tPM number: "
311                             <<PMNumber<<"\tADC: "<< ADC << "\tTime: "<< Time << endl;} 
312         if(fVerbose == 2) {
313             ftxt<<"DDL: "<<fileName<< "\tdigit number: "<< k<<"\tPM number: "
314                            <<PMNumber<<"\tADC: "<< ADC << "\tTime: "<< Time << endl;          
315         }
316         buffer->WriteBinary(PMNumber, ADC, Time);
317     }
318   if(fVerbose==2) ftxt.close();
319   }
320
321   delete buffer;
322   fLoader->UnloadDigits();
323 }
324
325