]> git.uio.no Git - u/mrichter/AliRoot.git/blob - VZERO/AliVZERO.cxx
cover case for AOD analysis
[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 #include "AliRawReader.h"
57 #include "AliVZERORawStream.h"
58
59 ClassImp(AliVZERO)
60  //__________________________________________________________________
61 AliVZERO::AliVZERO(): AliDetector(),
62           fIdSens1(0),
63           fThickness(0.),
64           fThickness1(0.),
65           fMaxStepQua(0.),
66           fMaxStepAlu(0.),
67           fMaxDestepQua(0.),
68           fMaxDestepAlu(0.)
69 {
70 /// Default Constructor
71     
72     AliDebug(1,Form("default (empty) ctor this=%p",this));
73     fIshunt          = 0;         
74 }
75 //_____________________________________________________________________________
76 AliVZERO::AliVZERO(const char *name, const char *title)
77        : AliDetector(name,title),
78          fIdSens1(0),
79          fThickness(4.4),
80          fThickness1(2.0),
81          fMaxStepQua(0.05),
82          fMaxStepAlu(0.01),
83          fMaxDestepQua(-1.0),
84          fMaxDestepAlu(-1.0)
85 {
86   
87   // Standard constructor for VZERO Detector
88   
89   AliDebug(1,Form("ctor this=%p",this));
90   
91   //  fIshunt       =  1;  // All hits are associated with primary particles  
92    
93   fHits         =  new TClonesArray("AliVZEROhit", 400);
94   fDigits       =  new TClonesArray("AliVZEROdigit",400); 
95    
96   gAlice->GetMCApp()->AddHitList(fHits);
97
98 //   fThickness    =  4.4;   // total thickness of the V0R box in cm
99 //   fThickness1   =  2.0;   // thickness of scintillating cells in cm
100 //   
101 //   fMaxStepQua   =  0.05; 
102 //   fMaxStepAlu   =  0.01; 
103 //   
104 //   fMaxDestepQua =  -1.0;
105 //   fMaxDestepAlu =  -1.0;
106   
107   
108 }
109
110 //_____________________________________________________________________________
111 AliVZERO::~AliVZERO()
112 {
113   //
114   // Default destructor for VZERO Detector
115   //
116   
117     if (fHits) {
118         fHits->Delete();
119         delete fHits;
120         fHits=0; }
121     
122     if (fDigits) {
123         fDigits->Delete();
124         delete fDigits;
125         fDigits=0; }
126 }
127
128 //_____________________________________________________________________________
129 void AliVZERO::CreateGeometry()
130 {
131   //
132   // Builds simple Geant3 geometry 
133   //
134 }
135 //_____________________________________________________________________________
136 void AliVZERO::CreateMaterials()
137 {
138   //
139   // Creates materials used for Geant3 geometry 
140   //
141 }
142
143 //_____________________________________________________________________________
144 void AliVZERO::Init()
145 {
146   //
147   // Initialises the VZERO  class after it has been built
148   //
149 }
150
151
152 //_____________________________________________________________________________
153 void AliVZERO::SetMaxStepQua(Float_t p1)
154 {
155   //
156   // Possible parametrisation of steps in active materials
157   //
158      fMaxStepQua = p1;
159 }
160
161 //_____________________________________________________________________________
162 void AliVZERO::SetMaxStepAlu(Float_t p1)
163 {
164   //
165   // Possible parametrisation of steps in Aluminum foils (not used in 
166   // version v2)
167   //
168     fMaxStepAlu = p1;
169 }
170
171 //_____________________________________________________________________________
172 void AliVZERO::SetMaxDestepQua(Float_t p1)
173 {
174   //
175   // Possible parametrisation of steps in active materials (quartz)
176   //
177     fMaxDestepQua = p1;
178 }
179
180 //_____________________________________________________________________________
181 void AliVZERO::SetMaxDestepAlu(Float_t p1)
182 {
183   //
184   // Possible parametrisation of steps in Aluminum (not used in 
185   // version v2)
186   //
187     fMaxDestepAlu = p1;
188 }
189
190 //_____________________________________________________________________________
191 AliLoader* AliVZERO::MakeLoader(const char* topfoldername)
192
193   //
194   // Builds VZEROgetter (AliLoader type)
195   // if detector wants to use customized getter, it must overload this method
196   //
197 //  Info("MakeLoader","Creating AliVZEROLoader. Top folder is %s.",topfoldername); 
198  
199   AliDebug(1,Form("Creating AliVZEROLoader, Top folder is %s ",topfoldername));
200   fLoader = new AliVZEROLoader(GetName(),topfoldername);
201   return fLoader;
202 }
203
204 //_____________________________________________________________________________
205 void AliVZERO::SetTreeAddress()
206 {
207   //
208   // Sets tree address for hits.
209   //
210   if (fLoader->TreeH() && (fHits == 0x0))
211     fHits = new  TClonesArray("AliVZEROhit", 400);
212
213   AliDetector::SetTreeAddress();
214 }
215
216 //_____________________________________________________________________________
217 AliDigitizer* AliVZERO::CreateDigitizer(AliRunDigitizer* manager) const
218 {
219   //
220   // Creates a digitizer for VZERO
221   //
222   return new AliVZERODigitizer(manager);
223 }
224
225 //_____________________________________________________________________________
226 void AliVZERO::Hits2Digits(){
227   //
228   // Converts hits to digits of the current event
229   //
230   // Inputs file name
231   const char *alifile = "galice.root";   
232
233   // Create the run digitizer 
234   AliRunDigitizer* manager = new AliRunDigitizer(1, 1);
235   manager->SetInputStream(0, alifile);
236   manager->SetOutputFile("H2Dfile");
237
238   // Creates the VZERO digitizer 
239   AliVZERODigitizer* dig = new AliVZERODigitizer(manager);
240
241   // Creates the digits
242   dig->Exec("");
243
244 }
245 //_____________________________________________________________________________
246 void AliVZERO::Digits2Raw()
247 {
248    //
249    //  Converts digits of the current event to raw data
250    //
251   
252    AliVZERO *fVZERO = (AliVZERO*)gAlice->GetDetector("VZERO");
253    fLoader->LoadDigits();
254    TTree* digits = fLoader->TreeD();
255    if (!digits) {
256       Error("Digits2Raw", "no digits tree");
257       return;
258    }
259    TClonesArray * VZEROdigits = new TClonesArray("AliVZEROdigit",1000);
260    fVZERO->SetTreeAddress();            
261    digits->GetBranch("VZERODigit")->SetAddress(&VZEROdigits); 
262   
263    const char *fileName    = AliDAQ::DdlFileName("VZERO",0);
264    AliVZEROBuffer* buffer  = new AliVZEROBuffer(fileName);
265   
266    //  Verbose level
267    //  0: Silent
268    //  1: cout messages
269    //  2: txt files with digits 
270    //  BE CAREFUL, verbose level 2 MUST be used only for debugging and
271    //  it is highly suggested to use this mode only for debugging digits files
272    //  reasonably small, because otherwise the size of the txt files can reach
273    //  quickly several MB wasting time and disk space.
274   
275    ofstream ftxt;
276    buffer->SetVerbose(0);
277    Int_t verbose = buffer->GetVerbose();
278
279    // Get Trigger information first
280    // Read trigger inputs from trigger-detector object
281    AliDataLoader * dataLoader = fLoader->GetDigitsDataLoader();
282    if( !dataLoader->IsFileOpen() ) 
283         dataLoader->OpenFile( "READ" );
284    AliTriggerDetector* trgdet = (AliTriggerDetector*)dataLoader->GetDirectory()->Get( "Trigger" );
285    UInt_t triggerInfo = 0;
286    if(trgdet) {
287       triggerInfo = trgdet->GetMask() & 0xffff;
288    }
289    else {
290       AliError(Form("There is no trigger object for %s",fLoader->GetName()));
291    }
292
293    buffer->WriteTriggerInfo((UInt_t)triggerInfo); 
294    buffer->WriteTriggerScalers(); 
295    buffer->WriteBunchNumbers(); 
296   
297    // Now retrieve the channel information: charge+time and 
298    // dump it into ADC and Time arrays
299    // We assume here an ordered (by PMNumber) array of
300    // digits!!
301
302    Int_t nEntries = Int_t(digits->GetEntries());
303    Float_t ADC[64];
304    Int_t PMNumber[64];
305    Float_t Time[64];
306    Bool_t Integrator[64];
307   
308    for (Int_t i = 0; i < nEntries; i++) {
309      fVZERO->ResetDigits();
310      digits->GetEvent(i);
311      Int_t ndig = VZEROdigits->GetEntriesFast(); 
312    
313      if(ndig == 0) continue;
314      if(verbose == 2) {ftxt.open("VZEROdigits.txt",ios::app);}
315      for(Int_t k=0; k<ndig; k++){
316          AliVZEROdigit* fVZERODigit = (AliVZEROdigit*) VZEROdigits->At(k);
317          // Convert aliroot channel k into FEE channel iChannel before writing data
318          Int_t iChannel      = buffer->GetOnlineChannel(k);             
319          ADC[iChannel]       = fVZERODigit->ADC();
320          PMNumber[iChannel]  = fVZERODigit->PMNumber();
321          Time[iChannel]      = fVZERODigit->Time();
322          Integrator[iChannel]= fVZERODigit->Integrator(); 
323          if(verbose == 1) { cout <<"DDL: "<<fileName<< "\tdigit number: "<< k<<"\tPM number: "
324                             <<PMNumber[k]<<"\tADC: "<< ADC[k] << "\tTime: "<< Time[k] << endl;} 
325          if(verbose == 2) {
326               ftxt<<"DDL: "<<fileName<< "\tdigit number: "<< k<<"\tPM number: "
327                            <<PMNumber[k]<<"\tADC: "<< ADC[k] << "\tTime: "<< Time[k] << endl;}        
328 //       printf("DDL: %s, channel: %d, PM: %d, ADC: %f, Time: %f \n", 
329 //                  fileName,k,PMNumber[k],ADC[k],Time[k]); 
330      }        
331    if(verbose==2) ftxt.close();
332    }
333
334    // Now fill raw data
335           
336    for (Int_t  iCIU = 0; iCIU < 8; iCIU++) { 
337  
338    // decoding of one Channel Interface Unit numbered iCIU - there are 8 channels per CIU (and 8 CIUs) :
339   
340       for(Int_t iChannel_Offset = iCIU*8; iChannel_Offset < (iCIU*8)+8; iChannel_Offset=iChannel_Offset+4) { 
341          for(Int_t iChannel = iChannel_Offset; iChannel < iChannel_Offset+4; iChannel++) {
342              buffer->WriteChannel(iChannel, (Int_t) ADC[iChannel], Time[iChannel], Integrator[iChannel]);       
343          }
344          buffer->WriteBeamFlags(); 
345          buffer->WriteMBInfo(); 
346          buffer->WriteMBFlags();   
347          buffer->WriteBeamScalers(); 
348       } 
349 //      for(Int_t iChannel=0; iChannel < 8; iChannel++) {
350       for(Int_t iChannel=7; iChannel >= 0; iChannel--) {
351           buffer->WriteTiming(iChannel, (Int_t) ADC[iChannel], Time[iChannel]); 
352       }
353
354     // End of decoding of one CIU card
355     
356   } // end of decoding the eight CIUs
357      
358   delete buffer;
359   fLoader->UnloadDigits();  
360 }
361
362 //_____________________________________________________________________________
363 Bool_t AliVZERO::Raw2SDigits(AliRawReader* rawReader){
364   // Converts the VZERO raw data into digits
365   // The method is used for merging simulated and
366   // real data events
367   TStopwatch timer;
368   timer.Start();
369
370   if(!fLoader) {
371     AliError("no VZERO loader found");
372     return kFALSE; }
373
374   TTree* treeD  = fLoader->TreeD();
375   if(!treeD) {
376       fLoader->MakeTree("D");
377       treeD = fLoader->TreeD(); }
378         
379   AliVZEROdigit  digit;
380   AliVZEROdigit* pdigit = &digit;
381   const Int_t kBufferSize = 4000;
382    
383   treeD->Branch("VZERO", "AliVZEROdigit",  &pdigit, kBufferSize);
384
385   rawReader->Reset();
386   AliVZERORawStream* rawStream  = new AliVZERORawStream(rawReader);    
387      
388   if (!rawStream->Next()) return kFALSE; // No VZERO data found
389   
390   for(Int_t i=0; i<64; i++) {
391       new(pdigit) AliVZEROdigit(i, rawStream->GetADC(i), rawStream->GetTime(i)); 
392       treeD->Fill();
393   }
394  
395 // Checks if everything is OK by printing results 
396
397 //   for(int i=0;i<64;i++) {
398 //      printf("Channel %d : %f %f \n",i,rawStream->GetADC(i),rawStream->GetTime(i)); }
399 //   treeD->Print(); printf(" \n"); 
400         
401   fLoader->WriteDigits("OVERWRITE");
402   fLoader->UnloadDigits();      
403         
404   delete rawStream;
405
406   timer.Stop();
407   timer.Print();
408   return kTRUE;
409 }
410
411