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