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