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