]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCAL.cxx
default geant cuts for EMCAL's material: cutele=cutgam=100KeV
[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, "Scintillator$", 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$", 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   Float_t cutgam=10.e-5; // 100 kev;
162   Float_t cutele=10.e-5; // 100 kev;
163   TString ntmp(GetTitle()); 
164   ntmp.ToUpper();
165   if(ntmp.Contains("10KEV")) {
166     cutele = cutgam = 1.e-5;
167   } else if(ntmp.Contains("50KEV")) {
168     cutele = cutgam = 5.e-5;
169   } else if(ntmp.Contains("100KEV")) {
170     cutele = cutgam = 1.e-4;
171   } else if(ntmp.Contains("200KEV")) {
172     cutele = cutgam = 2.e-4;
173   } else if(ntmp.Contains("500KEV")) {
174     cutele = cutgam = 5.e-4;
175   }
176
177   gMC->Gstpar(idtmed[1600],"CUTGAM", cutgam);
178   gMC->Gstpar(idtmed[1600],"CUTELE", cutele); // 1MEV -> 0.1MEV; 15-aug-05
179   gMC->Gstpar(idtmed[1600],"BCUTE",  cutgam);  // BCUTE and BCUTM start from GUTGUM
180   gMC->Gstpar(idtmed[1600],"BCUTM",  cutgam);  // BCUTE and BCUTM start from GUTGUM
181   // --- Generate explicitly delta rays in Lead ---
182   gMC->Gstpar(idtmed[1600], "LOSS",3.) ;
183   gMC->Gstpar(idtmed[1600], "DRAY",1.) ;
184   gMC->Gstpar(idtmed[1600], "DCUTE", cutele) ;
185   gMC->Gstpar(idtmed[1600], "DCUTM", cutele) ;
186
187 // --- in aluminium parts ---
188   gMC->Gstpar(idtmed[1602],"CUTGAM", cutgam) ;
189   gMC->Gstpar(idtmed[1602],"CUTELE", cutele) ;
190   gMC->Gstpar(idtmed[1602],"BCUTE",  cutgam);  // BCUTE and BCUTM start from GUTGUM
191   gMC->Gstpar(idtmed[1602],"BCUTM",  cutgam);  // BCUTE and BCUTM start from GUTGUM
192   gMC->Gstpar(idtmed[1602], "LOSS",3.) ;
193   gMC->Gstpar(idtmed[1602], "DRAY",1.) ;
194   gMC->Gstpar(idtmed[1602], "DCUTE", cutele) ;
195   gMC->Gstpar(idtmed[1602], "DCUTM", cutele) ;
196
197 // --- and finally thresholds for photons and electrons in the scintillator ---
198   gMC->Gstpar(idtmed[1601],"CUTGAM", cutgam) ;
199   gMC->Gstpar(idtmed[1601],"CUTELE", cutele) ;// 1MEV -> 0.1MEV; 15-aug-05
200   gMC->Gstpar(idtmed[1601],"BCUTE",  cutgam);  // BCUTE and BCUTM start from GUTGUM
201   gMC->Gstpar(idtmed[1601],"BCUTM",  cutgam);  // BCUTE and BCUTM start from GUTGUM
202   gMC->Gstpar(idtmed[1601], "LOSS",3.) ; // generate delta rays 
203   gMC->Gstpar(idtmed[1601], "DRAY",1.) ;
204   gMC->Gstpar(idtmed[1601], "DCUTE", cutele) ;
205   gMC->Gstpar(idtmed[1601], "DCUTM", cutele) ;
206
207   // S steel - 
208   gMC->Gstpar(idtmed[1603],"CUTGAM", cutgam);
209   gMC->Gstpar(idtmed[1603],"CUTELE", cutele);
210   gMC->Gstpar(idtmed[1603],"BCUTE",  cutgam);  // BCUTE and BCUTM start from GUTGUM
211   gMC->Gstpar(idtmed[1603],"BCUTM",  cutgam);  // BCUTE and BCUTM start from GUTGUM
212   // --- Generate explicitly delta rays 
213   gMC->Gstpar(idtmed[1603], "LOSS",3.);
214   gMC->Gstpar(idtmed[1603], "DRAY",1.);
215   gMC->Gstpar(idtmed[1603], "DCUTE", cutele) ;
216   gMC->Gstpar(idtmed[1603], "DCUTM", cutele) ;
217
218   //set constants for Birk's Law implentation
219   fBirkC0 =  1;
220   fBirkC1 =  0.013/dP;
221   fBirkC2 =  9.6e-6/(dP * dP);
222
223 }
224       
225 //____________________________________________________________________________
226 void AliEMCAL::Digits2Raw()
227 {
228   // convert digits of the current event to raw data
229   AliEMCALLoader * loader = dynamic_cast<AliEMCALLoader*>(fLoader) ; 
230
231   // get the digits
232   loader->LoadDigits();
233   TClonesArray* digits = loader->Digits() ;
234
235   if (!digits) {
236     Error("Digits2Raw", "no digits found !");
237     return;
238   }
239
240   // get the digitizer 
241   loader->LoadDigitizer();
242   AliEMCALDigitizer * digitizer = dynamic_cast<AliEMCALDigitizer *>(loader->Digitizer())  ; 
243   
244   // get the geometry
245   AliEMCALGeometry* geom = GetGeometry();
246   if (!geom) {
247     Error("Digits2Raw", "no geometry found !");
248     return;
249   }
250
251   // some digitization constants
252   const Int_t    kDDLOffset = 0x800;
253   const Int_t    kThreshold = 1;
254   const Int_t    kChannelsperDDL = 897 ; 
255   AliAltroBuffer* buffer = NULL;
256   Int_t prevDDL = -1;
257   Int_t adcValuesLow[fkTimeBins];
258   Int_t adcValuesHigh[fkTimeBins];
259   
260   // loop over digits (assume ordered digits)
261   for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
262     AliEMCALDigit* digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit)) ;
263     if (digit->GetAmp() < kThreshold) 
264       continue;
265     Int_t iDDL = digit->GetId() / kChannelsperDDL ;
266     // for each DDL id is numbered from 1 to  kChannelsperDDL -1 
267     Int_t idDDL = digit->GetId() - iDDL * ( kChannelsperDDL - 1 ) ;  
268     // new DDL
269     if (iDDL != prevDDL) {
270       // write real header and close previous file
271       if (buffer) {
272         buffer->Flush();
273         buffer->WriteDataHeader(kFALSE, kFALSE);
274         delete buffer;
275       }
276
277       // open new file and write dummy header
278       TString fileName("EMCAL_") ;
279       fileName += (iDDL + kDDLOffset) ; 
280       fileName += ".ddl" ; 
281       buffer = new AliAltroBuffer(fileName.Data(), 1);
282       buffer->WriteDataHeader(kTRUE, kFALSE);  //Dummy;
283
284       prevDDL = iDDL;
285     }
286
287     // out of time range signal (?)
288     if (digit->GetTimeR() > GetRawFormatTimeMax() ) {
289       buffer->FillBuffer(digit->GetAmp());
290       buffer->FillBuffer(GetRawFormatTimeBins() );  // time bin
291       buffer->FillBuffer(3);          // bunch length
292       buffer->WriteTrailer(3, idDDL, 0, 0);  // trailer
293
294       // calculate the time response function
295     } else {
296       Double_t energy = 0 ;  
297       energy = digit->GetAmp() * digitizer->GetECAchannel() + digitizer->GetECApedestal() ; 
298       
299       Bool_t lowgain = RawSampledResponse(digit->GetTimeR(), energy, adcValuesHigh, adcValuesLow) ; 
300       
301       if (lowgain) 
302         buffer->WriteChannel(iDDL, 0, fLowGainOffset, 
303                              GetRawFormatTimeBins(), adcValuesLow, kThreshold);
304       else 
305         buffer->WriteChannel(iDDL, 0, 0, 
306                              GetRawFormatTimeBins(), adcValuesHigh, kThreshold);
307       
308     }
309   }
310   
311   // write real header and close last file
312   if (buffer) {
313     buffer->Flush();
314     buffer->WriteDataHeader(kFALSE, kFALSE);
315     delete buffer;
316   }
317
318   loader->UnloadDigits();
319 }
320
321 //____________________________________________________________________________
322 void AliEMCAL::Hits2SDigits()  
323
324 // create summable digits
325
326   AliEMCALSDigitizer emcalDigitizer(fLoader->GetRunLoader()->GetFileName().Data()) ;
327   emcalDigitizer.SetEventRange(0, -1) ; // do all the events
328   emcalDigitizer.ExecuteTask() ;
329 }
330
331 //____________________________________________________________________________
332 AliLoader* AliEMCAL::MakeLoader(const char* topfoldername)
333 {
334 //different behaviour than standard (singleton getter)
335 // --> to be discussed and made eventually coherent
336  fLoader = new AliEMCALLoader(GetName(),topfoldername);
337  return fLoader;
338 }
339
340 //__________________________________________________________________
341 Double_t AliEMCAL::RawResponseFunction(Double_t *x, Double_t *par)
342 {
343   // Shape of the electronics raw reponse:
344   // It is a semi-gaussian, 2nd order Gamma function of the general form
345   // v(t) = n**n * Q * A**n / C *(t/tp)**n * exp(-n * t/tp) with 
346   // tp : peaking time par[0]
347   // n  : order of the function
348   // C  : integrating capacitor in the preamplifier
349   // A  : open loop gain of the preamplifier
350   // Q  : the total APD charge to be measured Q = C * energy
351   
352   Double_t signal ;
353   Double_t xx = x[0] - ( fgTimeTrigger + par[3] ) ; 
354   
355   if (xx < 0 || xx > fgTimeMax) 
356     signal = 0. ;  
357   else { 
358     Double_t fac = par[0] * TMath::Power(fgOrder, fgOrder) * TMath::Power(par[1], fgOrder) / fgCapa ; 
359     signal = fac * par[2] * TMath::Power(xx / fgTimePeak, fgOrder) * TMath::Exp(-fgOrder * (xx / fgTimePeak)) ; 
360   }
361   return signal ;  
362 }
363
364 //__________________________________________________________________
365 Double_t AliEMCAL::RawResponseFunctionMax(Double_t charge, Double_t gain) 
366 {
367   return ( charge * TMath::Power(fgOrder, fgOrder) * TMath::Power(gain, fgOrder) 
368      / ( fgCapa * TMath::Exp(fgOrder) ) );  
369
370 }
371 //__________________________________________________________________
372 Bool_t AliEMCAL::RawSampledResponse(
373 const Double_t dtime, const Double_t damp, Int_t * adcH, Int_t * adcL) const 
374 {
375   // for a start time dtime and an amplitude damp given by digit, 
376   // calculates the raw sampled response AliEMCAL::RawResponseFunction
377
378   const Int_t kRawSignalOverflow = 0x3FF ; 
379   Bool_t lowGain = kFALSE ; 
380
381   TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 4);
382
383   for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
384     signalF.SetParameter(0, GetRawFormatHighCharge() ) ; 
385     signalF.SetParameter(1, GetRawFormatHighGain() ) ; 
386     signalF.SetParameter(2, damp) ; 
387     signalF.SetParameter(3, dtime) ; 
388     Double_t time = iTime * GetRawFormatTimeMax() / GetRawFormatTimeBins() ;
389     Double_t signal = signalF.Eval(time) ;     
390     if ( static_cast<Int_t>(signal+0.5) > kRawSignalOverflow ){  // larger than 10 bits 
391       signal = kRawSignalOverflow ;
392       lowGain = kTRUE ; 
393     }
394     adcH[iTime] =  static_cast<Int_t>(signal + 0.5) ;
395
396     signalF.SetParameter(0, GetRawFormatLowCharge() ) ;     
397     signalF.SetParameter(1, GetRawFormatLowGain() ) ; 
398     signal = signalF.Eval(time) ;  
399     if ( static_cast<Int_t>(signal+0.5) > kRawSignalOverflow)  // larger than 10 bits 
400       signal = kRawSignalOverflow ;
401     adcL[iTime] = static_cast<Int_t>(0.5 + signal ) ; 
402
403   }
404   return lowGain ; 
405 }
406
407 //____________________________________________________________________________
408 void AliEMCAL::SetTreeAddress()
409
410   // Linking Hits in Tree to Hits array
411   TBranch *branch;
412   //  char branchname[20];
413   //  sprintf(branchname,"%s",GetName());
414   // Branch address for hit tree
415   TTree *treeH = TreeH();
416   if (treeH) {
417     //    treeH->Print();
418     branch = treeH->GetBranch(GetName());
419     if (branch) { 
420         if (fHits == 0x0) 
421         fHits= new TClonesArray("AliEMCALHit",1000);
422         branch->SetAddress(&fHits);
423     } else {
424       Warning("SetTreeAddress","<%s> Failed",GetName());
425     }
426   } else {
427     //    Warning("SetTreeAddress"," no treeH ");
428   }
429 }