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