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