]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSv1.cxx
#include "AliHeader.h" needed
[u/mrichter/AliRoot.git] / PHOS / AliPHOSv1.cxx
CommitLineData
7587f5a5 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
b2a60966 16/* $Id$ */
5f20d3fb 17
7587f5a5 18//_________________________________________________________________________
5f20d3fb 19// Implementation version v1 of PHOS Manager class
a3dfe79c 20//---
21// Layout EMC + PPSD has name GPS2:
ed4205d8 22// Produces cumulated hits
a3dfe79c 23//---
24// Layout EMC + CPV has name IHEP:
ed4205d8 25// Produces hits for CPV, cumulated hits
26//---
27// Layout EMC + CPV + PPSD has name GPS:
28// Produces hits for CPV, cumulated hits
29//---
5f20d3fb 30//*-- Author: Yves Schutz (SUBATECH)
b2a60966 31
7587f5a5 32
33// --- ROOT system ---
bea63bea 34
35#include "TBRIK.h"
36#include "TNode.h"
7587f5a5 37#include "TRandom.h"
94de3818 38#include "TTree.h"
7587f5a5 39
5f20d3fb 40
7587f5a5 41// --- Standard library ---
42
de9ec31b 43#include <stdio.h>
44#include <string.h>
45#include <stdlib.h>
46#include <strstream.h>
7587f5a5 47
48// --- AliRoot header files ---
49
50#include "AliPHOSv1.h"
51#include "AliPHOSHit.h"
97cee223 52#include "AliPHOSCPVDigit.h"
7587f5a5 53#include "AliRun.h"
54#include "AliConst.h"
94de3818 55#include "AliMC.h"
97cee223 56#include "AliPHOSGeometry.h"
7587f5a5 57
58ClassImp(AliPHOSv1)
59
bea63bea 60//____________________________________________________________________________
02ab1add 61AliPHOSv1::AliPHOSv1():
62AliPHOSv0()
bea63bea 63{
5f20d3fb 64 // ctor
bea63bea 65}
66
7587f5a5 67//____________________________________________________________________________
68AliPHOSv1::AliPHOSv1(const char *name, const char *title):
e04976bd 69AliPHOSv0(name,title)
7587f5a5 70{
5f20d3fb 71 // ctor : title is used to identify the layout
fa412d9b 72 // GPS2 = 5 modules (EMC + PPSD)
73 // IHEP = 5 modules (EMC + CPV )
ed4205d8 74 // MIXT = 4 modules (EMC + CPV ) and 1 module (EMC + PPSD)
5f20d3fb 75 //
ed4205d8 76 // We store hits :
5f20d3fb 77 // - fHits (the "normal" one), which retains the hits associated with
78 // the current primary particle being tracked
79 // (this array is reset after each primary has been tracked).
80 //
fa412d9b 81
037cc66d 82
5f20d3fb 83
84 // We do not want to save in TreeH the raw hits
85 // But save the cumulated hits instead (need to create the branch myself)
86 // It is put in the Digit Tree because the TreeH is filled after each primary
fa412d9b 87 // and the TreeD at the end of the event (branch is set in FinishEvent() ).
5f20d3fb 88
ed4205d8 89 fHits= new TClonesArray("AliPHOSHit",1000) ;
5f20d3fb 90
ed4205d8 91 fNhits = 0 ;
5f20d3fb 92
5f20d3fb 93 fIshunt = 1 ; // All hits are associated with primary particles
fa412d9b 94
5f20d3fb 95}
96
7587f5a5 97//____________________________________________________________________________
bea63bea 98AliPHOSv1::~AliPHOSv1()
b2a60966 99{
bea63bea 100 // dtor
5f20d3fb 101
ed4205d8 102 if ( fHits) {
103 fHits->Delete() ;
104 delete fHits ;
105 fHits = 0 ;
8dfa469d 106 }
7587f5a5 107}
108
7587f5a5 109//____________________________________________________________________________
b37750a6 110void AliPHOSv1::AddHit(Int_t shunt, Int_t primary, Int_t tracknumber, Int_t Id, Float_t * hits)
bea63bea 111{
112 // Add a hit to the hit list.
5f20d3fb 113 // A PHOS hit is the sum of all hits in a single crystal
114 // or in a single PPSD gas cell
bea63bea 115
5f20d3fb 116 Int_t hitCounter ;
bea63bea 117 AliPHOSHit *newHit ;
5f20d3fb 118 AliPHOSHit *curHit ;
119 Bool_t deja = kFALSE ;
bea63bea 120
b37750a6 121 newHit = new AliPHOSHit(shunt, primary, tracknumber, Id, hits) ;
bea63bea 122
7854a24a 123 for ( hitCounter = fNhits-1 ; hitCounter >= 0 && !deja ; hitCounter-- ) {
ed4205d8 124 curHit = (AliPHOSHit*) (*fHits)[hitCounter] ;
4e669ac8 125 if(curHit->GetPrimary() != primary) break ; // We add hits with the same primary, while GEANT treats primaries succesively
ed4205d8 126 if( *curHit == *newHit ) {
127 *curHit = *curHit + *newHit ;
128 deja = kTRUE ;
5f20d3fb 129 }
130 }
131
132 if ( !deja ) {
ed4205d8 133 new((*fHits)[fNhits]) AliPHOSHit(*newHit) ;
134 fNhits++ ;
5f20d3fb 135 }
136
bea63bea 137 delete newHit;
bea63bea 138}
139
5f20d3fb 140//____________________________________________________________________________
7587f5a5 141void AliPHOSv1::StepManager(void)
142{
b2a60966 143 // Accumulates hits as long as the track stays in a single crystal or PPSD gas Cell
b2a60966 144
4f5bbbd4 145 Int_t relid[4] ; // (box, layer, row, column) indices
146 Int_t absid ; // absolute cell ID number
2168f43b 147 Float_t xyze[4]={0,0,0,0} ; // position wrt MRS and energy deposited
4f5bbbd4 148 TLorentzVector pos ; // Lorentz vector of the track current position
fa412d9b 149 Int_t copy ;
7587f5a5 150
bea63bea 151 Int_t tracknumber = gAlice->CurrentTrack() ;
fa412d9b 152 Int_t primary = gAlice->GetPrimary( gAlice->CurrentTrack() );
153 TString name = fGeom->GetName() ;
037cc66d 154
155
ed4205d8 156 if ( name == "GPS2" || name == "MIXT" ) { // ======> CPV is a GPS' PPSD
fa412d9b 157
dc999bc0 158 if( gMC->CurrentVolID(copy) == gMC->VolId("PPCE") ) // We are inside a gas cell
7587f5a5 159 {
160 gMC->TrackPosition(pos) ;
161 xyze[0] = pos[0] ;
162 xyze[1] = pos[1] ;
163 xyze[2] = pos[2] ;
bea63bea 164 xyze[3] = gMC->Edep() ;
7587f5a5 165
b37750a6 166 if ( xyze[3] != 0 ) { // there is deposited energy
7587f5a5 167 gMC->CurrentVolOffID(5, relid[0]) ; // get the PHOS Module number
fad3e5b9 168 if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(5),"PHO1") == 0 ){
ed4205d8 169 relid[0] += fGeom->GetNModules() - fGeom->GetNPPSDModules();
fad3e5b9 170 }
7587f5a5 171 gMC->CurrentVolOffID(3, relid[1]) ; // get the Micromegas Module number
ed4205d8 172 // 1-> fGeom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() upper
173 // > fGeom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() lower
7587f5a5 174 gMC->CurrentVolOffID(1, relid[2]) ; // get the row number of the cell
175 gMC->CurrentVolID(relid[3]) ; // get the column number
176
177 // get the absolute Id number
178
bea63bea 179 fGeom->RelToAbsNumbering(relid, absid) ;
7587f5a5 180
bea63bea 181 // add current hit to the hit list
b37750a6 182 AddHit(fIshunt, primary, tracknumber, absid, xyze);
037cc66d 183
7587f5a5 184
185 } // there is deposited energy
fa412d9b 186 } // We are inside the gas of the CPV
187 } // GPS2 configuration
188
ed4205d8 189 if ( name == "IHEP" || name == "MIXT" ) { // ======> CPV is a IHEP's one
fa412d9b 190
191 // Yuri Kharlov, 28 September 2000
192
97cee223 193 if( gMC->CurrentVolID(copy) == gMC->VolId("PCPQ") &&
b37750a6 194 (gMC->IsTrackEntering() ) &&
037cc66d 195 gMC->TrackCharge() != 0) {
fa412d9b 196
b37750a6 197 gMC -> TrackPosition(pos);
198 Float_t xyzm[3], xyzd[3] ;
199 Int_t i;
200 for (i=0; i<3; i++) xyzm[i] = pos[i];
201 gMC -> Gmtod (xyzm, xyzd, 1); // transform coordinate from master to daughter system
202
203 Float_t xyd[3]={0,0,0} ; //local posiiton of the entering
204 xyd[0] = xyzd[0];
205 xyd[1] =-xyzd[1];
206 xyd[2] =-xyzd[2];
207
208
209 // Current momentum of the hit's track in the local ref. system
210 TLorentzVector pmom ; //momentum of the particle initiated hit
211 gMC -> TrackMomentum(pmom);
212 Float_t pm[3], pd[3];
213 for (i=0; i<3; i++) pm[i] = pmom[i];
214 gMC -> Gmtod (pm, pd, 2); // transform 3-momentum from master to daughter system
215 pmom[0] = pd[0];
216 pmom[1] =-pd[1];
217 pmom[2] =-pd[2];
218
037cc66d 219 // Digitize the current CPV hit:
220
221 // 1. find pad response and
fa412d9b 222
a3dfe79c 223 Int_t moduleNumber;
224 gMC->CurrentVolOffID(3,moduleNumber);
225 moduleNumber--;
fa412d9b 226
fa412d9b 227
3d402178 228 TClonesArray *cpvDigits = new TClonesArray("AliPHOSCPVDigit",0); // array of digits for current hit
a3dfe79c 229 CPVDigitize(pmom,xyd,moduleNumber,cpvDigits);
fa412d9b 230
231 Float_t xmean = 0;
232 Float_t zmean = 0;
233 Float_t qsum = 0;
cd461ab8 234 Int_t idigit,ndigits;
fa412d9b 235
236 // 2. go through the current digit list and sum digits in pads
237
238 ndigits = cpvDigits->GetEntriesFast();
cd461ab8 239 for (idigit=0; idigit<ndigits-1; idigit++) {
3d402178 240 AliPHOSCPVDigit *cpvDigit1 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit);
fa412d9b 241 Float_t x1 = cpvDigit1->GetXpad() ;
242 Float_t z1 = cpvDigit1->GetYpad() ;
243 for (Int_t jdigit=idigit+1; jdigit<ndigits; jdigit++) {
3d402178 244 AliPHOSCPVDigit *cpvDigit2 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(jdigit);
fa412d9b 245 Float_t x2 = cpvDigit2->GetXpad() ;
246 Float_t z2 = cpvDigit2->GetYpad() ;
247 if (x1==x2 && z1==z2) {
248 Float_t qsum = cpvDigit1->GetQpad() + cpvDigit2->GetQpad() ;
249 cpvDigit2->SetQpad(qsum) ;
250 cpvDigits->RemoveAt(idigit) ;
251 }
252 }
253 }
254 cpvDigits->Compress() ;
255
256 // 3. add digits to temporary hit list fTmpHits
257
258 ndigits = cpvDigits->GetEntriesFast();
cd461ab8 259 for (idigit=0; idigit<ndigits; idigit++) {
3d402178 260 AliPHOSCPVDigit *cpvDigit = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit);
a3dfe79c 261 relid[0] = moduleNumber + 1 ; // CPV (or PHOS) module number
fa412d9b 262 relid[1] =-1 ; // means CPV
263 relid[2] = cpvDigit->GetXpad() ; // column number of a pad
264 relid[3] = cpvDigit->GetYpad() ; // row number of a pad
265
266 // get the absolute Id number
267 fGeom->RelToAbsNumbering(relid, absid) ;
268
269 // add current digit to the temporary hit list
270 xyze[0] = 0. ;
271 xyze[1] = 0. ;
272 xyze[2] = 0. ;
273 xyze[3] = cpvDigit->GetQpad() ; // amplitude in a pad
274 primary = -1; // No need in primary for CPV
b37750a6 275 AddHit(fIshunt, primary, tracknumber, absid, xyze);
fa412d9b 276
277 if (cpvDigit->GetQpad() > 0.02) {
278 xmean += cpvDigit->GetQpad() * (cpvDigit->GetXpad() + 0.5);
279 zmean += cpvDigit->GetQpad() * (cpvDigit->GetYpad() + 0.5);
280 qsum += cpvDigit->GetQpad();
281 }
282 }
283 delete cpvDigits;
284 }
285 } // end of IHEP configuration
7587f5a5 286
037cc66d 287
fa412d9b 288 if(gMC->CurrentVolID(copy) == gMC->VolId("PXTL") ) { // We are inside a PBWO crystal
289 gMC->TrackPosition(pos) ;
290 xyze[0] = pos[0] ;
291 xyze[1] = pos[1] ;
292 xyze[2] = pos[2] ;
293 xyze[3] = gMC->Edep() ;
ed4205d8 294
037cc66d 295
b37750a6 296 if ( xyze[3] != 0 ) { // Track is inside the crystal and deposits some energy
ed4205d8 297
fa412d9b 298 gMC->CurrentVolOffID(10, relid[0]) ; // get the PHOS module number ;
037cc66d 299
fad3e5b9 300 if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(10),"PHO1") == 0 )
301 relid[0] += fGeom->GetNModules() - fGeom->GetNPPSDModules();
302
fa412d9b 303 relid[1] = 0 ; // means PBW04
304 gMC->CurrentVolOffID(4, relid[2]) ; // get the row number inside the module
305 gMC->CurrentVolOffID(3, relid[3]) ; // get the cell number inside the module
306
307 // get the absolute Id number
fa412d9b 308 fGeom->RelToAbsNumbering(relid, absid) ;
ed4205d8 309
fa412d9b 310 // add current hit to the hit list
b37750a6 311 AddHit(fIshunt, primary,tracknumber, absid, xyze);
037cc66d 312
ed4205d8 313
fa412d9b 314 } // there is deposited energy
315 } // we are inside a PHOS Xtal
316}
317
318//____________________________________________________________________________
319void AliPHOSv1::CPVDigitize (TLorentzVector p, Float_t *zxhit, Int_t moduleNumber, TClonesArray *cpvDigits)
320{
321 // ------------------------------------------------------------------------
322 // Digitize one CPV hit:
323 // On input take exact 4-momentum p and position zxhit of the hit,
324 // find the pad response around this hit and
325 // put the amplitudes in the pads into array digits
326 //
327 // Author: Yuri Kharlov (after Serguei Sadovsky)
328 // 2 October 2000
329 // ------------------------------------------------------------------------
330
a3dfe79c 331 const Float_t kCelWr = fGeom->GetPadSizePhi()/2; // Distance between wires (2 wires above 1 pad)
332 const Float_t kDetR = 0.1; // Relative energy fluctuation in track for 100 e-
333 const Float_t kdEdx = 4.0; // Average energy loss in CPV;
334 const Int_t kNgamz = 5; // Ionization size in Z
335 const Int_t kNgamx = 9; // Ionization size in Phi
336 const Float_t kNoise = 0.03; // charge noise in one pad
fa412d9b 337
338 Float_t rnor1,rnor2;
339
340 // Just a reminder on axes notation in the CPV module:
341 // axis Z goes along the beam
342 // axis X goes across the beam in the module plane
343 // axis Y is a normal to the module plane showing from the IP
344
345 Float_t hitX = zxhit[0];
346 Float_t hitZ =-zxhit[1];
347 Float_t pX = p.Px();
348 Float_t pZ =-p.Pz();
349 Float_t pNorm = p.Py();
a3dfe79c 350 Float_t eloss = kdEdx;
3d402178 351
fa412d9b 352 Float_t dZY = pZ/pNorm * fGeom->GetCPVGasThickness();
353 Float_t dXY = pX/pNorm * fGeom->GetCPVGasThickness();
354 gRandom->Rannor(rnor1,rnor2);
a3dfe79c 355 eloss *= (1 + kDetR*rnor1) *
356 TMath::Sqrt((1 + ( pow(dZY,2) + pow(dXY,2) ) / pow(fGeom->GetCPVGasThickness(),2)));
fa412d9b 357 Float_t zhit1 = hitZ + fGeom->GetCPVActiveSize(1)/2 - dZY/2;
358 Float_t xhit1 = hitX + fGeom->GetCPVActiveSize(0)/2 - dXY/2;
359 Float_t zhit2 = zhit1 + dZY;
360 Float_t xhit2 = xhit1 + dXY;
361
a3dfe79c 362 Int_t iwht1 = (Int_t) (xhit1 / kCelWr); // wire (x) coordinate "in"
363 Int_t iwht2 = (Int_t) (xhit2 / kCelWr); // wire (x) coordinate "out"
fa412d9b 364
365 Int_t nIter;
366 Float_t zxe[3][5];
367 if (iwht1==iwht2) { // incline 1-wire hit
368 nIter = 2;
369 zxe[0][0] = (zhit1 + zhit2 - dZY*0.57735) / 2;
a3dfe79c 370 zxe[1][0] = (iwht1 + 0.5) * kCelWr;
371 zxe[2][0] = eloss/2;
fa412d9b 372 zxe[0][1] = (zhit1 + zhit2 + dZY*0.57735) / 2;
a3dfe79c 373 zxe[1][1] = (iwht1 + 0.5) * kCelWr;
374 zxe[2][1] = eloss/2;
fa412d9b 375 }
376 else if (TMath::Abs(iwht1-iwht2) != 1) { // incline 3-wire hit
377 nIter = 3;
378 Int_t iwht3 = (iwht1 + iwht2) / 2;
a3dfe79c 379 Float_t xwht1 = (iwht1 + 0.5) * kCelWr; // wire 1
380 Float_t xwht2 = (iwht2 + 0.5) * kCelWr; // wire 2
381 Float_t xwht3 = (iwht3 + 0.5) * kCelWr; // wire 3
fa412d9b 382 Float_t xwr13 = (xwht1 + xwht3) / 2; // center 13
383 Float_t xwr23 = (xwht2 + xwht3) / 2; // center 23
384 Float_t dxw1 = xhit1 - xwr13;
385 Float_t dxw2 = xhit2 - xwr23;
a3dfe79c 386 Float_t egm1 = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
387 Float_t egm2 = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
388 Float_t egm3 = kCelWr / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
fa412d9b 389 zxe[0][0] = (dXY*(xwr13-xwht1)/dXY + zhit1 + zhit1) / 2;
390 zxe[1][0] = xwht1;
a3dfe79c 391 zxe[2][0] = eloss * egm1;
fa412d9b 392 zxe[0][1] = (dXY*(xwr23-xwht1)/dXY + zhit1 + zhit2) / 2;
393 zxe[1][1] = xwht2;
a3dfe79c 394 zxe[2][1] = eloss * egm2;
fa412d9b 395 zxe[0][2] = dXY*(xwht3-xwht1)/dXY + zhit1;
396 zxe[1][2] = xwht3;
a3dfe79c 397 zxe[2][2] = eloss * egm3;
fa412d9b 398 }
399 else { // incline 2-wire hit
400 nIter = 2;
a3dfe79c 401 Float_t xwht1 = (iwht1 + 0.5) * kCelWr;
402 Float_t xwht2 = (iwht2 + 0.5) * kCelWr;
fa412d9b 403 Float_t xwr12 = (xwht1 + xwht2) / 2;
404 Float_t dxw1 = xhit1 - xwr12;
405 Float_t dxw2 = xhit2 - xwr12;
406 Float_t egm1 = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) );
407 Float_t egm2 = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) );
408 zxe[0][0] = (zhit1 + zhit2 - dZY*egm1) / 2;
409 zxe[1][0] = xwht1;
a3dfe79c 410 zxe[2][0] = eloss * egm1;
fa412d9b 411 zxe[0][1] = (zhit1 + zhit2 + dZY*egm2) / 2;
412 zxe[1][1] = xwht2;
a3dfe79c 413 zxe[2][1] = eloss * egm2;
fa412d9b 414 }
bea63bea 415
fa412d9b 416 // Finite size of ionization region
417
ed4205d8 418 Int_t nCellZ = fGeom->GetNumberOfCPVPadsZ();
419 Int_t nCellX = fGeom->GetNumberOfCPVPadsPhi();
a3dfe79c 420 Int_t nz3 = (kNgamz+1)/2;
421 Int_t nx3 = (kNgamx+1)/2;
422 cpvDigits->Expand(nIter*kNgamx*kNgamz);
fa412d9b 423 TClonesArray &ldigits = *(TClonesArray *)cpvDigits;
424
425 for (Int_t iter=0; iter<nIter; iter++) {
426
427 Float_t zhit = zxe[0][iter];
428 Float_t xhit = zxe[1][iter];
429 Float_t qhit = zxe[2][iter];
430 Float_t zcell = zhit / fGeom->GetPadSizeZ();
431 Float_t xcell = xhit / fGeom->GetPadSizePhi();
432 if ( zcell<=0 || xcell<=0 ||
433 zcell>=nCellZ || xcell>=nCellX) return;
434 Int_t izcell = (Int_t) zcell;
435 Int_t ixcell = (Int_t) xcell;
436 Float_t zc = zcell - izcell - 0.5;
437 Float_t xc = xcell - ixcell - 0.5;
a3dfe79c 438 for (Int_t iz=1; iz<=kNgamz; iz++) {
fa412d9b 439 Int_t kzg = izcell + iz - nz3;
440 if (kzg<=0 || kzg>nCellZ) continue;
441 Float_t zg = (Float_t)(iz-nz3) - zc;
a3dfe79c 442 for (Int_t ix=1; ix<=kNgamx; ix++) {
fa412d9b 443 Int_t kxg = ixcell + ix - nx3;
444 if (kxg<=0 || kxg>nCellX) continue;
445 Float_t xg = (Float_t)(ix-nx3) - xc;
446
447 // Now calculate pad response
448 Float_t qpad = CPVPadResponseFunction(qhit,zg,xg);
a3dfe79c 449 qpad += kNoise*rnor2;
fa412d9b 450 if (qpad<0) continue;
451
452 // Fill the array with pad response ID and amplitude
3d402178 453 new(ldigits[cpvDigits->GetEntriesFast()]) AliPHOSCPVDigit(kxg,kzg,qpad);
fa412d9b 454 }
fa412d9b 455 }
fa412d9b 456 }
457}
458
459//____________________________________________________________________________
460Float_t AliPHOSv1::CPVPadResponseFunction(Float_t qhit, Float_t zhit, Float_t xhit) {
461 // ------------------------------------------------------------------------
462 // Calculate the amplitude in one CPV pad using the
463 // cumulative pad response function
464 // Author: Yuri Kharlov (after Serguei Sadovski)
465 // 3 October 2000
466 // ------------------------------------------------------------------------
467
468 Double_t dz = fGeom->GetPadSizeZ() / 2;
469 Double_t dx = fGeom->GetPadSizePhi() / 2;
470 Double_t z = zhit * fGeom->GetPadSizeZ();
471 Double_t x = xhit * fGeom->GetPadSizePhi();
472 Double_t amplitude = qhit *
473 (CPVCumulPadResponse(z+dz,x+dx) - CPVCumulPadResponse(z+dz,x-dx) -
474 CPVCumulPadResponse(z-dz,x+dx) + CPVCumulPadResponse(z-dz,x-dx));
475 return (Float_t)amplitude;
7587f5a5 476}
477
fa412d9b 478//____________________________________________________________________________
479Double_t AliPHOSv1::CPVCumulPadResponse(Double_t x, Double_t y) {
480 // ------------------------------------------------------------------------
481 // Cumulative pad response function
482 // It includes several terms from the CF decomposition in electrostatics
483 // Note: this cumulative function is wrong since omits some terms
484 // but the cell amplitude obtained with it is correct because
485 // these omitting terms cancel
486 // Author: Yuri Kharlov (after Serguei Sadovski)
487 // 3 October 2000
488 // ------------------------------------------------------------------------
489
a3dfe79c 490 const Double_t kA=1.0;
491 const Double_t kB=0.7;
fa412d9b 492
493 Double_t r2 = x*x + y*y;
494 Double_t xy = x*y;
495 Double_t cumulPRF = 0;
496 for (Int_t i=0; i<=4; i++) {
a3dfe79c 497 Double_t b1 = (2*i + 1) * kB;
fa412d9b 498 cumulPRF += TMath::Power(-1,i) * TMath::ATan( xy / (b1*TMath::Sqrt(b1*b1 + r2)) );
499 }
a3dfe79c 500 cumulPRF *= kA/(2*TMath::Pi());
fa412d9b 501 return cumulPRF;
502}
7eb9d12d 503