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