]> git.uio.no Git - u/mrichter/AliRoot.git/blob - VZERO/AliVZERO.cxx
Major fixes to the raw-data simulation. The ADC samples and signal width are added...
[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 of the current event
231   //
232   // Inputs file name
233   const char *alifile = "galice.root";   
234
235   // Create the run digitizer 
236   AliRunDigitizer* manager = new AliRunDigitizer(1, 1);
237   manager->SetInputStream(0, alifile);
238   manager->SetOutputFile("H2Dfile");
239
240   // Creates the VZERO digitizer 
241   AliVZERODigitizer* dig = new AliVZERODigitizer(manager);
242
243   // Creates the digits
244   dig->Exec("");
245
246 }
247 //_____________________________________________________________________________
248 void AliVZERO::Digits2Raw()
249 {
250    //
251    //  Converts digits of the current event to raw data
252    //
253   
254    AliVZERO *fVZERO = (AliVZERO*)gAlice->GetDetector("VZERO");
255    fLoader->LoadDigits();
256    TTree* digits = fLoader->TreeD();
257    if (!digits) {
258       Error("Digits2Raw", "no digits tree");
259       return;
260    }
261    TClonesArray * VZEROdigits = new TClonesArray("AliVZEROdigit",1000);
262    fVZERO->SetTreeAddress();            
263    digits->GetBranch("VZERODigit")->SetAddress(&VZEROdigits); 
264   
265    const char *fileName    = AliDAQ::DdlFileName("VZERO",0);
266    AliVZEROBuffer* buffer  = new AliVZEROBuffer(fileName);
267   
268    // Get Trigger information first
269    // Read trigger inputs from trigger-detector object
270    AliDataLoader * dataLoader = fLoader->GetDigitsDataLoader();
271    if( !dataLoader->IsFileOpen() ) 
272         dataLoader->OpenFile( "READ" );
273    AliTriggerDetector* trgdet = (AliTriggerDetector*)dataLoader->GetDirectory()->Get( "Trigger" );
274    UInt_t triggerInfo = 0;
275    if(trgdet) {
276       triggerInfo = trgdet->GetMask() & 0xffff;
277    }
278    else {
279       AliError(Form("There is no trigger object for %s",fLoader->GetName()));
280    }
281
282    buffer->WriteTriggerInfo((UInt_t)triggerInfo); 
283    buffer->WriteTriggerScalers(); 
284    buffer->WriteBunchNumbers(); 
285   
286    Int_t aBBflagsV0A = 0;
287    Int_t aBBflagsV0C = 0;
288    Int_t aBGflagsV0A = 0;
289    Int_t aBGflagsV0C = 0;
290
291    if (digits->GetUserInfo()->FindObject("BBflagsV0A")) {
292      aBBflagsV0A = ((TParameter<int>*)digits->GetUserInfo()->FindObject("BBflagsV0A"))->GetVal();
293    }
294    else
295      AliWarning("V0A beam-beam flags not found in digits tree UserInfo! The flags will not be written to the raw-data stream!");
296
297    if (digits->GetUserInfo()->FindObject("BBflagsV0C")) {
298      aBBflagsV0C = ((TParameter<int>*)digits->GetUserInfo()->FindObject("BBflagsV0C"))->GetVal();
299    }
300    else
301      AliWarning("V0C 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("BGflagsV0A")) {
304      aBGflagsV0A = ((TParameter<int>*)digits->GetUserInfo()->FindObject("BGflagsV0A"))->GetVal();
305    }
306    else
307      AliWarning("V0A beam-gas flags not found in digits tree UserInfo! The flags will not be written to the raw-data stream!");
308
309    if (digits->GetUserInfo()->FindObject("BGflagsV0C")) {
310      aBGflagsV0C = ((TParameter<int>*)digits->GetUserInfo()->FindObject("BGflagsV0C"))->GetVal();
311    }
312    else
313      AliWarning("V0C beam-gas flags not found in digits tree UserInfo! The flags will not be written to the raw-data stream!");
314
315    // Now retrieve the channel information: charge smaples + time and 
316    // dump it into ADC and Time arrays
317    Int_t nEntries = Int_t(digits->GetEntries());
318    Short_t aADC[64][AliVZEROdigit::kNClocks];
319    Float_t aTime[64];
320    Float_t aWidth[64];
321    Bool_t  aIntegrator[64];
322    Bool_t  aBBflag[64];
323    Bool_t  aBGflag[64];
324   
325    for (Int_t i = 0; i < nEntries; i++) {
326      fVZERO->ResetDigits();
327      digits->GetEvent(i);
328      Int_t ndig = VZEROdigits->GetEntriesFast(); 
329    
330      if(ndig == 0) continue;
331      for(Int_t k=0; k<ndig; k++){
332          AliVZEROdigit* fVZERODigit = (AliVZEROdigit*) VZEROdigits->At(k);
333          // Convert aliroot channel k into FEE channel iChannel before writing data
334          Int_t iChannel       = AliVZEROCalibData::GetBoardNumber(fVZERODigit->PMNumber()) * 8 +
335            AliVZEROCalibData::GetFEEChannelNumber(fVZERODigit->PMNumber());
336          for(Int_t iClock = 0; iClock < AliVZEROdigit::kNClocks; ++iClock) aADC[iChannel][iClock] = fVZERODigit->ChargeADC(iClock);
337          aTime[iChannel]      = fVZERODigit->Time();
338          aWidth[iChannel]     = fVZERODigit->Width();
339          aIntegrator[iChannel]= fVZERODigit->Integrator();
340          if(fVZERODigit->PMNumber() < 32) {
341            aBBflag[iChannel] = (aBBflagsV0C >> fVZERODigit->PMNumber()) & 0x1;
342            aBGflag[iChannel] = (aBGflagsV0C >> fVZERODigit->PMNumber()) & 0x1;
343          }
344          else {
345            aBBflag[iChannel] = (aBBflagsV0A >> (fVZERODigit->PMNumber()-32)) & 0x1;
346            aBGflag[iChannel] = (aBGflagsV0A >> (fVZERODigit->PMNumber()-32)) & 0x1;
347          }
348          AliDebug(1,Form("DDL: %s\tdigit number: %d\tPM number: %d\tADC: %f\tTime: %f",
349                          fileName,k,fVZERODigit->PMNumber(),aADC[k],aTime[k])); 
350      }        
351    }
352
353    // Now fill raw data
354           
355    for (Int_t  iCIU = 0; iCIU < 8; iCIU++) { 
356  
357    // decoding of one Channel Interface Unit numbered iCIU - there are 8 channels per CIU (and 8 CIUs) :
358   
359       for(Int_t iChannel_Offset = iCIU*8; iChannel_Offset < (iCIU*8)+8; iChannel_Offset=iChannel_Offset+4) { 
360          for(Int_t iChannel = iChannel_Offset; iChannel < iChannel_Offset+4; iChannel++) {
361              buffer->WriteChannel(iChannel, aADC[iChannel], aIntegrator[iChannel]);       
362          }
363          buffer->WriteBeamFlags(&aBBflag[iChannel_Offset],&aBGflag[iChannel_Offset]); 
364          buffer->WriteMBInfo(); 
365          buffer->WriteMBFlags();   
366          buffer->WriteBeamScalers(); 
367       } 
368
369       for(Int_t iChannel = iCIU*8 + 7; iChannel >= iCIU*8; iChannel--) {
370           buffer->WriteTiming(aTime[iChannel], aWidth[iChannel]); 
371       }
372
373     // End of decoding of one CIU card
374     
375   } // end of decoding the eight CIUs
376      
377   delete buffer;
378   fLoader->UnloadDigits();  
379 }
380
381 //_____________________________________________________________________________
382 Bool_t AliVZERO::Raw2SDigits(AliRawReader* rawReader){
383   // Converts the VZERO raw data into digits
384   // The method is used for merging simulated and
385   // real data events
386   TStopwatch timer;
387   timer.Start();
388
389   if(!fLoader) {
390     AliError("no VZERO loader found");
391     return kFALSE; }
392
393   TTree* treeD  = fLoader->TreeD();
394   if(!treeD) {
395       fLoader->MakeTree("D");
396       treeD = fLoader->TreeD(); }
397         
398   AliVZEROdigit  digit;
399   AliVZEROdigit* pdigit = &digit;
400   const Int_t kBufferSize = 4000;
401    
402   treeD->Branch("VZERO", "AliVZEROdigit",  &pdigit, kBufferSize);
403
404   rawReader->Reset();
405   AliVZERORawStream* rawStream  = new AliVZERORawStream(rawReader);    
406      
407   if (!rawStream->Next()) return kFALSE; // No VZERO data found
408   
409   for(Int_t i=0; i<64; i++) {
410       new(pdigit) AliVZEROdigit(i, rawStream->GetADC(i), rawStream->GetTime(i)); 
411       treeD->Fill();
412   }
413  
414 // Checks if everything is OK by printing results 
415
416 //   for(int i=0;i<64;i++) {
417 //      printf("Channel %d : %f %f \n",i,rawStream->GetADC(i),rawStream->GetTime(i)); }
418 //   treeD->Print(); printf(" \n"); 
419         
420   fLoader->WriteDigits("OVERWRITE");
421   fLoader->UnloadDigits();      
422         
423   delete rawStream;
424
425   timer.Stop();
426   timer.Print();
427   return kTRUE;
428 }
429
430