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