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