]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCAL.cxx
Adding required tracker object for dHLT to work in AliReconstruction framework.
[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 #include <TH1F.h> 
34 #include <TF1.h> 
35 #include <TRandom.h> 
36
37 // --- Standard library ---
38
39 // --- AliRoot header files ---
40 #include "AliMagF.h"
41 #include "AliEMCAL.h"
42 #include "AliEMCALGetter.h"
43 #include "AliRun.h"
44 #include "AliEMCALSDigitizer.h"
45 #include "AliEMCALDigitizer.h"
46 #include "AliAltroBuffer.h"
47
48 ClassImp(AliEMCAL)
49 Double_t AliEMCAL::fgCapa        = 1.;        // 1pF 
50 Int_t    AliEMCAL::fgOrder       = 2 ;
51 Double_t AliEMCAL::fgTimeMax     = 2.56E-5 ;  // each sample is over 100 ns fTimeMax/fTimeBins
52 Double_t AliEMCAL::fgTimePeak    = 4.1E-6 ;   // 4 micro seconds
53 Double_t AliEMCAL::fgTimeTrigger = 100E-9 ;      // 100ns, just for a reference
54  
55 //____________________________________________________________________________
56 AliEMCAL::AliEMCAL():AliDetector()
57 {
58   // Default ctor 
59   fName = "EMCAL" ;
60 }
61
62 //____________________________________________________________________________
63 AliEMCAL::AliEMCAL(const char* name, const char* title): AliDetector(name,title)
64 {
65   //   ctor : title is used to identify the layout
66
67   fHighCharge        = 8.2 ;          // adjusted for a high gain range of 5.12 GeV (10 bits)
68   fHighGain          = 6.64 ; 
69   fHighLowGainFactor = 16. ;          // adjusted for a low gain range of 82 GeV (10 bits) 
70   fLowGainOffset     = 1 ;            // offset added to the module id to distinguish high and low gain data
71 }
72
73 //____________________________________________________________________________
74 AliEMCAL::~AliEMCAL()
75 {
76
77 }
78
79 //____________________________________________________________________________
80 void AliEMCAL::Copy(AliEMCAL & emcal)
81 {
82   TObject::Copy(emcal) ; 
83   emcal.fHighCharge        = fHighCharge ;
84   emcal.fHighGain          = fHighGain ; 
85   emcal.fHighLowGainFactor = fHighLowGainFactor ;  
86   emcal.fLowGainOffset     = fLowGainOffset;   
87 }
88
89 //____________________________________________________________________________
90 AliDigitizer* AliEMCAL::CreateDigitizer(AliRunDigitizer* manager) const
91 {
92   return new AliEMCALDigitizer(manager);
93 }
94
95 //____________________________________________________________________________
96 void AliEMCAL::CreateMaterials()
97 {
98   // Definitions of materials to build EMCAL and associated tracking media.
99   // media number in idtmed are 1599 to 1698.
100
101   // --- Air ---               
102   Float_t aAir[4]={12.0107,14.0067,15.9994,39.948};
103   Float_t zAir[4]={6.,7.,8.,18.};
104   Float_t wAir[4]={0.000124,0.755267,0.231781,0.012827};
105   Float_t dAir = 1.20479E-3;
106   AliMixture(0, "Air$", aAir, zAir, dAir, 4, wAir) ;
107
108   // --- Lead ---                                                                     
109   AliMaterial(1, "Pb$", 207.2, 82, 11.35, 0.56, 0., 0, 0) ;
110
111
112   // --- The polysterene scintillator (CH) ---
113   Float_t aP[2] = {12.011, 1.00794} ;
114   Float_t zP[2] = {6.0, 1.0} ;
115   Float_t wP[2] = {1.0, 1.0} ;
116   Float_t dP = 1.032 ;
117
118   AliMixture(2, "Polystyrene$", aP, zP, dP, -2, wP) ;
119
120   // --- Aluminium ---
121   AliMaterial(3, "Al$", 26.98, 13., 2.7, 8.9, 999., 0, 0) ;
122   // ---         Absorption length is ignored ^
123
124   // 25-aug-04 by PAI - see  PMD/AliPMDv0.cxx for STEEL definition
125   Float_t asteel[4] = { 55.847,51.9961,58.6934,28.0855 };
126   Float_t zsteel[4] = { 26.,24.,28.,14. };
127   Float_t wsteel[4] = { .715,.18,.1,.005 };
128   AliMixture(4, "STAINLESS STEEL$", asteel, zsteel, 7.88, 4, wsteel);
129
130   // DEFINITION OF THE TRACKING MEDIA
131
132   // for EMCAL: idtmed[1599->1698] equivalent to fIdtmed[0->100]
133   Int_t * idtmed = fIdtmed->GetArray() - 1599 ; 
134   Int_t   isxfld = gAlice->Field()->Integ() ;
135   Float_t sxmgmx = gAlice->Field()->Max() ;
136
137   // Air                                                                         -> idtmed[1599]
138  AliMedium(0, "Air          $", 0, 0,
139              isxfld, sxmgmx, 10.0, 1.0, 0.1, 0.1, 10.0, 0, 0) ;
140
141   // The Lead                                                                      -> idtmed[1600]
142  
143   AliMedium(1, "Lead      $", 1, 0,
144              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
145
146  // The scintillator of the CPV made of Polystyrene scintillator                   -> idtmed[1601]
147   AliMedium(2, "CPV scint.   $", 2, 1,
148             isxfld, sxmgmx, 10.0, 0.001, 0.1, 0.001, 0.001, 0, 0) ;
149
150   // Various Aluminium parts made of Al                                            -> idtmed[1602]
151   AliMedium(3, "Al parts     $", 3, 0,
152              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.001, 0.001, 0, 0) ;
153
154   // 25-aug-04 by PAI : see  PMD/AliPMDv0.cxx for STEEL definition                 -> idtmed[1603]
155   AliMedium(4, "S steel$", 4, 0, 
156              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.001, 0.001, 0, 0) ;
157
158 // --- Set decent energy thresholds for gamma and electron tracking
159
160   // Tracking threshold for photons and electrons in Lead 
161   gMC->Gstpar(idtmed[1600],"CUTGAM",0.00008) ;
162   gMC->Gstpar(idtmed[1600],"CUTELE",0.001) ;
163   gMC->Gstpar(idtmed[1600],"BCUTE",0.0001) ;
164
165   // --- Generate explicitly delta rays in Lead ---
166   gMC->Gstpar(idtmed[1600], "LOSS",3.) ;
167   gMC->Gstpar(idtmed[1600], "DRAY",1.) ;
168   gMC->Gstpar(idtmed[1600], "DCUTE",0.00001) ;
169   gMC->Gstpar(idtmed[1600], "DCUTM",0.00001) ;
170
171 // --- in aluminium parts ---
172   gMC->Gstpar(idtmed[1602], "LOSS",3.) ;
173   gMC->Gstpar(idtmed[1602], "DRAY",1.) ;
174   gMC->Gstpar(idtmed[1602], "DCUTE",0.00001) ;
175   gMC->Gstpar(idtmed[1602], "DCUTM",0.00001) ;
176
177 // --- and finally thresholds for photons and electrons in the scintillator ---
178   gMC->Gstpar(idtmed[1601],"CUTGAM",0.00008) ;
179   gMC->Gstpar(idtmed[1601],"CUTELE",0.001) ;
180   gMC->Gstpar(idtmed[1601],"BCUTE",0.0001) ;
181
182   // the same parameters as for Pb - 10-sep-04
183   gMC->Gstpar(idtmed[1603],"CUTGAM",0.00008) ;
184   gMC->Gstpar(idtmed[1603],"CUTELE",0.001) ;
185   gMC->Gstpar(idtmed[1603],"BCUTE", 0.0001); 
186   // --- Generate explicitly delta rays in Lead ---
187   gMC->Gstpar(idtmed[1603], "LOSS",3.);
188   gMC->Gstpar(idtmed[1603], "DRAY",1.);
189   gMC->Gstpar(idtmed[1603], "DCUTE",0.00001);
190   gMC->Gstpar(idtmed[1603], "DCUTM",0.00001);
191
192   //set constants for Birk's Law implentation
193   fBirkC0 =  1;
194   fBirkC1 =  0.013/dP;
195   fBirkC2 =  9.6e-6/(dP * dP);
196
197 }
198       
199 //____________________________________________________________________________
200 void AliEMCAL::Digits2Raw()
201 {
202   // convert digits of the current event to raw data
203   AliEMCALLoader * loader = dynamic_cast<AliEMCALLoader*>(fLoader) ; 
204
205   // get the digits
206   loader->LoadDigits();
207   TClonesArray* digits = loader->Digits() ;
208
209   if (!digits) {
210     Error("Digits2Raw", "no digits found !");
211     return;
212   }
213
214   // get the digitizer 
215   loader->LoadDigitizer();
216   AliEMCALDigitizer * digitizer = dynamic_cast<AliEMCALDigitizer *>(loader->Digitizer())  ; 
217   
218   // get the geometry
219   AliEMCALGeometry* geom = GetGeometry();
220   if (!geom) {
221     Error("Digits2Raw", "no geometry found !");
222     return;
223   }
224
225   // some digitization constants
226   const Int_t    kDDLOffset = 0x800;
227   const Int_t    kThreshold = 1;
228   const Int_t    kChannelsperDDL = 897 ; 
229   AliAltroBuffer* buffer = NULL;
230   Int_t prevDDL = -1;
231   Int_t adcValuesLow[fkTimeBins];
232   Int_t adcValuesHigh[fkTimeBins];
233   
234   // loop over digits (assume ordered digits)
235   for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
236     AliEMCALDigit* digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit)) ;
237     if (digit->GetAmp() < kThreshold) 
238       continue;
239     Int_t iDDL = digit->GetId() / kChannelsperDDL ;
240     // for each DDL id is numbered from 1 to  kChannelsperDDL -1 
241     Int_t idDDL = digit->GetId() - iDDL * ( kChannelsperDDL - 1 ) ;  
242     // new DDL
243     if (iDDL != prevDDL) {
244       // write real header and close previous file
245       if (buffer) {
246         buffer->Flush();
247         buffer->WriteDataHeader(kFALSE, kFALSE);
248         delete buffer;
249       }
250
251       // open new file and write dummy header
252       TString fileName("EMCAL_") ;
253       fileName += (iDDL + kDDLOffset) ; 
254       fileName += ".ddl" ; 
255       buffer = new AliAltroBuffer(fileName.Data(), 1);
256       buffer->WriteDataHeader(kTRUE, kFALSE);  //Dummy;
257
258       prevDDL = iDDL;
259     }
260
261     // out of time range signal (?)
262     if (digit->GetTimeR() > GetRawFormatTimeMax() ) {
263       buffer->FillBuffer(digit->GetAmp());
264       buffer->FillBuffer(GetRawFormatTimeBins() );  // time bin
265       buffer->FillBuffer(3);          // bunch length
266       buffer->WriteTrailer(3, idDDL, 0, 0);  // trailer
267
268       // calculate the time response function
269     } else {
270       Double_t energy = 0 ;  
271       energy = digit->GetAmp() * digitizer->GetECAchannel() + digitizer->GetECApedestal() ; 
272       
273       Bool_t lowgain = RawSampledResponse(digit->GetTimeR(), energy, adcValuesHigh, adcValuesLow) ; 
274       
275       if (lowgain) 
276         buffer->WriteChannel(iDDL, 0, fLowGainOffset, 
277                              GetRawFormatTimeBins(), adcValuesLow, kThreshold);
278       else 
279         buffer->WriteChannel(iDDL, 0, 0, 
280                              GetRawFormatTimeBins(), adcValuesHigh, kThreshold);
281       
282     }
283   }
284   
285   // write real header and close last file
286   if (buffer) {
287     buffer->Flush();
288     buffer->WriteDataHeader(kFALSE, kFALSE);
289     delete buffer;
290   }
291
292   loader->UnloadDigits();
293 }
294
295 //____________________________________________________________________________
296 void AliEMCAL::Hits2SDigits()  
297
298 // create summable digits
299
300   AliEMCALSDigitizer emcalDigitizer(fLoader->GetRunLoader()->GetFileName().Data()) ;
301   emcalDigitizer.SetEventRange(0, -1) ; // do all the events
302   emcalDigitizer.ExecuteTask() ;
303 }
304
305 //____________________________________________________________________________
306 AliLoader* AliEMCAL::MakeLoader(const char* topfoldername)
307 {
308 //different behaviour than standard (singleton getter)
309 // --> to be discussed and made eventually coherent
310  fLoader = new AliEMCALLoader(GetName(),topfoldername);
311  return fLoader;
312 }
313
314 //__________________________________________________________________
315 Double_t AliEMCAL::RawResponseFunction(Double_t *x, Double_t *par)
316 {
317   // Shape of the electronics raw reponse:
318   // It is a semi-gaussian, 2nd order Gamma function of the general form
319   // v(t) = n**n * Q * A**n / C *(t/tp)**n * exp(-n * t/tp) with 
320   // tp : peaking time par[0]
321   // n  : order of the function
322   // C  : integrating capacitor in the preamplifier
323   // A  : open loop gain of the preamplifier
324   // Q  : the total APD charge to be measured Q = C * energy
325   
326   Double_t signal ;
327   Double_t xx = x[0] - ( fgTimeTrigger + par[3] ) ; 
328   
329   if (xx < 0 || xx > fgTimeMax) 
330     signal = 0. ;  
331   else { 
332     Double_t fac = par[0] * TMath::Power(fgOrder, fgOrder) * TMath::Power(par[1], fgOrder) / fgCapa ; 
333     signal = fac * par[2] * TMath::Power(xx / fgTimePeak, fgOrder) * TMath::Exp(-fgOrder * (xx / fgTimePeak)) ; 
334   }
335   return signal ;  
336 }
337
338 //__________________________________________________________________
339 Double_t AliEMCAL::RawResponseFunctionMax(Double_t charge, Double_t gain) 
340 {
341   return ( charge * TMath::Power(fgOrder, fgOrder) * TMath::Power(gain, fgOrder) 
342      / ( fgCapa * TMath::Exp(fgOrder) ) );  
343
344 }
345 //__________________________________________________________________
346 Bool_t AliEMCAL::RawSampledResponse(
347 const Double_t dtime, const Double_t damp, Int_t * adcH, Int_t * adcL) const 
348 {
349   // for a start time dtime and an amplitude damp given by digit, 
350   // calculates the raw sampled response AliEMCAL::RawResponseFunction
351
352   const Int_t kRawSignalOverflow = 0x3FF ; 
353   Bool_t lowGain = kFALSE ; 
354
355   TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 4);
356
357   for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
358     signalF.SetParameter(0, GetRawFormatHighCharge() ) ; 
359     signalF.SetParameter(1, GetRawFormatHighGain() ) ; 
360     signalF.SetParameter(2, damp) ; 
361     signalF.SetParameter(3, dtime) ; 
362     Double_t time = iTime * GetRawFormatTimeMax() / GetRawFormatTimeBins() ;
363     Double_t signal = signalF.Eval(time) ;     
364     if ( static_cast<Int_t>(signal+0.5) > kRawSignalOverflow ){  // larger than 10 bits 
365       signal = kRawSignalOverflow ;
366       lowGain = kTRUE ; 
367     }
368     adcH[iTime] =  static_cast<Int_t>(signal + 0.5) ;
369
370     signalF.SetParameter(0, GetRawFormatLowCharge() ) ;     
371     signalF.SetParameter(1, GetRawFormatLowGain() ) ; 
372     signal = signalF.Eval(time) ;  
373     if ( static_cast<Int_t>(signal+0.5) > kRawSignalOverflow)  // larger than 10 bits 
374       signal = kRawSignalOverflow ;
375     adcL[iTime] = static_cast<Int_t>(0.5 + signal ) ; 
376
377   }
378   return lowGain ; 
379 }
380
381 //____________________________________________________________________________
382 void AliEMCAL::SetTreeAddress()
383
384   // Linking Hits in Tree to Hits array
385   TBranch *branch;
386   //  char branchname[20];
387   //  sprintf(branchname,"%s",GetName());
388   // Branch address for hit tree
389   TTree *treeH = TreeH();
390   if (treeH) {
391     //    treeH->Print();
392     branch = treeH->GetBranch(GetName());
393     if (branch) { 
394         if (fHits == 0x0) 
395         fHits= new TClonesArray("AliEMCALHit",1000);
396         branch->SetAddress(&fHits);
397     } else {
398       Warning("SetTreeAddress","<%s> Failed",GetName());
399     }
400   } else {
401     //    Warning("SetTreeAddress"," no treeH ");
402   }
403 }