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