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