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