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