]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOS.cxx
small fix
[u/mrichter/AliRoot.git] / PHOS / AliPHOS.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 /* $Id$ */
16 /* History of cvs commits:
17  *
18  * $Log$
19  * Revision 1.97  2006/04/22 10:30:17  hristov
20  * Add fEnergy to AliPHOSDigit and operate with EMC amplitude in energy units (Yu.Kharlov)
21  *
22  * Revision 1.96  2006/04/07 08:41:59  hristov
23  * Follow AliAlignObj framework and remove AliPHOSAlignData (Yu.Kharlov)
24  *
25  * Revision 1.95  2006/03/14 19:40:41  kharlov
26  * Remove De-digitizing of raw data and digitizing the raw data fit
27  *
28  * Revision 1.94  2006/03/07 18:56:25  kharlov
29  * CDB is passed via environment variable
30  *
31  * Revision 1.93  2005/11/22 08:45:11  kharlov
32  * Calibration is read from CDB if any (Boris Polichtchouk)
33  *
34  * Revision 1.92  2005/11/03 13:09:19  hristov
35  * Removing meaningless const declarations (linuxicc)
36  *
37  * Revision 1.91  2005/07/27 15:08:53  kharlov
38  * Mixture ArCO2 is corrected
39  *
40  * Revision 1.90  2005/06/17 07:39:07  hristov
41  * Removing GetDebug and SetDebug from AliRun and AliModule. Using AliLog for the messages
42  *
43  * Revision 1.89  2005/05/28 12:10:07  schutz
44  * Copy constructor is corrected (by T.P.)
45  *
46  */
47
48 //_________________________________________________________________________
49 // Base Class for PHOS description:
50 //   PHOS consists of a PbWO4 calorimeter (EMCA) and a gazeous charged 
51 //    particles detector (CPV or PPSD).
52 //   The only provided method here is CreateMaterials, 
53 //    which defines the materials common to all PHOS versions.   
54 // 
55 //*-- Author: Laurent Aphecetche & Yves Schutz (SUBATECH) 
56 //////////////////////////////////////////////////////////////////////////////
57
58
59 // --- ROOT system ---
60 class TFile;
61 #include <TFolder.h> 
62 #include <TTree.h>
63 #include <TVirtualMC.h> 
64 #include <TH1F.h> 
65 #include <TF1.h> 
66 #include <TRandom.h> 
67
68 // --- Standard library ---
69
70 // --- AliRoot header files ---
71 #include "AliMagF.h"
72 #include "AliPHOS.h"
73 #include "AliPHOSGetter.h"
74 #include "AliRun.h"
75 #include "AliPHOSDigitizer.h"
76 #include "AliPHOSSDigitizer.h"
77 #include "AliPHOSDigit.h"
78 #include "AliAltroBuffer.h"
79 #include "AliLog.h"
80 #include "AliCDBManager.h"
81 #include "AliCDBEntry.h"
82 #include "AliCDBStorage.h"
83 #include "AliPHOSCalibData.h"
84
85 ClassImp(AliPHOS)
86
87 Double_t AliPHOS::fgCapa        = 1.;        // 1pF 
88 Int_t    AliPHOS::fgOrder       = 2 ;
89 Double_t AliPHOS::fgTimeMax     = 2.56E-5 ;  // each sample is over 100 ns fTimeMax/fTimeBins
90 Double_t AliPHOS::fgTimePeak    = 4.1E-6 ;   // 4 micro seconds
91 Double_t AliPHOS::fgTimeTrigger = 100E-9 ;      // 100ns, just for a reference
92
93 //____________________________________________________________________________
94   AliPHOS:: AliPHOS() : AliDetector()
95 {
96   // Default ctor
97   fName   = "PHOS" ;
98
99 }
100
101 //____________________________________________________________________________
102 AliPHOS::AliPHOS(const char* name, const char* title): AliDetector(name, title)
103 {
104   //   ctor : title is used to identify the layout
105
106   fHighCharge        = 8.2 ;          // adjusted for a high gain range of 5.12 GeV (10 bits)
107   fHighGain          = 6.64 ; 
108   fHighLowGainFactor = 16. ;          // adjusted for a low gain range of 82 GeV (10 bits) 
109   fLowGainOffset     = GetGeometry()->GetNModules() + 1 ;   
110     // offset added to the module id to distinguish high and low gain data
111 }
112
113 //____________________________________________________________________________
114 AliPHOS::~AliPHOS() 
115 {  
116 }
117
118 //____________________________________________________________________________
119 void AliPHOS::Copy(TObject &obj)const
120 {
121   // copy method to be used by the cpy ctor
122   TObject::Copy(obj);
123   
124   AliPHOS &phos = static_cast<AliPHOS &>(obj); 
125   
126   phos.fHighCharge        = fHighCharge ;
127   phos.fHighGain          = fHighGain ; 
128   phos.fHighLowGainFactor = fHighLowGainFactor ;  
129   phos.fLowGainOffset     = fLowGainOffset;   
130 }
131
132 //____________________________________________________________________________
133 AliDigitizer* AliPHOS::CreateDigitizer(AliRunDigitizer* manager) const
134 {
135   return new AliPHOSDigitizer(manager);
136 }
137
138 //____________________________________________________________________________
139 void AliPHOS::CreateMaterials()
140 {
141   // Definitions of materials to build PHOS and associated tracking media.
142   // media number in idtmed are 699 to 798.
143
144   // --- The PbWO4 crystals ---
145   Float_t aX[3] = {207.19, 183.85, 16.0} ;
146   Float_t zX[3] = {82.0, 74.0, 8.0} ;
147   Float_t wX[3] = {1.0, 1.0, 4.0} ;
148   Float_t dX = 8.28 ;
149
150   AliMixture(0, "PbWO4$", aX, zX, dX, -3, wX) ;
151
152
153   // --- The polysterene scintillator (CH) ---
154   Float_t aP[2] = {12.011, 1.00794} ;
155   Float_t zP[2] = {6.0, 1.0} ;
156   Float_t wP[2] = {1.0, 1.0} ;
157   Float_t dP = 1.032 ;
158
159   AliMixture(1, "Polystyrene$", aP, zP, dP, -2, wP) ;
160
161   // --- Aluminium ---
162   AliMaterial(2, "Al$", 26.98, 13., 2.7, 8.9, 999., 0, 0) ;
163   // ---         Absorption length is ignored ^
164
165  // --- Tyvek (CnH2n) ---
166   Float_t aT[2] = {12.011, 1.00794} ;
167   Float_t zT[2] = {6.0, 1.0} ;
168   Float_t wT[2] = {1.0, 2.0} ;
169   Float_t dT = 0.331 ;
170
171   AliMixture(3, "Tyvek$", aT, zT, dT, -2, wT) ;
172
173   // --- Polystyrene foam ---
174   Float_t aF[2] = {12.011, 1.00794} ;
175   Float_t zF[2] = {6.0, 1.0} ;
176   Float_t wF[2] = {1.0, 1.0} ;
177   Float_t dF = 0.12 ;
178
179   AliMixture(4, "Foam$", aF, zF, dF, -2, wF) ;
180
181  // --- Titanium ---
182   Float_t aTIT[3] = {47.88, 26.98, 54.94} ;
183   Float_t zTIT[3] = {22.0, 13.0, 25.0} ;
184   Float_t wTIT[3] = {69.0, 6.0, 1.0} ;
185   Float_t dTIT = 4.5 ;
186
187   AliMixture(5, "Titanium$", aTIT, zTIT, dTIT, -3, wTIT);
188
189  // --- Silicon ---
190   AliMaterial(6, "Si$", 28.0855, 14., 2.33, 9.36, 42.3, 0, 0) ;
191
192
193
194   // --- Foam thermo insulation ---
195   Float_t aTI[2] = {12.011, 1.00794} ;
196   Float_t zTI[2] = {6.0, 1.0} ;
197   Float_t wTI[2] = {1.0, 1.0} ;
198   Float_t dTI = 0.04 ;
199
200   AliMixture(7, "Thermo Insul.$", aTI, zTI, dTI, -2, wTI) ;
201
202   // --- Textolith ---
203   Float_t aTX[4] = {16.0, 28.09, 12.011, 1.00794} ;
204   Float_t zTX[4] = {8.0, 14.0, 6.0, 1.0} ;
205   Float_t wTX[4] = {292.0, 68.0, 462.0, 736.0} ;
206   Float_t dTX    = 1.75 ;
207
208   AliMixture(8, "Textolit$", aTX, zTX, dTX, -4, wTX) ;
209
210   //--- FR4  ---
211   Float_t aFR[4] = {16.0, 28.09, 12.011, 1.00794} ;
212   Float_t zFR[4] = {8.0, 14.0, 6.0, 1.0} ;
213   Float_t wFR[4] = {292.0, 68.0, 462.0, 736.0} ;
214   Float_t dFR = 1.8 ; 
215
216   AliMixture(9, "FR4$", aFR, zFR, dFR, -4, wFR) ;
217
218   // --- The Composite Material for  micromegas (so far polyetylene) ---                                       
219   Float_t aCM[2] = {12.01, 1.} ; 
220   Float_t zCM[2] = {6., 1.} ; 
221   Float_t wCM[2] = {1., 2.} ; 
222   Float_t dCM = 0.935 ; 
223
224   AliMixture(10, "Compo Mat$", aCM, zCM, dCM, -2, wCM) ;
225
226   // --- Copper ---                                                                    
227   AliMaterial(11, "Cu$", 63.546, 29, 8.96, 1.43, 14.8, 0, 0) ;
228  
229   // --- G10 : Printed Circuit material ---                                                  
230   Float_t aG10[4] = { 12., 1., 16., 28.} ;
231   Float_t zG10[4] = { 6., 1., 8., 14.} ;
232   Float_t wG10[4] = { .259, .288, .248, .205} ;
233   Float_t dG10  = 1.7 ;
234   
235   AliMixture(12, "G10$", aG10, zG10, dG10, -4, wG10);
236
237   // --- Lead ---                                                                     
238   AliMaterial(13, "Pb$", 207.2, 82, 11.35, 0.56, 0., 0, 0) ;
239
240  // --- The gas mixture ---                                                                
241  // Co2
242   Float_t aCO[2] = {12.0, 16.0} ; 
243   Float_t zCO[2] = {6.0, 8.0} ; 
244   Float_t wCO[2] = {1.0, 2.0} ; 
245   Float_t dCO = 0.001977 ; 
246
247   AliMixture(14, "CO2$", aCO, zCO, dCO, -2, wCO);
248
249  // Ar
250   Float_t dAr = 0.001782 ; 
251   AliMaterial(15, "Ar$", 39.948, 18.0, dAr, 14.0, 0., 0, 0) ;   
252  
253   // Ar+CO2 Mixture (80% / 20%)
254   Float_t arContent = 0.80 ;  // Ar-content of the ArCO2-mixture
255   Float_t aArCO[3]  = {39.948, 12.0, 16.0} ;
256   Float_t zArCO[3]  = {18.0  ,  6.0,  8.0} ;
257   Float_t wArCO[3];
258   wArCO[0] = arContent;
259   wArCO[1] = (1-arContent)*1;
260   wArCO[2] = (1-arContent)*2;
261   Float_t dArCO = arContent*dAr + (1-arContent)*dCO ;
262   AliMixture(16, "ArCO2$", aArCO, zArCO, dArCO,  -3, wArCO) ;
263
264   // --- Stainless steel (let it be pure iron) ---
265   AliMaterial(17, "Steel$", 55.845, 26, 7.87, 1.76, 0., 0, 0) ;
266
267
268   // --- Fiberglass ---
269   Float_t aFG[4] = {16.0, 28.09, 12.011, 1.00794} ;
270   Float_t zFG[4] = {8.0, 14.0, 6.0, 1.0} ;
271   Float_t wFG[4] = {292.0, 68.0, 462.0, 736.0} ;
272   Float_t dFG    = 1.9 ;
273
274   AliMixture(18, "Fibergla$", aFG, zFG, dFG, -4, wFG) ;
275
276   // --- Cables in Air box  ---
277   // SERVICES
278
279   Float_t aCA[4] = { 1.,12.,55.8,63.5 };
280   Float_t zCA[4] = { 1.,6.,26.,29. }; 
281   Float_t wCA[4] = { .014,.086,.42,.48 };
282   Float_t dCA    = 0.8 ;  //this density is raw estimation, if you know better - correct
283
284   AliMixture(19, "Cables  $", aCA, zCA, dCA, -4, wCA) ;
285
286
287   // --- Air ---
288   Float_t aAir[4]={12.0107,14.0067,15.9994,39.948};
289   Float_t zAir[4]={6.,7.,8.,18.};
290   Float_t wAir[4]={0.000124,0.755267,0.231781,0.012827};
291   Float_t dAir = 1.20479E-3;
292  
293   AliMixture(99, "Air$", aAir, zAir, dAir, 4, wAir) ;
294
295   // DEFINITION OF THE TRACKING MEDIA
296
297   // for PHOS: idtmed[699->798] equivalent to fIdtmed[0->100]
298   Int_t * idtmed = fIdtmed->GetArray() - 699 ; 
299   Int_t   isxfld = gAlice->Field()->Integ() ;
300   Float_t sxmgmx = gAlice->Field()->Max() ;
301
302   // The scintillator of the calorimeter made of PBW04                              -> idtmed[699]
303   AliMedium(0, "PHOS Xtal    $", 0, 1,
304             isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
305
306   // The scintillator of the CPV made of Polystyrene scintillator                   -> idtmed[700]
307   AliMedium(1, "CPV scint.   $", 1, 1,
308             isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
309
310   // Various Aluminium parts made of Al                                             -> idtmed[701]
311   AliMedium(2, "Al parts     $", 2, 0,
312              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.001, 0.001, 0, 0) ;
313
314   // The Tywek which wraps the calorimeter crystals                                 -> idtmed[702]
315   AliMedium(3, "Tyvek wrapper$", 3, 0,
316              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.001, 0.001, 0, 0) ;
317
318   // The Polystyrene foam around the calorimeter module                             -> idtmed[703]
319   AliMedium(4, "Polyst. foam $", 4, 0,
320              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
321
322   // The Titanium around the calorimeter crystal                                    -> idtmed[704]
323   AliMedium(5, "Titan. cover $", 5, 0,
324              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.0001, 0.0001, 0, 0) ;
325
326   // The Silicon of the pin diode to read out the calorimeter crystal               -> idtmed[705] 
327  AliMedium(6, "Si PIN       $", 6, 0,
328              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.01, 0.01, 0, 0) ;
329
330  // The thermo insulating material of the box which contains the calorimeter module -> idtmed[706]
331   AliMedium(7, "Thermo Insul.$", 7, 0,
332              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
333
334   // The Textolit which makes up the box which contains the calorimeter module      -> idtmed[707]
335   AliMedium(8, "Textolit     $", 8, 0,
336              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
337
338   // FR4: The Plastic which makes up the frame of micromegas                        -> idtmed[708]
339   AliMedium(9, "FR4 $", 9, 0,
340              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.0001, 0, 0) ; 
341
342
343   // The Composite Material for  micromegas                                         -> idtmed[709]
344   AliMedium(10, "CompoMat   $", 10, 0,
345              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
346
347   // Copper                                                                         -> idtmed[710]
348   AliMedium(11, "Copper     $", 11, 0,
349              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.0001, 0, 0) ;
350
351   // G10: Printed Circuit material                                                  -> idtmed[711]
352  
353   AliMedium(12, "G10        $", 12, 0,
354              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.01, 0, 0) ;
355
356   // The Lead                                                                       -> idtmed[712]
357  
358   AliMedium(13, "Lead      $", 13, 0,
359              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
360
361   // The gas mixture: ArCo2                                                         -> idtmed[715]
362  
363   AliMedium(16, "ArCo2      $", 16, 1,
364              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.01, 0, 0) ;
365  
366   // Stainless steel                                                                -> idtmed[716]
367   AliMedium(17, "Steel     $", 17, 0,
368              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.0001, 0, 0) ;
369
370   // Fibergalss                                                                     -> idtmed[717]
371   AliMedium(18, "Fiberglass$", 18, 0,
372              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
373
374   // Cables in air                                                                  -> idtmed[718]
375   AliMedium(19, "Cables    $", 19, 0,
376              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
377
378   // Air                                                                            -> idtmed[798] 
379   AliMedium(99, "Air          $", 99, 0,
380              isxfld, sxmgmx, 10.0, 1.0, 0.1, 0.1, 10.0, 0, 0) ;
381
382   // --- Set decent energy thresholds for gamma and electron tracking
383
384   // Tracking threshold for photons and electrons in the scintillator crystal 
385   gMC->Gstpar(idtmed[699], "CUTGAM",0.5E-4) ; 
386   gMC->Gstpar(idtmed[699], "CUTELE",1.0E-4) ;
387  
388   // --- Generate explicitly delta rays in the titan cover ---
389   gMC->Gstpar(idtmed[704], "LOSS",3.) ;
390   gMC->Gstpar(idtmed[704], "DRAY",1.) ;
391   // --- and in aluminium parts ---
392   gMC->Gstpar(idtmed[701], "LOSS",3.) ;
393   gMC->Gstpar(idtmed[701], "DRAY",1.) ;
394   // --- and in PIN diode
395   gMC->Gstpar(idtmed[705], "LOSS",3) ;
396   gMC->Gstpar(idtmed[705], "DRAY",1) ;
397   // --- and in the passive convertor
398   gMC->Gstpar(idtmed[712], "LOSS",3) ;
399   gMC->Gstpar(idtmed[712], "DRAY",1) ;
400   // Tracking threshold for photons and electrons in the gas ArC02 
401   gMC->Gstpar(idtmed[715], "CUTGAM",1.E-5) ; 
402   gMC->Gstpar(idtmed[715], "CUTELE",1.E-5) ;
403   gMC->Gstpar(idtmed[715], "CUTNEU",1.E-5) ;
404   gMC->Gstpar(idtmed[715], "CUTHAD",1.E-5) ;
405   gMC->Gstpar(idtmed[715], "CUTMUO",1.E-5) ;
406   gMC->Gstpar(idtmed[715], "BCUTE",1.E-5) ;
407   gMC->Gstpar(idtmed[715], "BCUTM",1.E-5) ;
408   gMC->Gstpar(idtmed[715], "DCUTE",1.E-5) ;
409   gMC->Gstpar(idtmed[715], "DCUTM",1.E-5) ;
410   gMC->Gstpar(idtmed[715], "PPCUTM",1.E-5) ;
411   gMC->Gstpar(idtmed[715], "LOSS",2.) ;
412   gMC->Gstpar(idtmed[715], "DRAY",0.) ;
413   gMC->Gstpar(idtmed[715], "STRA",2.) ;
414
415 }
416
417 //____________________________________________________________________________
418 void AliPHOS::Digits2Raw()
419 {
420 // convert digits of the current event to raw data
421   
422   AliPHOSLoader * loader = dynamic_cast<AliPHOSLoader*>(fLoader) ; 
423
424   // get the digits
425   loader->LoadDigits();
426   TClonesArray* digits = loader->Digits() ;
427
428   if (!digits) {
429     AliError(Form("No digits found !"));
430     return;
431   }
432
433   // get the geometry
434   AliPHOSGeometry* geom = GetGeometry();
435   if (!geom) {
436     AliError(Form("No geometry found !"));
437     return;
438   }
439
440   // some digitization constants
441   const Int_t    kDDLOffset = 0x600; // assigned to PHOS
442 //   const Int_t    kThreshold = 1; // skip digits below this threshold // YVK
443   const Float_t    kThreshold = 0.001; // skip digits below 1 MeV
444   const Int_t      kAdcThreshold = 1;  // Lower ADC threshold to write to raw data
445
446   AliAltroBuffer* buffer = NULL;
447   Int_t prevDDL = -1;
448   Int_t adcValuesLow[fkTimeBins];
449   Int_t adcValuesHigh[fkTimeBins];
450
451   // loop over digits (assume ordered digits)
452   for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
453     AliPHOSDigit* digit = dynamic_cast<AliPHOSDigit *>(digits->At(iDigit)) ;
454     if (digit->GetEnergy() < kThreshold) 
455       continue;
456     Int_t relId[4];
457     geom->AbsToRelNumbering(digit->GetId(), relId);
458     Int_t module = relId[0];
459  
460    // Begin FIXME 
461     if (relId[1] != 0) 
462       continue;    // ignore digits from CPV
463    // End FIXME 
464
465     // PHOS EMCA has 4 DDL per module. Splitting is done based on the row number
466     Int_t iDDL = 4 * (module - 1) + (4 * (relId[2] - 1)) / geom->GetNPhi();
467
468     // new DDL
469     if (iDDL != prevDDL) {
470       // write real header and close previous file
471       if (buffer) {
472         buffer->Flush();
473         buffer->WriteDataHeader(kFALSE, kFALSE);
474         delete buffer;
475       }
476
477       // open new file and write dummy header
478       TString fileName("PHOS_") ;
479       fileName += (iDDL + kDDLOffset) ; 
480       fileName += ".ddl" ; 
481       buffer = new AliAltroBuffer(fileName.Data());
482       buffer->WriteDataHeader(kTRUE, kFALSE);  //Dummy;
483
484       prevDDL = iDDL;
485     }
486
487     // out of time range signal (?)
488     if (digit->GetTimeR() > GetRawFormatTimeMax() ) {
489       buffer->FillBuffer((Int_t)digit->GetEnergy());
490       buffer->FillBuffer(GetRawFormatTimeBins() );  // time bin
491       buffer->FillBuffer(3);          // bunch length      
492       buffer->WriteTrailer(3, relId[3], relId[2], module);  // trailer
493       
494     // calculate the time response function
495     } else {
496       Double_t energy = 0 ;
497       Int_t   module = relId[0];
498       if ( digit->GetId() <= geom->GetNModules() *  geom->GetNCristalsInModule()) {
499         energy=digit->GetEnergy();
500       }
501       else {
502 //      energy = digit->GetAmp()*digitizer->GetCPVchannel()+digitizer->GetCPVpedestal();
503         energy = 0; // CPV raw data format is now know yet
504       }        
505       Bool_t lowgain = RawSampledResponse(digit->GetTimeR(), energy, adcValuesHigh, adcValuesLow) ; 
506       
507       if (lowgain) 
508         buffer->WriteChannel(relId[3], relId[2], module + fLowGainOffset, 
509                              GetRawFormatTimeBins(), adcValuesLow , kAdcThreshold);
510       else 
511         buffer->WriteChannel(relId[3], relId[2], module, 
512                              GetRawFormatTimeBins(), adcValuesHigh, kAdcThreshold);
513       
514     }
515   }
516   
517   // write real header and close last file
518   if (buffer) {
519     buffer->Flush();
520     buffer->WriteDataHeader(kFALSE, kFALSE);
521     delete buffer;
522   }
523   
524   loader->UnloadDigits();
525 }
526
527 //____________________________________________________________________________
528 void AliPHOS::Hits2SDigits()  
529
530 // create summable digits
531
532   AliPHOSSDigitizer phosDigitizer(fLoader->GetRunLoader()->GetFileName().Data()) ;
533   phosDigitizer.SetEventRange(0, -1) ; // do all the events
534   phosDigitizer.ExecuteTask("all") ; 
535 }
536
537 //____________________________________________________________________________
538 AliLoader* AliPHOS::MakeLoader(const char* topfoldername)
539 {
540 //different behaviour than standard (singleton getter)
541 // --> to be discussed and made eventually coherent
542  fLoader = new AliPHOSLoader(GetName(),topfoldername);
543  return fLoader;
544 }
545
546 //__________________________________________________________________
547 Double_t AliPHOS::RawResponseFunction(Double_t *x, Double_t *par) 
548 {
549   // Shape of the electronics raw reponse:
550   // It is a semi-gaussian, 2nd order Gamma function of the general form
551   // v(t) = n**n * Q * A**n / C *(t/tp)**n * exp(-n * t/tp) with 
552   // tp : peaking time par[0]
553   // n  : order of the function
554   // C  : integrating capacitor in the preamplifier
555   // A  : open loop gain of the preamplifier
556   // Q  : the total APD charge to be measured Q = C * energy
557   
558   Double_t signal ;
559   Double_t xx = x[0] - ( fgTimeTrigger + par[3] ) ; 
560
561   if (xx < 0 || xx > fgTimeMax) 
562     signal = 0. ;  
563   else { 
564     Double_t fac = par[0] * TMath::Power(fgOrder, fgOrder) * TMath::Power(par[1], fgOrder) / fgCapa ; 
565     signal = fac * par[2] * TMath::Power(xx / fgTimePeak, fgOrder) * TMath::Exp(-fgOrder * (xx / fgTimePeak)) ; 
566   }
567   return signal ;  
568 }
569
570 //__________________________________________________________________
571 Double_t AliPHOS::RawResponseFunctionMax(Double_t charge, Double_t gain) 
572 {
573   return ( charge * TMath::Power(fgOrder, fgOrder) * TMath::Power(gain, fgOrder) 
574      / ( fgCapa * TMath::Exp(fgOrder) ) );  
575
576 }
577
578 //__________________________________________________________________
579 Bool_t AliPHOS::RawSampledResponse(Double_t dtime, Double_t damp, Int_t * adcH, Int_t * adcL) const 
580 {
581   // for a start time dtime and an amplitude damp given by digit, 
582   // calculates the raw sampled response AliPHOS::RawResponseFunction
583   // Input: dtime - signal start time
584   //        damp  - signal amplitude (energy)
585   // Output: adcH - array[fkTimeBins] of 10-bit samples for high-gain channel
586   //         adcL - array[fkTimeBins] of 10-bit samples for low-gain channel
587
588   const Int_t kRawSignalOverflow = 0x3FF ; 
589   Bool_t lowGain = kFALSE ; 
590
591   TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 4);
592
593   for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
594     signalF.SetParameter(0, GetRawFormatHighCharge() ) ; 
595     signalF.SetParameter(1, GetRawFormatHighGain() ) ; 
596     signalF.SetParameter(2, damp) ; 
597     signalF.SetParameter(3, dtime) ; 
598     Double_t time = iTime * GetRawFormatTimeMax() / GetRawFormatTimeBins() ;
599     Double_t signal = signalF.Eval(time) ;     
600     if ( static_cast<Int_t>(signal+0.5) > kRawSignalOverflow ){  // larger than 10 bits 
601       signal = kRawSignalOverflow ;
602       lowGain = kTRUE ; 
603     }
604     adcH[iTime] =  static_cast<Int_t>(signal + 0.5) ;
605
606     signalF.SetParameter(0, GetRawFormatLowCharge() ) ;     
607     signalF.SetParameter(1, GetRawFormatLowGain() ) ; 
608     signal = signalF.Eval(time) ;  
609     if ( static_cast<Int_t>(signal+0.5) > kRawSignalOverflow)  // larger than 10 bits 
610       signal = kRawSignalOverflow ;
611     adcL[iTime] = static_cast<Int_t>(0.5 + signal ) ; 
612
613   }
614   return lowGain ; 
615 }
616
617 //____________________________________________________________________________
618 void AliPHOS::SetTreeAddress()
619
620   // Links Hits in the Tree to Hits array
621   TBranch *branch;
622   char branchname[20];
623   sprintf(branchname,"%s",GetName());
624   // Branch address for hit tree
625     TTree *treeH = TreeH();
626   if (treeH) {
627     branch = treeH->GetBranch(branchname);
628     if (branch) 
629      { 
630        if (fHits == 0x0) fHits= new TClonesArray("AliPHOSHit",1000);
631        //AliInfo(Form("<%s> Setting Hits Address",GetName()));
632        branch->SetAddress(&fHits);
633      }
634   }
635 }
636