Changing name QualAss to QA
[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.115  2007/08/22 09:20:50  hristov
20  * Updated QA classes (Yves)
21  *
22  * Revision 1.114  2007/08/07 14:12:03  kharlov
23  * Quality assurance added (Yves Schutz)
24  *
25  * Revision 1.113  2007/07/18 16:29:54  policheh
26  * Raw Sdigits energy converted to GeV.
27  *
28  * Revision 1.112  2007/02/25 22:59:13  policheh
29  * Digits2Raw(): ALTRO buffer and mapping created per each DDL.
30  *
31  * Revision 1.111  2007/02/18 15:21:47  kharlov
32  * Corrections for memory leak in Digits2Raw due to AliAltroMapping
33  *
34  * Revision 1.110  2007/02/13 10:52:08  policheh
35  * Raw2SDigits() implemented
36  *
37  * Revision 1.109  2007/02/05 10:43:25  hristov
38  * Changes for correct initialization of Geant4 (Mihaela)
39  *
40  * Revision 1.108  2007/02/01 10:34:47  hristov
41  * Removing warnings on Solaris x86
42  *
43  * Revision 1.107  2007/01/29 16:29:37  kharlov
44  * Digits2Raw(): special workaround for digits with time out of range
45  *
46  * Revision 1.106  2007/01/17 17:28:56  kharlov
47  * Extract ALTRO sample generation to a separate class AliPHOSPulseGenerator
48  *
49  * Revision 1.105  2007/01/12 21:44:29  kharlov
50  * Simulate and reconstruct two gains simulaneouslsy
51  */
52
53 //_________________________________________________________________________
54 // Base Class for PHOS description:
55 //   PHOS consists of a PbWO4 calorimeter (EMCA) and a gazeous charged 
56 //    particles detector (CPV or PPSD).
57 //   The only provided method here is CreateMaterials, 
58 //    which defines the materials common to all PHOS versions.   
59 // 
60 //*-- Author: Laurent Aphecetche & Yves Schutz (SUBATECH) 
61 //////////////////////////////////////////////////////////////////////////////
62
63
64 // --- ROOT system ---
65 class TFile;
66 #include <TFolder.h> 
67 #include <TTree.h>
68 #include <TVirtualMC.h> 
69 #include <TH1F.h> 
70 #include <TF1.h> 
71 #include <TRandom.h> 
72
73 // --- Standard library ---
74
75 // --- AliRoot header files ---
76 #include "AliMagF.h"
77 #include "AliPHOS.h"
78 #include "AliPHOSGetter.h"
79 #include "AliRun.h"
80 #include "AliPHOSDigitizer.h"
81 #include "AliPHOSSDigitizer.h"
82 #include "AliPHOSDigit.h"
83 #include "AliAltroBuffer.h"
84 #include "AliAltroMapping.h"
85 #include "AliCaloAltroMapping.h"
86 #include "AliLog.h"
87 #include "AliCDBManager.h"
88 #include "AliCDBEntry.h"
89 #include "AliCDBStorage.h"
90 #include "AliPHOSCalibData.h"
91 #include "AliPHOSPulseGenerator.h"
92 #include "AliDAQ.h"
93 #include "AliPHOSRawDecoder.h"
94 #include "AliPHOSRawDigiProducer.h"
95 #include "AliPHOSQAChecker.h"
96
97 ClassImp(AliPHOS)
98
99 //____________________________________________________________________________
100   AliPHOS:: AliPHOS() : AliDetector()
101 {
102   // Default ctor
103   fName   = "PHOS" ;
104
105 }
106
107 //____________________________________________________________________________
108 AliPHOS::AliPHOS(const char* name, const char* title): AliDetector(name, title)
109 {
110   //   ctor : title is used to identify the layout
111 }
112
113 //____________________________________________________________________________
114 AliPHOS::~AliPHOS() 
115 {  
116 }
117
118 //____________________________________________________________________________
119 AliDigitizer* AliPHOS::CreateDigitizer(AliRunDigitizer* manager) const
120 {
121   return new AliPHOSDigitizer(manager);
122 }
123
124 //____________________________________________________________________________
125 void AliPHOS::CreateMaterials()
126 {
127   // Definitions of materials to build PHOS and associated tracking media.
128   // media number in idtmed are 699 to 798.
129
130   // --- The PbWO4 crystals ---
131   Float_t aX[3] = {207.19, 183.85, 16.0} ;
132   Float_t zX[3] = {82.0, 74.0, 8.0} ;
133   Float_t wX[3] = {1.0, 1.0, 4.0} ;
134   Float_t dX = 8.28 ;
135
136   AliMixture(0, "PbWO4$", aX, zX, dX, -3, wX) ;
137
138
139   // --- The polysterene scintillator (CH) ---
140   Float_t aP[2] = {12.011, 1.00794} ;
141   Float_t zP[2] = {6.0, 1.0} ;
142   Float_t wP[2] = {1.0, 1.0} ;
143   Float_t dP = 1.032 ;
144
145   AliMixture(1, "Polystyrene$", aP, zP, dP, -2, wP) ;
146
147   // --- Aluminium ---
148   AliMaterial(2, "Al$", 26.98, 13., 2.7, 8.9, 999., 0, 0) ;
149   // ---         Absorption length is ignored ^
150
151  // --- Tyvek (CnH2n) ---
152   Float_t aT[2] = {12.011, 1.00794} ;
153   Float_t zT[2] = {6.0, 1.0} ;
154   Float_t wT[2] = {1.0, 2.0} ;
155   Float_t dT = 0.331 ;
156
157   AliMixture(3, "Tyvek$", aT, zT, dT, -2, wT) ;
158
159   // --- Polystyrene foam ---
160   Float_t aF[2] = {12.011, 1.00794} ;
161   Float_t zF[2] = {6.0, 1.0} ;
162   Float_t wF[2] = {1.0, 1.0} ;
163   Float_t dF = 0.12 ;
164
165   AliMixture(4, "Foam$", aF, zF, dF, -2, wF) ;
166
167  // --- Titanium ---
168   Float_t aTIT[3] = {47.88, 26.98, 54.94} ;
169   Float_t zTIT[3] = {22.0, 13.0, 25.0} ;
170   Float_t wTIT[3] = {69.0, 6.0, 1.0} ;
171   Float_t dTIT = 4.5 ;
172
173   AliMixture(5, "Titanium$", aTIT, zTIT, dTIT, -3, wTIT);
174
175  // --- Silicon ---
176   AliMaterial(6, "Si$", 28.0855, 14., 2.33, 9.36, 42.3, 0, 0) ;
177
178
179
180   // --- Foam thermo insulation ---
181   Float_t aTI[2] = {12.011, 1.00794} ;
182   Float_t zTI[2] = {6.0, 1.0} ;
183   Float_t wTI[2] = {1.0, 1.0} ;
184   Float_t dTI = 0.04 ;
185
186   AliMixture(7, "Thermo Insul.$", aTI, zTI, dTI, -2, wTI) ;
187
188   // --- Textolith ---
189   Float_t aTX[4] = {16.0, 28.09, 12.011, 1.00794} ;
190   Float_t zTX[4] = {8.0, 14.0, 6.0, 1.0} ;
191   Float_t wTX[4] = {292.0, 68.0, 462.0, 736.0} ;
192   Float_t dTX    = 1.75 ;
193
194   AliMixture(8, "Textolit$", aTX, zTX, dTX, -4, wTX) ;
195
196   //--- FR4  ---
197   Float_t aFR[4] = {16.0, 28.09, 12.011, 1.00794} ;
198   Float_t zFR[4] = {8.0, 14.0, 6.0, 1.0} ;
199   Float_t wFR[4] = {292.0, 68.0, 462.0, 736.0} ;
200   Float_t dFR = 1.8 ; 
201
202   AliMixture(9, "FR4$", aFR, zFR, dFR, -4, wFR) ;
203
204   // --- The Composite Material for  micromegas (so far polyetylene) ---                                       
205   Float_t aCM[2] = {12.01, 1.} ; 
206   Float_t zCM[2] = {6., 1.} ; 
207   Float_t wCM[2] = {1., 2.} ; 
208   Float_t dCM = 0.935 ; 
209
210   AliMixture(10, "Compo Mat$", aCM, zCM, dCM, -2, wCM) ;
211
212   // --- Copper ---                                                                    
213   AliMaterial(11, "Cu$", 63.546, 29, 8.96, 1.43, 14.8, 0, 0) ;
214  
215   // --- G10 : Printed Circuit material ---                                                  
216   Float_t aG10[4] = { 12., 1., 16., 28.} ;
217   Float_t zG10[4] = { 6., 1., 8., 14.} ;
218   Float_t wG10[4] = { .259, .288, .248, .205} ;
219   Float_t dG10  = 1.7 ; 
220
221   
222   AliMixture(12, "G10$", aG10, zG10, dG10, -4, wG10);
223
224   // --- Lead ---                                                                     
225   AliMaterial(13, "Pb$", 207.2, 82, 11.35, 0.56, 0., 0, 0) ;
226
227  // --- The gas mixture ---                                                                
228  // Co2
229   Float_t aCO[2] = {12.0, 16.0} ; 
230   Float_t zCO[2] = {6.0, 8.0} ; 
231   Float_t wCO[2] = {1.0, 2.0} ; 
232   Float_t dCO = 0.001977 ; 
233
234   AliMixture(14, "CO2$", aCO, zCO, dCO, -2, wCO);
235
236  // Ar
237   Float_t dAr = 0.001782 ; 
238   AliMaterial(15, "Ar$", 39.948, 18.0, dAr, 14.0, 0., 0, 0) ;   
239  
240   // Ar+CO2 Mixture (80% / 20%)
241   Float_t arContent = 0.80 ;  // Ar-content of the ArCO2-mixture
242   Float_t aArCO[3]  = {39.948, 12.0, 16.0} ;
243   Float_t zArCO[3]  = {18.0  ,  6.0,  8.0} ;
244   Float_t wArCO[3];
245   wArCO[0] = arContent;
246   wArCO[1] = (1-arContent)*1;
247   wArCO[2] = (1-arContent)*2;
248   Float_t dArCO = arContent*dAr + (1-arContent)*dCO ;
249   AliMixture(16, "ArCO2$", aArCO, zArCO, dArCO,  -3, wArCO) ; 
250
251
252   // --- Stainless steel (let it be pure iron) ---
253   AliMaterial(17, "Steel$", 55.845, 26, 7.87, 1.76, 0., 0, 0) ;
254
255
256   // --- Fiberglass ---
257   Float_t aFG[4] = {16.0, 28.09, 12.011, 1.00794} ;
258   Float_t zFG[4] = {8.0, 14.0, 6.0, 1.0} ;
259   Float_t wFG[4] = {292.0, 68.0, 462.0, 736.0} ;
260   Float_t dFG    = 1.9 ;
261
262   AliMixture(18, "Fibergla$", aFG, zFG, dFG, -4, wFG) ;
263
264   // --- Cables in Air box  ---
265   // SERVICES
266
267   Float_t aCA[4] = { 1.,12.,55.8,63.5 };
268   Float_t zCA[4] = { 1.,6.,26.,29. }; 
269   Float_t wCA[4] = { .014,.086,.42,.48 };
270   Float_t dCA    = 0.8 ;  //this density is raw estimation, if you know better - correct
271
272   AliMixture(19, "Cables  $", aCA, zCA, dCA, -4, wCA) ;
273
274
275   // --- Air ---
276   Float_t aAir[4]={12.0107,14.0067,15.9994,39.948};
277   Float_t zAir[4]={6.,7.,8.,18.};
278   Float_t wAir[4]={0.000124,0.755267,0.231781,0.012827};
279   Float_t dAir = 1.20479E-3;
280  
281   AliMixture(99, "Air$", aAir, zAir, dAir, 4, wAir) ;
282
283   // DEFINITION OF THE TRACKING MEDIA
284
285   // for PHOS: idtmed[699->798] equivalent to fIdtmed[0->100]
286   Int_t   isxfld = gAlice->Field()->Integ() ;
287   Float_t sxmgmx = gAlice->Field()->Max() ;
288
289   // The scintillator of the calorimeter made of PBW04                              -> idtmed[699]
290   AliMedium(0, "PHOS Xtal    $", 0, 1,
291             isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
292
293   // The scintillator of the CPV made of Polystyrene scintillator                   -> idtmed[700]
294   AliMedium(1, "CPV scint.   $", 1, 1,
295             isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
296
297   // Various Aluminium parts made of Al                                             -> idtmed[701]
298   AliMedium(2, "Al parts     $", 2, 0,
299              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.001, 0.001, 0, 0) ;
300
301   // The Tywek which wraps the calorimeter crystals                                 -> idtmed[702]
302   AliMedium(3, "Tyvek wrapper$", 3, 0,
303              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.001, 0.001, 0, 0) ;
304
305   // The Polystyrene foam around the calorimeter module                             -> idtmed[703]
306   AliMedium(4, "Polyst. foam $", 4, 0,
307              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
308
309   // The Titanium around the calorimeter crystal                                    -> idtmed[704]
310   AliMedium(5, "Titan. cover $", 5, 0,
311              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.0001, 0.0001, 0, 0) ;
312
313   // The Silicon of the pin diode to read out the calorimeter crystal               -> idtmed[705] 
314  AliMedium(6, "Si PIN       $", 6, 0,
315              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.01, 0.01, 0, 0) ;
316
317  // The thermo insulating material of the box which contains the calorimeter module -> idtmed[706]
318   AliMedium(7, "Thermo Insul.$", 7, 0,
319              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
320
321   // The Textolit which makes up the box which contains the calorimeter module      -> idtmed[707]
322   AliMedium(8, "Textolit     $", 8, 0,
323              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
324
325   // FR4: The Plastic which makes up the frame of micromegas                        -> idtmed[708]
326   AliMedium(9, "FR4 $", 9, 0,
327              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.0001, 0, 0) ; 
328
329
330   // The Composite Material for  micromegas                                         -> idtmed[709]
331   AliMedium(10, "CompoMat   $", 10, 0,
332              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
333
334   // Copper                                                                         -> idtmed[710]
335   AliMedium(11, "Copper     $", 11, 0,
336              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.0001, 0, 0) ;
337
338   // G10: Printed Circuit material                                                  -> idtmed[711]
339  
340   AliMedium(12, "G10        $", 12, 0,
341              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.01, 0, 0) ;
342
343   // The Lead                                                                       -> idtmed[712]
344  
345   AliMedium(13, "Lead      $", 13, 0,
346              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
347
348   // The gas mixture: ArCo2                                                         -> idtmed[715]
349  
350   AliMedium(16, "ArCo2      $", 16, 1,
351              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.01, 0, 0) ;
352  
353   // Stainless steel                                                                -> idtmed[716]
354   AliMedium(17, "Steel     $", 17, 0,
355              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.0001, 0, 0) ;
356
357   // Fibergalss                                                                     -> idtmed[717]
358   AliMedium(18, "Fiberglass$", 18, 0,
359              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
360
361   // Cables in air                                                                  -> idtmed[718]
362   AliMedium(19, "Cables    $", 19, 0,
363              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
364
365   // Air                                                                            -> idtmed[798] 
366   AliMedium(99, "Air          $", 99, 0,
367              isxfld, sxmgmx, 10.0, 1.0, 0.1, 0.1, 10.0, 0, 0) ;
368 }
369
370 //_____________________________________________________________________________
371 void AliPHOS::Init()
372 {
373   //
374   // Initialises cuts for PHOS
375   //
376   // --- Set decent energy thresholds for gamma and electron tracking
377   Int_t * idtmed = fIdtmed->GetArray() - 699 ; 
378
379   // Tracking threshold for photons and electrons in the scintillator crystal 
380   gMC->Gstpar(idtmed[699], "CUTGAM",0.5E-4) ; 
381   gMC->Gstpar(idtmed[699], "CUTELE",1.0E-4) ;
382  
383   // --- Generate explicitly delta rays in the titan cover ---
384   gMC->Gstpar(idtmed[704], "LOSS",3.) ;
385   gMC->Gstpar(idtmed[704], "DRAY",1.) ;
386   // --- and in aluminium parts ---
387   gMC->Gstpar(idtmed[701], "LOSS",3.) ;
388   gMC->Gstpar(idtmed[701], "DRAY",1.) ;
389   // --- and in PIN diode
390   gMC->Gstpar(idtmed[705], "LOSS",3) ;
391   gMC->Gstpar(idtmed[705], "DRAY",1) ;
392   // --- and in the passive convertor
393   gMC->Gstpar(idtmed[712], "LOSS",3) ;
394   gMC->Gstpar(idtmed[712], "DRAY",1) ;
395   // Tracking threshold for photons and electrons in the gas ArC02 
396   gMC->Gstpar(idtmed[715], "CUTGAM",1.E-5) ; 
397   gMC->Gstpar(idtmed[715], "CUTELE",1.E-5) ;
398   gMC->Gstpar(idtmed[715], "CUTNEU",1.E-5) ;
399   gMC->Gstpar(idtmed[715], "CUTHAD",1.E-5) ;
400   gMC->Gstpar(idtmed[715], "CUTMUO",1.E-5) ;
401   gMC->Gstpar(idtmed[715], "BCUTE",1.E-5) ;
402   gMC->Gstpar(idtmed[715], "BCUTM",1.E-5) ;
403   gMC->Gstpar(idtmed[715], "DCUTE",1.E-5) ;
404   gMC->Gstpar(idtmed[715], "DCUTM",1.E-5) ;
405   gMC->Gstpar(idtmed[715], "PPCUTM",1.E-5) ;
406   gMC->Gstpar(idtmed[715], "LOSS",2.) ;
407   gMC->Gstpar(idtmed[715], "DRAY",0.) ;
408   gMC->Gstpar(idtmed[715], "STRA",2.) ;
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 geometry
428   AliPHOSGeometry* geom = GetGeometry();
429   if (!geom) {
430     AliError(Form("No geometry found !"));
431     return;
432   }
433
434   // some digitization constants
435   const Float_t    kThreshold = 0.001; // skip digits below 1 MeV
436   const Int_t      kAdcThreshold = 1;  // Lower ADC threshold to write to raw data
437
438   Int_t prevDDL = -1;
439
440   // Create a shaper pulse object
441   AliPHOSPulseGenerator *pulse = new AliPHOSPulseGenerator();
442   
443   Int_t *adcValuesLow = new Int_t[pulse->GetRawFormatTimeBins()];
444   Int_t *adcValuesHigh= new Int_t[pulse->GetRawFormatTimeBins()];
445   
446   const Int_t maxDDL = 20;
447   AliAltroBuffer  *buffer[maxDDL];
448   AliAltroMapping *mapping[maxDDL];
449
450   for(Int_t jDDL=0; jDDL<maxDDL; jDDL++) {
451     buffer[jDDL]=0;
452     mapping[jDDL]=0;
453   }
454
455   //!!!!for debug!!!
456   Int_t modMax=-111;
457   Int_t colMax=-111;
458   Int_t rowMax=-111;
459   Float_t eMax=-333;
460   //!!!for debug!!!
461
462   // loop over digits (assume ordered digits)
463   for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
464     AliPHOSDigit* digit = dynamic_cast<AliPHOSDigit *>(digits->At(iDigit)) ;
465     if (digit->GetEnergy() < kThreshold) 
466       continue;
467     Int_t relId[4];
468     geom->AbsToRelNumbering(digit->GetId(), relId);
469     Int_t module = relId[0];
470  
471     // Begin FIXME 
472     if (relId[1] != 0) 
473       continue;    // ignore digits from CPV
474     // End FIXME 
475
476     Int_t row = relId[2]-1;
477     Int_t col = relId[3]-1;
478
479     Int_t iRCU = -111;
480
481     //RCU0
482     if(0<=row&&row<32 && 0<=col&&col<28) iRCU=0;
483
484     //RCU1
485     if(0<=row&&row<32 && 28<=col&&col<56) iRCU=1;
486
487     //RCU2
488     if(32<=row&&row<64 && 0<=col&&col<28) iRCU=2;
489
490     //RCU3
491     if(32<=row&&row<64 && 28<=col&&col<56) iRCU=3;
492
493
494     // PHOS EMCA has 4 DDL per module. Splitting is based on the (row,column) numbers.
495     // PHOS internal convention: 1<module<5.
496     Int_t iDDL = 4 * (module - 1) + iRCU;
497
498     // new DDL
499     if (iDDL != prevDDL) {
500       if (buffer[iDDL] == 0) {
501         // open new file and write dummy header
502         TString fileName = AliDAQ::DdlFileName("PHOS",iDDL);
503
504         TString path = gSystem->Getenv("ALICE_ROOT");
505         path += "/PHOS/mapping/RCU";
506         path += iRCU;
507         path += ".data";
508
509         mapping[iDDL] = new AliCaloAltroMapping(path.Data());
510         buffer[iDDL]  = new AliAltroBuffer(fileName.Data(),mapping[iDDL]);
511         buffer[iDDL]->WriteDataHeader(kTRUE, kFALSE);  //Dummy;
512       }
513       prevDDL = iDDL;
514     }
515
516     AliDebug(2,Form("digit E=%.4f GeV, t=%g s, (mod,col,row)=(%d,%d,%d)\n",
517                     digit->GetEnergy(),digit->GetTimeR(),
518                     relId[0]-1,relId[3]-1,relId[2]-1));
519     // if a signal is out of time range, write only trailer
520     if (digit->GetTimeR() > pulse->GetRawFormatTimeMax()*0.5 ) {
521       AliInfo("Signal is out of time range.\n");
522       buffer[iDDL]->FillBuffer(0);
523       buffer[iDDL]->FillBuffer(pulse->GetRawFormatTimeBins() );  // time bin
524       buffer[iDDL]->FillBuffer(3);                               // bunch length
525       buffer[iDDL]->WriteTrailer(3, relId[3]-1, relId[2]-1, 0);  // trailer
526       
527     // calculate the time response function
528     } else {
529       Double_t energy = 0 ;
530       module = relId[0];
531       if (digit->GetId() <= geom->GetNModules() * geom->GetNCristalsInModule()) {
532         energy=digit->GetEnergy();
533         if(energy>eMax) {eMax=energy; modMax=module; colMax=col; rowMax=row;}
534       }
535       else {
536         energy = 0; // CPV raw data format is now know yet
537       }
538       pulse->SetAmplitude(energy);
539       pulse->SetTZero(digit->GetTimeR());
540       pulse->MakeSamples();
541       pulse->GetSamples(adcValuesHigh, adcValuesLow) ; 
542       buffer[iDDL]->WriteChannel(relId[3]-1, relId[2]-1, 0, 
543                            pulse->GetRawFormatTimeBins(), adcValuesLow , kAdcThreshold);
544       buffer[iDDL]->WriteChannel(relId[3]-1, relId[2]-1, 1, 
545                            pulse->GetRawFormatTimeBins(), adcValuesHigh, kAdcThreshold);
546     }
547   }
548   
549   // write real header and close last file
550   for (Int_t iDDL=0; iDDL<maxDDL; iDDL++) {
551     if (buffer[iDDL]) {
552       buffer[iDDL]->Flush();
553       buffer[iDDL]->WriteDataHeader(kFALSE, kFALSE);
554       delete buffer[iDDL];
555       if (mapping[iDDL]) delete mapping[iDDL];
556     }
557   }
558   
559   AliDebug(1,Form("Digit with max. energy:  modMax %d colMax %d rowMax %d  eMax %f\n",
560          modMax,colMax,rowMax,eMax));
561
562   delete pulse;
563   loader->UnloadDigits();
564 }
565
566 //____________________________________________________________________________
567 void AliPHOS::Hits2SDigits()  
568
569 // create summable digits
570
571   AliPHOSSDigitizer phosDigitizer(fLoader->GetRunLoader()->GetFileName().Data()) ;
572   phosDigitizer.SetEventRange(0, -1) ; // do all the events
573  
574   phosDigitizer.ExecuteTask("all") ; 
575 }
576
577
578 //____________________________________________________________________________
579 AliLoader* AliPHOS::MakeLoader(const char* topfoldername)
580 {
581 //different behaviour than standard (singleton getter)
582 // --> to be discussed and made eventually coherent
583  fLoader = new AliPHOSLoader(GetName(),topfoldername);
584  return fLoader;
585 }
586
587 //____________________________________________________________________________
588 void AliPHOS::SetTreeAddress()
589
590   // Links Hits in the Tree to Hits array
591   TBranch *branch;
592   char branchname[20];
593   sprintf(branchname,"%s",GetName());
594   // Branch address for hit tree
595     TTree *treeH = TreeH();
596   if (treeH) {
597     branch = treeH->GetBranch(branchname);
598     if (branch) 
599      { 
600        if (fHits == 0x0) fHits= new TClonesArray("AliPHOSHit",1000);
601        //AliInfo(Form("<%s> Setting Hits Address",GetName()));
602        branch->SetAddress(&fHits);
603      }
604   }
605 }
606
607 //____________________________________________________________________________
608 Bool_t AliPHOS::Raw2SDigits(AliRawReader* rawReader)
609 {
610
611   AliPHOSLoader * loader = dynamic_cast<AliPHOSLoader*>(fLoader) ;
612
613   TTree * tree = 0 ;
614   tree = loader->TreeS() ;
615   if ( !tree ) {
616     loader->MakeTree("S");
617     tree = loader->TreeS() ;
618   }
619
620   TClonesArray * sdigits = loader->SDigits() ;
621   if(!sdigits) { 
622     loader->MakeSDigitsArray();
623     sdigits = loader->SDigits();
624   }
625   sdigits->Clear();
626
627   AliPHOSRawDecoder dc(rawReader);
628   AliPHOSRawDigiProducer pr;
629   pr.MakeDigits(sdigits,&dc);
630
631   Int_t bufferSize = 32000 ;
632   // TBranch * sdigitsBranch = tree->Branch("PHOS",&sdigits,bufferSize);
633   tree->Branch("PHOS",&sdigits,bufferSize);
634   tree->Fill();
635
636   fLoader->WriteSDigits("OVERWRITE");
637   return kTRUE;
638     
639 }