]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCAL.cxx
Introducing the raw data format
[u/mrichter/AliRoot.git] / EMCAL / AliEMCAL.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 // Base Class for EMCAL description:
20 // This class contains material definitions    
21 // for the EMCAL - It does not place the detector in Alice
22 //*-- Author: Yves Schutz (SUBATECH) 
23 //
24 //*-- Additional Contributions: Sahal Yacoob (LBNL/UCT)
25 //
26 //////////////////////////////////////////////////////////////////////////////
27
28 // --- ROOT system ---
29 class TFile;
30 #include <TFolder.h> 
31 #include <TTree.h>
32 #include <TVirtualMC.h> 
33
34 // --- Standard library ---
35
36 // --- AliRoot header files ---
37 #include "AliMagF.h"
38 #include "AliEMCAL.h"
39 #include "AliEMCALGetter.h"
40 #include "AliRun.h"
41 #include "AliEMCALSDigitizer.h"
42 #include "AliEMCALDigitizer.h"
43 #include "AliAltroBuffer.h"
44
45 ClassImp(AliEMCAL)
46 //____________________________________________________________________________
47 AliEMCAL::AliEMCAL():AliDetector()
48 {
49   // Default ctor 
50   fName="EMCAL";
51   //  fGeom = 0 ; 
52   fRan    = new TRandom(123456) ; 
53 }
54
55 //____________________________________________________________________________
56 AliEMCAL::AliEMCAL(const char* name, const char* title): AliDetector(name,title)
57 {
58   //   ctor : title is used to identify the layout
59   // fGeom = 0;
60   fRan    = new TRandom(123456) ; 
61 }
62
63 //____________________________________________________________________________
64 AliEMCAL::~AliEMCAL()
65 {
66
67 }
68
69 //____________________________________________________________________________
70 void AliEMCAL::Copy(AliEMCAL & emcal)
71 {
72   TObject::Copy(emcal) ; 
73 }
74
75 //____________________________________________________________________________
76 AliDigitizer* AliEMCAL::CreateDigitizer(AliRunDigitizer* manager) const
77 {
78   return new AliEMCALDigitizer(manager);
79 }
80
81 //____________________________________________________________________________
82 void AliEMCAL::CreateMaterials()
83 {
84   // Definitions of materials to build EMCAL and associated tracking media.
85   // media number in idtmed are 1599 to 1698.
86
87   // --- Air ---               
88   Float_t aAir[4]={12.0107,14.0067,15.9994,39.948};
89   Float_t zAir[4]={6.,7.,8.,18.};
90   Float_t wAir[4]={0.000124,0.755267,0.231781,0.012827};
91   Float_t dAir = 1.20479E-3;
92   AliMixture(0, "Air$", aAir, zAir, dAir, 4, wAir) ;
93
94   // --- Lead ---                                                                     
95   AliMaterial(1, "Pb$", 207.2, 82, 11.35, 0.56, 0., 0, 0) ;
96
97
98   // --- The polysterene scintillator (CH) ---
99   Float_t aP[2] = {12.011, 1.00794} ;
100   Float_t zP[2] = {6.0, 1.0} ;
101   Float_t wP[2] = {1.0, 1.0} ;
102   Float_t dP = 1.032 ;
103
104   AliMixture(2, "Polystyrene$", aP, zP, dP, -2, wP) ;
105
106   // --- Aluminium ---
107   AliMaterial(3, "Al$", 26.98, 13., 2.7, 8.9, 999., 0, 0) ;
108   // ---         Absorption length is ignored ^
109
110   // DEFINITION OF THE TRACKING MEDIA
111
112   // for EMCAL: idtmed[1599->1698] equivalent to fIdtmed[0->100]
113   Int_t * idtmed = fIdtmed->GetArray() - 1599 ; 
114   Int_t   isxfld = gAlice->Field()->Integ() ;
115   Float_t sxmgmx = gAlice->Field()->Max() ;
116
117   // Air                                                                         -> idtmed[1599]
118  AliMedium(0, "Air          $", 0, 0,
119              isxfld, sxmgmx, 10.0, 1.0, 0.1, 0.1, 10.0, 0, 0) ;
120
121   // The Lead                                                                      -> idtmed[1600]
122  
123   AliMedium(1, "Lead      $", 1, 0,
124              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
125
126  // The scintillator of the CPV made of Polystyrene scintillator                   -> idtmed[1601]
127   AliMedium(2, "CPV scint.   $", 2, 1,
128             isxfld, sxmgmx, 10.0, 0.001, 0.1, 0.001, 0.001, 0, 0) ;
129
130   // Various Aluminium parts made of Al                                            -> idtmed[1602]
131   AliMedium(3, "Al parts     $", 3, 0,
132              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.001, 0.001, 0, 0) ;
133
134
135 // --- Set decent energy thresholds for gamma and electron tracking
136
137   // Tracking threshold for photons and electrons in Lead 
138   gMC->Gstpar(idtmed[1600],"CUTGAM",0.00008) ;
139   gMC->Gstpar(idtmed[1600],"CUTELE",0.001) ;
140   gMC->Gstpar(idtmed[1600],"BCUTE",0.0001) ;
141
142   // --- Generate explicitly delta rays in Lead ---
143   gMC->Gstpar(idtmed[1600], "LOSS",3.) ;
144   gMC->Gstpar(idtmed[1600], "DRAY",1.) ;
145   gMC->Gstpar(idtmed[1600], "DCUTE",0.00001) ;
146   gMC->Gstpar(idtmed[1600], "DCUTM",0.00001) ;
147
148 // --- in aluminium parts ---
149   gMC->Gstpar(idtmed[1602], "LOSS",3.) ;
150   gMC->Gstpar(idtmed[1602], "DRAY",1.) ;
151   gMC->Gstpar(idtmed[1602], "DCUTE",0.00001) ;
152   gMC->Gstpar(idtmed[1602], "DCUTM",0.00001) ;
153
154 // --- and finally thresholds for photons and electrons in the scintillator ---
155   gMC->Gstpar(idtmed[1601],"CUTGAM",0.00008) ;
156   gMC->Gstpar(idtmed[1601],"CUTELE",0.001) ;
157   gMC->Gstpar(idtmed[1601],"BCUTE",0.0001) ;
158
159   //set constants for Birk's Law implentation
160   fBirkC0 =  1;
161   fBirkC1 =  0.013/dP;
162   fBirkC2 =  9.6e-6/(dP * dP);
163
164
165 }
166       
167 //____________________________________________________________________________
168 void AliEMCAL::Digits2Raw()
169 {
170 // convert digits of the current event to raw data
171
172   // get the digits
173   AliEMCALGetter * gime = AliEMCALGetter::Instance(AliRunLoader::GetGAliceName()) ; 
174   if (!gime) {
175     Error("Digits2Raw", "EMCAL Getter not instantiated") ;
176     return ; 
177   }
178   gime->Event(gime->EventNumber(), "D") ; 
179   TClonesArray* digits = gime->Digits() ;
180
181   if (!digits) {
182     Error("Digits2Raw", "no digits found !");
183     return;
184   }
185
186   // get the geometry
187   AliEMCALGeometry* geom = gime->EMCALGeometry();
188   if (!geom) {
189     Error("Digits2Raw", "no geometry found !");
190     return;
191   }
192
193   // some digitization constants
194   const Int_t    kDDLOffset = 0x800;
195   const Double_t kTimeMax = 1.28E-5;
196   const Int_t    kTimeBins = 256;
197   const Double_t kTimePeak = 2.0E-6;
198   const Double_t kTimeRes = 1.5E-6;
199   const Int_t    kThreshold = 3;
200   const Int_t    kHighGainFactor = 40;
201   const Int_t    kHighGainOffset = 0x200;
202   // PHOS has 4 DDL per module; I assume therefore that kChannelsperDDL=896+1 EMCAL channel go to one DDL
203   const Int_t    kChannelsperDDL = 897 ; 
204
205   AliAltroBuffer* buffer = NULL;
206   Int_t prevDDL = -1;
207   Int_t adcValuesLow[kTimeBins];
208   Int_t adcValuesHigh[kTimeBins];
209
210   // loop over digits (assume ordered digits)
211   for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
212     AliEMCALDigit* digit = gime->Digit(iDigit);
213     if (digit->GetAmp() < kThreshold) 
214       continue;
215     Int_t iDDL = digit->GetId() / kChannelsperDDL ;
216     // for each DDL id is numbered from 1 to  kChannelsperDDL -1 
217     Int_t idDDL = digit->GetId() - iDDL * ( kChannelsperDDL - 1 ) ;  
218     // new DDL
219     if (iDDL != prevDDL) {
220       // write real header and close previous file
221       if (buffer) {
222         buffer->Flush();
223         buffer->WriteDataHeader(kFALSE, kFALSE);
224         delete buffer;
225       }
226
227       // open new file and write dummy header
228       TString fileName("EMCAL_") ;
229       fileName += (iDDL + kDDLOffset) ; 
230       fileName += ".ddl" ; 
231       buffer = new AliAltroBuffer(fileName.Data(), 1);
232       buffer->WriteDataHeader(kTRUE, kFALSE);  //Dummy;
233
234       prevDDL = iDDL;
235     }
236
237     // out of time range signal (?)
238     if (digit->GetTime() > kTimeMax) {
239       buffer->FillBuffer(digit->GetAmp());
240       buffer->FillBuffer(kTimeBins);  // time bin
241       buffer->FillBuffer(3);          // bunch length
242       buffer->WriteTrailer(3, idDDL, 0, 0);  // trailer
243       
244     // simulate linear rise and gaussian decay of signal
245     } else {
246       Bool_t highGain = kFALSE;
247
248       // fill time bin values :
249       // 1. the signal starts at the time given by the digit
250       // 2. the rise is linear and the maximum is reached kTimePeak after start
251       // 3. the decay is gaussian with a sigma of kTimeRes
252       // 4. the signal is binned into kTimeBins bins 
253       for (Int_t iTime = 0; iTime < kTimeBins; iTime++) {
254         Double_t time = iTime * kTimeMax/kTimeBins;
255         Int_t signal = 0;
256         if (time < digit->GetTime() + kTimePeak) {      // signal is rising
257           signal = static_cast<Int_t>(fRan->Rndm() + digit->GetAmp() * 
258                          (time - digit->GetTime()) / kTimePeak);
259         } else {                                        // signal is decaying
260           signal = static_cast<Int_t>(fRan->Rndm() + digit->GetAmp() * 
261                  TMath::Gaus(time, digit->GetTime() + kTimePeak, kTimeRes));
262         }
263         if (signal < 0) 
264           signal = 0;
265         adcValuesLow[iTime] = signal;
266         if (signal > 0x3FF) // larger than 10 bits 
267           adcValuesLow[iTime] = 0x3FF;
268         adcValuesHigh[iTime] = signal / kHighGainFactor;
269         if (adcValuesHigh[iTime] > 0) 
270           highGain = kTRUE;
271       }
272
273       // write low and eventually high gain channel
274       buffer->WriteChannel(idDDL, 0, 0, 
275                            kTimeBins, adcValuesLow, kThreshold);
276       if (highGain) {
277         buffer->WriteChannel(idDDL, 0, kHighGainOffset, 
278                              kTimeBins, adcValuesHigh, 1);
279       }
280     }
281   }
282
283   // write real header and close last file
284   if (buffer) {
285     buffer->Flush();
286     buffer->WriteDataHeader(kFALSE, kFALSE);
287     delete buffer;
288   }
289
290   gime->EmcalLoader()->UnloadDigits();
291 }
292
293 //____________________________________________________________________________
294 void AliEMCAL::Hits2SDigits()  
295
296 // create summable digits
297
298   AliEMCALSDigitizer* emcalDigitizer = 
299     new AliEMCALSDigitizer(fLoader->GetRunLoader()->GetFileName().Data()) ;
300   emcalDigitizer->SetEventRange(0, -1) ; // do all the events
301   emcalDigitizer->ExecuteTask() ;
302 }
303
304 //____________________________________________________________________________
305 AliLoader* AliEMCAL::MakeLoader(const char* topfoldername)
306 {
307 //different behaviour than standard (singleton getter)
308 // --> to be discussed and made eventually coherent
309  fLoader = new AliEMCALLoader(GetName(),topfoldername);
310  return fLoader;
311 }
312
313 //____________________________________________________________________________
314 void AliEMCAL::SetTreeAddress()
315
316   // Linking Hits in Tree to Hits array
317   TBranch *branch;
318   char branchname[20];
319   sprintf(branchname,"%s",GetName());
320   
321   // Branch address for hit tree
322   TTree *treeH = TreeH();
323   if (treeH) {
324     branch = treeH->GetBranch(branchname);
325     if (branch) 
326       { 
327         if (fHits == 0x0) 
328           fHits= new TClonesArray("AliEMCALHit",1000);
329         branch->SetAddress(&fHits);
330       }
331     else
332       {
333         Warning("SetTreeAddress","<%s> Failed",GetName());
334       }
335   }
336 }
337
338
339
340
341