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