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