]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSv1.cxx
6d221a3920bcbf51a3c14f427731542d1003a632
[u/mrichter/AliRoot.git] / PHOS / AliPHOSv1.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
16 /* $Id$ */
17
18 /* History of cvs commits:
19  *
20  * $Log$
21  * Revision 1.104  2005/05/28 14:19:05  schutz
22  * Compilation warnings fixed by T.P.
23  *
24  */
25
26 //_________________________________________________________________________
27 // Implementation version v1 of PHOS Manager class 
28 //---
29 //---
30 // Layout EMC + CPV  has name IHEP:
31 // Produces hits for CPV, cumulated hits
32 //---
33 //---
34 //*-- Author: Yves Schutz (SUBATECH)
35
36
37 // --- ROOT system ---
38 #include <TParticle.h>
39 #include <TVirtualMC.h>
40
41 // --- Standard library ---
42
43
44 // --- AliRoot header files ---
45 #include "AliPHOSCPVDigit.h"
46 #include "AliPHOSGeometry.h"
47 #include "AliPHOSHit.h"
48 #include "AliPHOSv1.h"
49 #include "AliRun.h"
50 #include "AliMC.h"
51
52 ClassImp(AliPHOSv1)
53
54 //____________________________________________________________________________
55 AliPHOSv1::AliPHOSv1():
56   fLightYieldMean(0.),
57   fIntrinsicPINEfficiency(0.),
58   fLightYieldAttenuation(0.),
59   fRecalibrationFactor(0.),
60   fElectronsPerGeV(0.),
61   fAPDGain(0.),
62   fLightFactor(0.),
63   fAPDFactor(0.)
64 {
65   //Def ctor.
66 }
67
68 //____________________________________________________________________________
69 AliPHOSv1::AliPHOSv1(const char *name, const char *title):
70   AliPHOSv0(name,title),
71   fLightYieldMean(0.),
72   fIntrinsicPINEfficiency(0.),
73   fLightYieldAttenuation(0.),
74   fRecalibrationFactor(0.),
75   fElectronsPerGeV(0.),
76   fAPDGain(0.),
77   fLightFactor(0.),
78   fAPDFactor(0.)
79 {
80   //
81   // We store hits :
82   //   - fHits (the "normal" one), which retains the hits associated with
83   //     the current primary particle being tracked
84   //     (this array is reset after each primary has been tracked).
85   //
86
87
88
89   // We do not want to save in TreeH the raw hits
90   // But save the cumulated hits instead (need to create the branch myself)
91   // It is put in the Digit Tree because the TreeH is filled after each primary
92   // and the TreeD at the end of the event (branch is set in FinishEvent() ). 
93   
94   fHits= new TClonesArray("AliPHOSHit",1000) ;
95   gAlice->GetMCApp()->AddHitList(fHits) ; 
96
97   fNhits = 0 ;
98
99   fIshunt     =  2 ; // All hits are associated with primary particles
100
101   //Photoelectron statistics:
102   // The light yield is a poissonian distribution of the number of
103   // photons created in the PbWo4 crystal, calculated using following formula
104   // NumberOfPhotons = EnergyLost * LightYieldMean* APDEfficiency *
105   //              exp (-LightYieldAttenuation * DistanceToPINdiodeFromTheHit);
106   // LightYieldMean is parameter calculated to be over 47000 photons per GeV
107   // APDEfficiency is 0.02655
108   // k_0 is 0.0045 from Valery Antonenko
109   // The number of electrons created in the APD is
110   // NumberOfElectrons = APDGain * LightYield
111   // The APD Gain is 300
112   fLightYieldMean = 47000;
113   fIntrinsicPINEfficiency = 0.02655 ; //APD= 0.1875/0.1271 * 0.018 (PIN)
114   fLightYieldAttenuation  = 0.0045 ; 
115   fRecalibrationFactor    = 13.418/ fLightYieldMean ;
116   fElectronsPerGeV        = 2.77e+8 ;
117   fAPDGain                = 300. ;
118   fLightFactor            = fLightYieldMean * fIntrinsicPINEfficiency ; 
119   fAPDFactor              = (fRecalibrationFactor/100.) * fAPDGain ;   
120 }
121
122 AliPHOSv1::AliPHOSv1(AliPHOSv1 & phos) :
123   AliPHOSv0(phos),
124   fLightYieldMean(0.),
125   fIntrinsicPINEfficiency(0.),
126   fLightYieldAttenuation(0.),
127   fRecalibrationFactor(0.),
128   fElectronsPerGeV(0.),
129   fAPDGain(0.),
130   fLightFactor(0.),
131   fAPDFactor(0.)
132 {
133   //Copy ctor. Can be wrong.
134   phos.Copy(*this) ; 
135 }
136
137
138 //____________________________________________________________________________
139 AliPHOSv1::~AliPHOSv1()
140 {
141   // dtor
142  if ( fHits) {
143     fHits->Delete() ; 
144     delete fHits ;
145     fHits = 0 ; 
146  }
147 }
148
149 //____________________________________________________________________________
150 void AliPHOSv1::Copy(TObject & base)const
151 {
152   TObject::Copy(base) ; 
153   AliPHOSv0::Copy(base) ;
154   AliPHOSv1 &phos = static_cast<AliPHOSv1 &>(base); 
155   phos.fLightYieldMean         = fLightYieldMean ; 
156   phos.fIntrinsicPINEfficiency = fIntrinsicPINEfficiency ; 
157   phos.fLightYieldAttenuation  = fLightYieldAttenuation ; 
158   phos.fRecalibrationFactor    = fRecalibrationFactor ; 
159   phos.fElectronsPerGeV        = fElectronsPerGeV ; 
160   phos.fAPDGain                = fAPDGain ; 
161   phos.fLightFactor            = fLightFactor ; 
162   phos.fAPDFactor              = fAPDFactor ; 
163 }
164
165 //____________________________________________________________________________
166 void AliPHOSv1::AddHit(Int_t shunt, Int_t primary, Int_t Id, Float_t * hits)
167 {
168   // Add a hit to the hit list.
169   // A PHOS hit is the sum of all hits in a single crystal from one primary and within some time gate
170
171   Int_t hitCounter ;
172   AliPHOSHit *newHit ;
173   AliPHOSHit *curHit ;
174   Bool_t deja = kFALSE ;
175   AliPHOSGeometry * geom = GetGeometry() ; 
176
177   newHit = new AliPHOSHit(shunt, primary, Id, hits) ;
178
179   for ( hitCounter = fNhits-1 ; hitCounter >= 0 && !deja ; hitCounter-- ) {
180     curHit = dynamic_cast<AliPHOSHit*>((*fHits)[hitCounter]) ;
181     if(curHit->GetPrimary() != primary) break ; 
182            // We add hits with the same primary, while GEANT treats primaries succesively 
183     if( *curHit == *newHit ) {
184       *curHit + *newHit ;
185       deja = kTRUE ;
186     }
187   }
188          
189   if ( !deja ) {
190     new((*fHits)[fNhits]) AliPHOSHit(*newHit) ;
191     // get the block Id number
192     Int_t relid[4] ;
193     geom->AbsToRelNumbering(Id, relid) ;
194
195     fNhits++ ;
196   }
197   
198   delete newHit;
199 }
200
201 //____________________________________________________________________________
202 void AliPHOSv1::FinishPrimary() 
203 {
204   // called at the end of each track (primary) by AliRun
205   // hits are reset for each new track
206   // accumulate the total hit-multiplicity
207
208 }
209
210 //____________________________________________________________________________
211 void AliPHOSv1::FinishEvent() 
212 {
213   // called at the end of each event by AliRun
214   // accumulate the hit-multiplicity and total energy per block 
215   // if the values have been updated check it
216   
217   AliDetector::FinishEvent(); 
218 }
219 //____________________________________________________________________________
220 void AliPHOSv1::StepManager(void)
221 {
222    // Accumulates hits as long as the track stays in a single crystal or CPV gas Cell
223
224   Int_t          relid[4] ;           // (box, layer, row, column) indices
225   Int_t          absid    ;           // absolute cell ID number
226   Float_t        xyzte[5]={-1000.,-1000.,-1000.,0.,0.}  ; // position wrt MRS, time and energy deposited
227   TLorentzVector pos      ;           // Lorentz vector of the track current position
228   Int_t          copy     ;
229
230   TString name      =  GetGeometry()->GetName() ; 
231
232   Int_t moduleNumber ;
233   
234   static Int_t idPCPQ = gMC->VolId("PCPQ");
235   if( gMC->CurrentVolID(copy) == idPCPQ &&
236       (gMC->IsTrackEntering() ) &&
237       gMC->TrackCharge() != 0) {      
238     
239     gMC -> TrackPosition(pos);
240     
241     Float_t xyzm[3], xyzd[3] ;
242     Int_t i;
243     for (i=0; i<3; i++) xyzm[i] = pos[i];
244     gMC -> Gmtod (xyzm, xyzd, 1);    // transform coordinate from master to daughter system
245     
246     Float_t        xyd[3]={0,0,0}   ;   //local position of the entering
247     xyd[0]  = xyzd[0];
248     xyd[1]  =-xyzd[2];
249     xyd[2]  =-xyzd[1];
250     
251     // Current momentum of the hit's track in the local ref. system
252     TLorentzVector pmom     ;        //momentum of the particle initiated hit
253     gMC -> TrackMomentum(pmom);
254     Float_t pm[3], pd[3];
255     for (i=0; i<3; i++)  
256       pm[i]   = pmom[i];
257     
258     gMC -> Gmtod (pm, pd, 2);        // transform 3-momentum from master to daughter system
259     pmom[0] = pd[0];
260     pmom[1] =-pd[1];
261     pmom[2] =-pd[2];
262
263     // Digitize the current CPV hit:
264     
265     // 1. find pad response and    
266     gMC->CurrentVolOffID(3,moduleNumber);
267     moduleNumber--;
268     
269     TClonesArray *cpvDigits = new TClonesArray("AliPHOSCPVDigit",0);   // array of digits for current hit
270     CPVDigitize(pmom,xyd,cpvDigits);
271       
272     Float_t xmean = 0;
273     Float_t zmean = 0;
274     Float_t qsum  = 0;
275     Int_t   idigit,ndigits;
276     
277     // 2. go through the current digit list and sum digits in pads
278     
279     ndigits = cpvDigits->GetEntriesFast();
280     for (idigit=0; idigit<ndigits-1; idigit++) {
281       AliPHOSCPVDigit  *cpvDigit1 = dynamic_cast<AliPHOSCPVDigit*>(cpvDigits->UncheckedAt(idigit));
282       Float_t x1 = cpvDigit1->GetXpad() ;
283       Float_t z1 = cpvDigit1->GetYpad() ;
284       for (Int_t jdigit=idigit+1; jdigit<ndigits; jdigit++) {
285         AliPHOSCPVDigit  *cpvDigit2 = dynamic_cast<AliPHOSCPVDigit*>(cpvDigits->UncheckedAt(jdigit));
286         Float_t x2 = cpvDigit2->GetXpad() ;
287         Float_t z2 = cpvDigit2->GetYpad() ;
288         if (x1==x2 && z1==z2) {
289           Float_t qsum = cpvDigit1->GetQpad() + cpvDigit2->GetQpad() ;
290           cpvDigit2->SetQpad(qsum) ;
291           cpvDigits->RemoveAt(idigit) ;
292         }
293       }
294     }
295     cpvDigits->Compress() ;
296     
297     // 3. add digits to temporary hit list fTmpHits
298     
299     ndigits = cpvDigits->GetEntriesFast();
300     for (idigit=0; idigit<ndigits; idigit++) {
301       AliPHOSCPVDigit  *cpvDigit = dynamic_cast<AliPHOSCPVDigit*>(cpvDigits->UncheckedAt(idigit));
302       relid[0] = moduleNumber + 1 ;                             // CPV (or PHOS) module number
303       relid[1] =-1 ;                                            // means CPV
304       relid[2] = cpvDigit->GetXpad() ;                          // column number of a pad
305       relid[3] = cpvDigit->GetYpad() ;                          // row    number of a pad
306       
307       // get the absolute Id number
308       GetGeometry()->RelToAbsNumbering(relid, absid) ; 
309       
310       // add current digit to the temporary hit list
311
312       xyzte[3] = gMC->TrackTime() ;
313       xyzte[4] = cpvDigit->GetQpad() ;                          // amplitude in a pad
314
315       Int_t primary  =  gAlice->GetMCApp()->GetPrimary( gAlice->GetMCApp()->GetCurrentTrackNumber() ); 
316       AddHit(fIshunt, primary, absid, xyzte);  
317       
318       if (cpvDigit->GetQpad() > 0.02) {
319         xmean += cpvDigit->GetQpad() * (cpvDigit->GetXpad() + 0.5);
320         zmean += cpvDigit->GetQpad() * (cpvDigit->GetYpad() + 0.5);
321         qsum  += cpvDigit->GetQpad();
322       }
323     }
324     if (cpvDigits) {
325       cpvDigits->Delete();
326       delete cpvDigits;
327       cpvDigits=0;
328     }
329   }
330
331  
332   static Int_t idPXTL = gMC->VolId("PXTL");  
333   if(gMC->CurrentVolID(copy) == idPXTL ) { //  We are inside a PBWO crystal
334
335     gMC->TrackPosition(pos) ;
336     xyzte[0] = pos[0] ;
337     xyzte[1] = pos[1] ;
338     xyzte[2] = pos[2] ;
339
340     Float_t global[3], local[3] ;
341     global[0] = pos[0] ;
342     global[1] = pos[1] ;
343     global[2] = pos[2] ;
344     Float_t lostenergy = gMC->Edep(); 
345     
346     //Put in the TreeK particle entering PHOS and all its parents
347     if ( gMC->IsTrackEntering() ){
348       Float_t xyzd[3] ;
349       gMC -> Gmtod (xyzte, xyzd, 1);    // transform coordinate from master to daughter system    
350       if (xyzd[1] < -GetGeometry()->GetCrystalSize(1)/2.+0.1){   //Entered close to forward surface  
351         Int_t parent = gAlice->GetMCApp()->GetCurrentTrackNumber() ; 
352         TParticle * part = gAlice->GetMCApp()->Particle(parent) ; 
353         Float_t vert[3],vertd[3] ;
354         vert[0]=part->Vx() ;
355         vert[1]=part->Vy() ;
356         vert[2]=part->Vz() ;
357         gMC -> Gmtod (vert, vertd, 1);    // transform coordinate from master to daughter system
358         if(vertd[1]<-GetGeometry()->GetCrystalSize(1)/2.-0.1){ //Particle is created in foront of PHOS 
359                                                                //0.1 to get rid of numerical errors 
360           part->SetBit(kKeepBit);
361           while ( parent != -1 ) {
362             part = gAlice->GetMCApp()->Particle(parent) ; 
363             part->SetBit(kKeepBit);
364             parent = part->GetFirstMother() ; 
365           }
366         }
367       }
368     }
369     if ( lostenergy != 0 ) {  // Track is inside the crystal and deposits some energy 
370       xyzte[3] = gMC->TrackTime() ;     
371       
372       gMC->CurrentVolOffID(10, moduleNumber) ; // get the PHOS module number ;
373       
374       Int_t strip ;
375       gMC->CurrentVolOffID(3, strip);
376       Int_t cell ;
377       gMC->CurrentVolOffID(2, cell);
378       
379       Int_t row = 1 + GetGeometry()->GetNZ() - strip % GetGeometry()->GetNZ() ;
380       Int_t col = (Int_t) TMath::Ceil((Double_t) strip/GetGeometry()->GetNZ()) -1 ;
381       
382       absid = (moduleNumber-1)*GetGeometry()->GetNCristalsInModule() + 
383         row + (col*GetGeometry()->GetEMCAGeometry()->GetNCellsInStrip() + cell-1)*GetGeometry()->GetNZ() ;
384       
385       gMC->Gmtod(global, local, 1) ;
386       
387       //Calculates the light yield, the number of photons produced in the
388       //crystal 
389       Float_t lightYield = gRandom->Poisson(fLightFactor * lostenergy *
390                                             exp(-fLightYieldAttenuation *
391                                                 (local[1]+GetGeometry()->GetCrystalSize(1)/2.0 ))
392                                             ) ;
393
394       //Calculates de energy deposited in the crystal  
395       xyzte[4] = fAPDFactor * lightYield  ;
396       
397       Int_t primary ;
398       if(fIshunt == 2){
399         primary = gAlice->GetMCApp()->GetCurrentTrackNumber() ;
400         TParticle * part = gAlice->GetMCApp()->Particle(primary) ;
401         while ( !part->TestBit(kKeepBit) ) {
402           primary = part->GetFirstMother() ;
403           if(primary == -1){        
404             primary  =  gAlice->GetMCApp()->GetPrimary( gAlice->GetMCApp()->GetCurrentTrackNumber() ); 
405             break ; //there is a possibility that particle passed e.g. thermal isulator and hits a side 
406           //surface of the crystal. In this case it may have no primary at all. 
407           //We can not easily separate this case from the case when this is part of the shower, 
408           //developed in the neighboring crystal.
409           }
410           part = gAlice->GetMCApp()->Particle(primary) ;
411         }
412       }
413       else
414         primary  =  gAlice->GetMCApp()->GetPrimary( gAlice->GetMCApp()->GetCurrentTrackNumber() ); 
415
416       
417       
418       // add current hit to the hit list
419       // Info("StepManager","%d %d", primary, tracknumber) ; 
420       AddHit(fIshunt, primary, absid, xyzte);
421         
422     } // there is deposited energy
423   } // we are inside a PHOS Xtal
424   
425 }
426
427 //____________________________________________________________________________
428 void AliPHOSv1::CPVDigitize (TLorentzVector p, Float_t *zxhit, TClonesArray *cpvDigits)
429 {
430   // ------------------------------------------------------------------------
431   // Digitize one CPV hit:
432   // On input take exact 4-momentum p and position zxhit of the hit,
433   // find the pad response around this hit and
434   // put the amplitudes in the pads into array digits
435   //
436   // Author: Yuri Kharlov (after Serguei Sadovsky)
437   // 2 October 2000
438   // ------------------------------------------------------------------------
439
440   const Float_t kCelWr  = GetGeometry()->GetPadSizePhi()/2;  // Distance between wires (2 wires above 1 pad)
441   const Float_t kDetR   = 0.1;     // Relative energy fluctuation in track for 100 e-
442   const Float_t kdEdx   = 4.0;     // Average energy loss in CPV;
443   const Int_t   kNgamz  = 5;       // Ionization size in Z
444   const Int_t   kNgamx  = 9;       // Ionization size in Phi
445   const Float_t kNoise = 0.03;    // charge noise in one pad
446
447   Float_t rnor1,rnor2;
448
449   // Just a reminder on axes notation in the CPV module:
450   // axis Z goes along the beam
451   // axis X goes across the beam in the module plane
452   // axis Y is a normal to the module plane showing from the IP
453
454   Float_t hitX  = zxhit[0];
455   Float_t hitZ  =-zxhit[1];
456   Float_t pX    = p.Px();
457   Float_t pZ    =-p.Pz();
458   Float_t pNorm = p.Py();
459   Float_t eloss = kdEdx;
460
461 //Info("CPVDigitize", "YVK : %f %f | %f %f %d", hitX, hitZ, pX, pZ, pNorm) ;
462
463   Float_t dZY   = pZ/pNorm * GetGeometry()->GetCPVGasThickness();
464   Float_t dXY   = pX/pNorm * GetGeometry()->GetCPVGasThickness();
465   gRandom->Rannor(rnor1,rnor2);
466   eloss *= (1 + kDetR*rnor1) *
467            TMath::Sqrt((1 + ( pow(dZY,2) + pow(dXY,2) ) / pow(GetGeometry()->GetCPVGasThickness(),2)));
468   Float_t zhit1 = hitZ + GetGeometry()->GetCPVActiveSize(1)/2 - dZY/2;
469   Float_t xhit1 = hitX + GetGeometry()->GetCPVActiveSize(0)/2 - dXY/2;
470   Float_t zhit2 = zhit1 + dZY;
471   Float_t xhit2 = xhit1 + dXY;
472
473   Int_t   iwht1 = (Int_t) (xhit1 / kCelWr);           // wire (x) coordinate "in"
474   Int_t   iwht2 = (Int_t) (xhit2 / kCelWr);           // wire (x) coordinate "out"
475
476   Int_t   nIter;
477   Float_t zxe[3][5];
478   if (iwht1==iwht2) {                      // incline 1-wire hit
479     nIter = 2;
480     zxe[0][0] = (zhit1 + zhit2 - dZY*0.57735) / 2;
481     zxe[1][0] = (iwht1 + 0.5) * kCelWr;
482     zxe[2][0] =  eloss/2;
483     zxe[0][1] = (zhit1 + zhit2 + dZY*0.57735) / 2;
484     zxe[1][1] = (iwht1 + 0.5) * kCelWr;
485     zxe[2][1] =  eloss/2;
486   }
487   else if (TMath::Abs(iwht1-iwht2) != 1) { // incline 3-wire hit
488     nIter = 3;
489     Int_t iwht3 = (iwht1 + iwht2) / 2;
490     Float_t xwht1 = (iwht1 + 0.5) * kCelWr; // wire 1
491     Float_t xwht2 = (iwht2 + 0.5) * kCelWr; // wire 2
492     Float_t xwht3 = (iwht3 + 0.5) * kCelWr; // wire 3
493     Float_t xwr13 = (xwht1 + xwht3) / 2;   // center 13
494     Float_t xwr23 = (xwht2 + xwht3) / 2;   // center 23
495     Float_t dxw1  = xhit1 - xwr13;
496     Float_t dxw2  = xhit2 - xwr23;
497     Float_t egm1  = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
498     Float_t egm2  = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
499     Float_t egm3  =           kCelWr / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
500     zxe[0][0] = (dXY*(xwr13-xwht1)/dXY + zhit1 + zhit1) / 2;
501     zxe[1][0] =  xwht1;
502     zxe[2][0] =  eloss * egm1;
503     zxe[0][1] = (dXY*(xwr23-xwht1)/dXY + zhit1 + zhit2) / 2;
504     zxe[1][1] =  xwht2;
505     zxe[2][1] =  eloss * egm2;
506     zxe[0][2] =  dXY*(xwht3-xwht1)/dXY + zhit1;
507     zxe[1][2] =  xwht3;
508     zxe[2][2] =  eloss * egm3;
509   }
510   else {                                   // incline 2-wire hit
511     nIter = 2;
512     Float_t xwht1 = (iwht1 + 0.5) * kCelWr;
513     Float_t xwht2 = (iwht2 + 0.5) * kCelWr;
514     Float_t xwr12 = (xwht1 + xwht2) / 2;
515     Float_t dxw1  = xhit1 - xwr12;
516     Float_t dxw2  = xhit2 - xwr12;
517     Float_t egm1  = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) );
518     Float_t egm2  = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) );
519     zxe[0][0] = (zhit1 + zhit2 - dZY*egm1) / 2;
520     zxe[1][0] =  xwht1;
521     zxe[2][0] =  eloss * egm1;
522     zxe[0][1] = (zhit1 + zhit2 + dZY*egm2) / 2;
523     zxe[1][1] =  xwht2;
524     zxe[2][1] =  eloss * egm2;
525   }
526
527   // Finite size of ionization region
528
529   Int_t nCellZ  = GetGeometry()->GetNumberOfCPVPadsZ();
530   Int_t nCellX  = GetGeometry()->GetNumberOfCPVPadsPhi();
531   Int_t nz3     = (kNgamz+1)/2;
532   Int_t nx3     = (kNgamx+1)/2;
533   cpvDigits->Expand(nIter*kNgamx*kNgamz);
534   TClonesArray &ldigits = *(static_cast<TClonesArray *>(cpvDigits));
535
536   for (Int_t iter=0; iter<nIter; iter++) {
537
538     Float_t zhit = zxe[0][iter];
539     Float_t xhit = zxe[1][iter];
540     Float_t qhit = zxe[2][iter];
541     Float_t zcell = zhit / GetGeometry()->GetPadSizeZ();
542     Float_t xcell = xhit / GetGeometry()->GetPadSizePhi();
543     if ( zcell<=0      || xcell<=0 ||
544          zcell>=nCellZ || xcell>=nCellX) return;
545     Int_t izcell = (Int_t) zcell;
546     Int_t ixcell = (Int_t) xcell;
547     Float_t zc = zcell - izcell - 0.5;
548     Float_t xc = xcell - ixcell - 0.5;
549     for (Int_t iz=1; iz<=kNgamz; iz++) {
550       Int_t kzg = izcell + iz - nz3;
551       if (kzg<=0 || kzg>nCellZ) continue;
552       Float_t zg = (Float_t)(iz-nz3) - zc;
553       for (Int_t ix=1; ix<=kNgamx; ix++) {
554         Int_t kxg = ixcell + ix - nx3;
555         if (kxg<=0 || kxg>nCellX) continue;
556         Float_t xg = (Float_t)(ix-nx3) - xc;
557         
558         // Now calculate pad response
559         Float_t qpad = CPVPadResponseFunction(qhit,zg,xg);
560         qpad += kNoise*rnor2;
561         if (qpad<0) continue;
562         
563         // Fill the array with pad response ID and amplitude
564         new(ldigits[cpvDigits->GetEntriesFast()]) AliPHOSCPVDigit(kxg,kzg,qpad);
565       }
566     }
567   }
568 }
569
570 //____________________________________________________________________________
571 Float_t AliPHOSv1::CPVPadResponseFunction(Float_t qhit, Float_t zhit, Float_t xhit) {
572   // ------------------------------------------------------------------------
573   // Calculate the amplitude in one CPV pad using the
574   // cumulative pad response function
575   // Author: Yuri Kharlov (after Serguei Sadovski)
576   // 3 October 2000
577   // ------------------------------------------------------------------------
578
579   Double_t dz = GetGeometry()->GetPadSizeZ()   / 2;
580   Double_t dx = GetGeometry()->GetPadSizePhi() / 2;
581   Double_t z  = zhit * GetGeometry()->GetPadSizeZ();
582   Double_t x  = xhit * GetGeometry()->GetPadSizePhi();
583   Double_t amplitude = qhit *
584     (CPVCumulPadResponse(z+dz,x+dx) - CPVCumulPadResponse(z+dz,x-dx) -
585      CPVCumulPadResponse(z-dz,x+dx) + CPVCumulPadResponse(z-dz,x-dx));
586   return (Float_t)amplitude;
587 }
588
589 //____________________________________________________________________________
590 Double_t AliPHOSv1::CPVCumulPadResponse(Double_t x, Double_t y) {
591   // ------------------------------------------------------------------------
592   // Cumulative pad response function
593   // It includes several terms from the CF decomposition in electrostatics
594   // Note: this cumulative function is wrong since omits some terms
595   //       but the cell amplitude obtained with it is correct because
596   //       these omitting terms cancel
597   // Author: Yuri Kharlov (after Serguei Sadovski)
598   // 3 October 2000
599   // ------------------------------------------------------------------------
600
601   const Double_t kA=1.0;
602   const Double_t kB=0.7;
603
604   Double_t r2       = x*x + y*y;
605   Double_t xy       = x*y;
606   Double_t cumulPRF = 0;
607   for (Int_t i=0; i<=4; i++) {
608     Double_t b1 = (2*i + 1) * kB;
609     cumulPRF += TMath::Power(-1,i) * TMath::ATan( xy / (b1*TMath::Sqrt(b1*b1 + r2)) );
610   }
611   cumulPRF *= kA/(2*TMath::Pi());
612   return cumulPRF;
613 }
614