]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOS.cxx
Additional protection. The call to RecWithStack has to be removed in case of raw...
[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.94  2006/03/07 18:56:25  kharlov
20  * CDB is passed via environment variable
21  *
22  * Revision 1.93  2005/11/22 08:45:11  kharlov
23  * Calibration is read from CDB if any (Boris Polichtchouk)
24  *
25  * Revision 1.92  2005/11/03 13:09:19  hristov
26  * Removing meaningless const declarations (linuxicc)
27  *
28  * Revision 1.91  2005/07/27 15:08:53  kharlov
29  * Mixture ArCO2 is corrected
30  *
31  * Revision 1.90  2005/06/17 07:39:07  hristov
32  * Removing GetDebug and SetDebug from AliRun and AliModule. Using AliLog for the messages
33  *
34  * Revision 1.89  2005/05/28 12:10:07  schutz
35  * Copy constructor is corrected (by T.P.)
36  *
37  */
38
39 //_________________________________________________________________________
40 // Base Class for PHOS description:
41 //   PHOS consists of a PbWO4 calorimeter (EMCA) and a gazeous charged 
42 //    particles detector (CPV or PPSD).
43 //   The only provided method here is CreateMaterials, 
44 //    which defines the materials common to all PHOS versions.   
45 // 
46 //*-- Author: Laurent Aphecetche & Yves Schutz (SUBATECH) 
47 //////////////////////////////////////////////////////////////////////////////
48
49
50 // --- ROOT system ---
51 class TFile;
52 #include <TFolder.h> 
53 #include <TTree.h>
54 #include <TVirtualMC.h> 
55 #include <TH1F.h> 
56 #include <TF1.h> 
57 #include <TRandom.h> 
58
59 // --- Standard library ---
60
61 // --- AliRoot header files ---
62 #include "AliMagF.h"
63 #include "AliPHOS.h"
64 #include "AliPHOSGetter.h"
65 #include "AliRun.h"
66 #include "AliPHOSDigitizer.h"
67 #include "AliPHOSSDigitizer.h"
68 #include "AliPHOSDigit.h"
69 #include "AliAltroBuffer.h"
70 #include "AliLog.h"
71 #include "AliCDBManager.h"
72 #include "AliCDBEntry.h"
73 #include "AliCDBStorage.h"
74 #include "AliPHOSCalibData.h"
75
76 ClassImp(AliPHOS)
77
78 Double_t AliPHOS::fgCapa        = 1.;        // 1pF 
79 Int_t    AliPHOS::fgOrder       = 2 ;
80 Double_t AliPHOS::fgTimeMax     = 2.56E-5 ;  // each sample is over 100 ns fTimeMax/fTimeBins
81 Double_t AliPHOS::fgTimePeak    = 4.1E-6 ;   // 4 micro seconds
82 Double_t AliPHOS::fgTimeTrigger = 100E-9 ;      // 100ns, just for a reference
83
84 //____________________________________________________________________________
85   AliPHOS:: AliPHOS() : AliDetector()
86 {
87   // Default ctor
88   fName   = "PHOS" ;
89
90 }
91
92 //____________________________________________________________________________
93 AliPHOS::AliPHOS(const char* name, const char* title): AliDetector(name, title)
94 {
95   //   ctor : title is used to identify the layout
96
97   // Check if CDB_PATH is defined and take alignment data from CDB
98   AliPHOSAlignData* alignda = 0;
99   if (gSystem->Getenv("CDB_PATH")) {
100     TString cdbPath = gSystem->Getenv("CDB_PATH");
101     AliCDBStorage *cdbStorage = AliCDBManager::Instance()->GetStorage(cdbPath);
102     if (cdbStorage != NULL) {
103       alignda  =
104         (AliPHOSAlignData*)(cdbStorage->Get("PHOS/Alignment/Geometry",0)->GetObject());
105       if(AliLog::GetGlobalDebugLevel()>0) alignda->Print();
106     }
107     else {
108       Fatal("AliPHOS", "No CDB storage at the path %s", cdbPath.Data()) ;
109     }
110   }
111   
112   fHighCharge        = 8.2 ;          // adjusted for a high gain range of 5.12 GeV (10 bits)
113   fHighGain          = 6.64 ; 
114   fHighLowGainFactor = 16. ;          // adjusted for a low gain range of 82 GeV (10 bits) 
115   fLowGainOffset     = GetGeometry(alignda)->GetNModules() + 1 ;   
116     // offset added to the module id to distinguish high and low gain data
117 }
118
119 //____________________________________________________________________________
120 AliPHOS::~AliPHOS() 
121 {  
122 }
123
124 //____________________________________________________________________________
125 void AliPHOS::Copy(TObject &obj)const
126 {
127   // copy method to be used by the cpy ctor
128   TObject::Copy(obj);
129   
130   AliPHOS &phos = static_cast<AliPHOS &>(obj); 
131   
132   phos.fHighCharge        = fHighCharge ;
133   phos.fHighGain          = fHighGain ; 
134   phos.fHighLowGainFactor = fHighLowGainFactor ;  
135   phos.fLowGainOffset     = fLowGainOffset;   
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 digitizer 
440   loader->LoadDigitizer();
441 //   AliPHOSDigitizer * digitizer = dynamic_cast<AliPHOSDigitizer *>(loader->Digitizer())  ; 
442   
443   // get the geometry
444   AliPHOSGeometry* geom = GetGeometry();
445   if (!geom) {
446     AliError(Form("No geometry found !"));
447     return;
448   }
449
450   // some digitization constants
451   const Int_t    kDDLOffset = 0x600; // assigned to PHOS
452   const Int_t    kThreshold = 1; // skip digits below this threshold
453
454   AliAltroBuffer* buffer = NULL;
455   Int_t prevDDL = -1;
456   Int_t adcValuesLow[fkTimeBins];
457   Int_t adcValuesHigh[fkTimeBins];
458
459 //   AliPHOSCalibData* calib=0;
460
461 //   //retrieve calibration database
462 //   if(AliCDBManager::Instance()->IsDefaultStorageSet()){
463 //     AliCDBEntry *entry = (AliCDBEntry*) AliCDBManager::Instance()->GetDefaultStorage()
464 //       ->Get("PHOS/GainFactors_and_Pedestals/Calibration",gAlice->GetRunNumber());
465 //     calib = (AliPHOSCalibData*) entry->GetObject();
466 //   }
467
468   // loop over digits (assume ordered digits)
469   for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
470     AliPHOSDigit* digit = dynamic_cast<AliPHOSDigit *>(digits->At(iDigit)) ;
471     if (digit->GetAmp() < kThreshold) 
472       continue;
473     Int_t relId[4];
474     geom->AbsToRelNumbering(digit->GetId(), relId);
475     Int_t module = relId[0];
476  
477    // Begin FIXME 
478     if (relId[1] != 0) 
479       continue;    // ignore digits from CPV
480    // End FIXME 
481
482     // PHOS EMCA has 4 DDL per module. Splitting is done based on the row number
483     Int_t iDDL = 4 * (module - 1) + (4 * (relId[2] - 1)) / geom->GetNPhi();
484
485     // new DDL
486     if (iDDL != prevDDL) {
487       // write real header and close previous file
488       if (buffer) {
489         buffer->Flush();
490         buffer->WriteDataHeader(kFALSE, kFALSE);
491         delete buffer;
492       }
493
494       // open new file and write dummy header
495       TString fileName("PHOS_") ;
496       fileName += (iDDL + kDDLOffset) ; 
497       fileName += ".ddl" ; 
498       buffer = new AliAltroBuffer(fileName.Data(), 1);
499       buffer->WriteDataHeader(kTRUE, kFALSE);  //Dummy;
500
501       prevDDL = iDDL;
502     }
503
504     // out of time range signal (?)
505     if (digit->GetTimeR() > GetRawFormatTimeMax() ) {
506       buffer->FillBuffer(digit->GetAmp());
507       buffer->FillBuffer(GetRawFormatTimeBins() );  // time bin
508       buffer->FillBuffer(3);          // bunch length      
509       buffer->WriteTrailer(3, relId[3], relId[2], module);  // trailer
510       
511     // calculate the time response function
512     } else {
513       Double_t energy = 0 ;
514       Int_t   module = relId[0];
515 //       Int_t   column = relId[3];
516 //       Int_t   row    = relId[2];
517       if ( digit->GetId() <= geom->GetNModules() *  geom->GetNCristalsInModule()) {
518 //      if(calib)
519 //        energy = digit->GetAmp()*calib->GetADCchannelEmc(module,column,row) + 
520 //          calib->GetADCpedestalEmc(module,column,row);
521 //      else
522 //        energy=digit->GetAmp()*digitizer->GetEMCchannel()+digitizer->GetEMCpedestal();
523 //       } 
524         energy=digit->GetAmp();
525       }
526       else {
527 //      energy = digit->GetAmp()*digitizer->GetCPVchannel()+digitizer->GetCPVpedestal();
528         energy = 0; // CPV raw data format is now know yet
529       }        
530       Bool_t lowgain = RawSampledResponse(digit->GetTimeR(), energy, adcValuesHigh, adcValuesLow) ; 
531       
532      if (lowgain) 
533         buffer->WriteChannel(relId[3], relId[2], module + fLowGainOffset, 
534                            GetRawFormatTimeBins(), adcValuesLow, kThreshold);
535       else 
536         buffer->WriteChannel(relId[3], relId[2], module, 
537                              GetRawFormatTimeBins(), adcValuesHigh, kThreshold);
538       
539     }
540   }
541   
542   // write real header and close last file
543   if (buffer) {
544     buffer->Flush();
545     buffer->WriteDataHeader(kFALSE, kFALSE);
546     delete buffer;
547   }
548   
549   loader->UnloadDigits();
550 }
551
552 //____________________________________________________________________________
553 void AliPHOS::Hits2SDigits()  
554
555 // create summable digits
556
557   AliPHOSSDigitizer phosDigitizer(fLoader->GetRunLoader()->GetFileName().Data()) ;
558   phosDigitizer.SetEventRange(0, -1) ; // do all the events
559   phosDigitizer.ExecuteTask("all") ; 
560 }
561
562 //____________________________________________________________________________
563 AliLoader* AliPHOS::MakeLoader(const char* topfoldername)
564 {
565 //different behaviour than standard (singleton getter)
566 // --> to be discussed and made eventually coherent
567  fLoader = new AliPHOSLoader(GetName(),topfoldername);
568  return fLoader;
569 }
570
571 //__________________________________________________________________
572 Double_t AliPHOS::RawResponseFunction(Double_t *x, Double_t *par) 
573 {
574   // Shape of the electronics raw reponse:
575   // It is a semi-gaussian, 2nd order Gamma function of the general form
576   // v(t) = n**n * Q * A**n / C *(t/tp)**n * exp(-n * t/tp) with 
577   // tp : peaking time par[0]
578   // n  : order of the function
579   // C  : integrating capacitor in the preamplifier
580   // A  : open loop gain of the preamplifier
581   // Q  : the total APD charge to be measured Q = C * energy
582   
583   Double_t signal ;
584   Double_t xx = x[0] - ( fgTimeTrigger + par[3] ) ; 
585
586   if (xx < 0 || xx > fgTimeMax) 
587     signal = 0. ;  
588   else { 
589     Double_t fac = par[0] * TMath::Power(fgOrder, fgOrder) * TMath::Power(par[1], fgOrder) / fgCapa ; 
590     signal = fac * par[2] * TMath::Power(xx / fgTimePeak, fgOrder) * TMath::Exp(-fgOrder * (xx / fgTimePeak)) ; 
591   }
592   return signal ;  
593 }
594
595 //__________________________________________________________________
596 Double_t AliPHOS::RawResponseFunctionMax(Double_t charge, Double_t gain) 
597 {
598   return ( charge * TMath::Power(fgOrder, fgOrder) * TMath::Power(gain, fgOrder) 
599      / ( fgCapa * TMath::Exp(fgOrder) ) );  
600
601 }
602
603 //__________________________________________________________________
604 Bool_t AliPHOS::RawSampledResponse(Double_t dtime, Double_t damp, Int_t * adcH, Int_t * adcL) const 
605 {
606   // for a start time dtime and an amplitude damp given by digit, 
607   // calculates the raw sampled response AliPHOS::RawResponseFunction
608
609   const Int_t kRawSignalOverflow = 0x3FF ; 
610   Bool_t lowGain = kFALSE ; 
611
612   TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 4);
613
614   for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
615     signalF.SetParameter(0, GetRawFormatHighCharge() ) ; 
616     signalF.SetParameter(1, GetRawFormatHighGain() ) ; 
617     signalF.SetParameter(2, damp) ; 
618     signalF.SetParameter(3, dtime) ; 
619     Double_t time = iTime * GetRawFormatTimeMax() / GetRawFormatTimeBins() ;
620     Double_t signal = signalF.Eval(time) ;     
621     if ( static_cast<Int_t>(signal+0.5) > kRawSignalOverflow ){  // larger than 10 bits 
622       signal = kRawSignalOverflow ;
623       lowGain = kTRUE ; 
624     }
625     adcH[iTime] =  static_cast<Int_t>(signal + 0.5) ;
626
627     signalF.SetParameter(0, GetRawFormatLowCharge() ) ;     
628     signalF.SetParameter(1, GetRawFormatLowGain() ) ; 
629     signal = signalF.Eval(time) ;  
630     if ( static_cast<Int_t>(signal+0.5) > kRawSignalOverflow)  // larger than 10 bits 
631       signal = kRawSignalOverflow ;
632     adcL[iTime] = static_cast<Int_t>(0.5 + signal ) ; 
633
634   }
635   return lowGain ; 
636 }
637
638 //____________________________________________________________________________
639 void AliPHOS::SetTreeAddress()
640
641   // Links Hits in the Tree to Hits array
642   TBranch *branch;
643   char branchname[20];
644   sprintf(branchname,"%s",GetName());
645   // Branch address for hit tree
646     TTree *treeH = TreeH();
647   if (treeH) {
648     branch = treeH->GetBranch(branchname);
649     if (branch) 
650      { 
651        if (fHits == 0x0) fHits= new TClonesArray("AliPHOSHit",1000);
652        //AliInfo(Form("<%s> Setting Hits Address",GetName()));
653        branch->SetAddress(&fHits);
654      }
655   }
656 }
657