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