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