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