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