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