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