These are the interface classes to Geant 3 to read the geometry information
[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"
52#include "AliPHOSDigit.h"
bea63bea 53#include "AliPHOSReconstructioner.h"
7587f5a5 54#include "AliRun.h"
55#include "AliConst.h"
94de3818 56#include "AliMC.h"
7587f5a5 57
58ClassImp(AliPHOSv1)
59
bea63bea 60//____________________________________________________________________________
61AliPHOSv1::AliPHOSv1()
62{
5f20d3fb 63 // ctor
fc879520 64
3d402178 65 // Create an empty array of AliPHOSCPVModule to satisfy
fa412d9b 66 // AliPHOSv1::Streamer when reading root file
ed4205d8 67
68 fReconstructioner = 0;
69 fTrackSegmentMaker = 0;
70
71 if ( 0==(fEMCModules=new TClonesArray("AliPHOSCPVModule",0)) ) {
72 Error("AliPHOSv1","Can not create array of EMC modules");
73 exit(1);
74 }
fa412d9b 75
ed4205d8 76 if ( 0==(fCPVModules=new TClonesArray("AliPHOSCPVModule",0)) ) {
fa412d9b 77 Error("AliPHOSv1","Can not create array of CPV modules");
78 exit(1);
79 }
fc879520 80
bea63bea 81}
82
7587f5a5 83//____________________________________________________________________________
84AliPHOSv1::AliPHOSv1(const char *name, const char *title):
e04976bd 85AliPHOSv0(name,title)
7587f5a5 86{
5f20d3fb 87 // ctor : title is used to identify the layout
fa412d9b 88 // GPS2 = 5 modules (EMC + PPSD)
89 // IHEP = 5 modules (EMC + CPV )
ed4205d8 90 // MIXT = 4 modules (EMC + CPV ) and 1 module (EMC + PPSD)
5f20d3fb 91 //
ed4205d8 92 // We store hits :
5f20d3fb 93 // - fHits (the "normal" one), which retains the hits associated with
94 // the current primary particle being tracked
95 // (this array is reset after each primary has been tracked).
96 //
fa412d9b 97
5f20d3fb 98 fPinElectronicNoise = 0.010 ;
b95514c0 99 fDigitThreshold = 0.1 ; // 1 GeV
5f20d3fb 100
101 // We do not want to save in TreeH the raw hits
102 // But save the cumulated hits instead (need to create the branch myself)
103 // It is put in the Digit Tree because the TreeH is filled after each primary
fa412d9b 104 // and the TreeD at the end of the event (branch is set in FinishEvent() ).
5f20d3fb 105
ed4205d8 106 fHits= new TClonesArray("AliPHOSHit",1000) ;
5f20d3fb 107
ed4205d8 108 fNhits = 0 ;
5f20d3fb 109
ed4205d8 110 fReconstructioner = 0;
111 fTrackSegmentMaker = 0;
5f20d3fb 112
5f20d3fb 113 fIshunt = 1 ; // All hits are associated with primary particles
fa412d9b 114
ed4205d8 115 // Create array of EMC modules of the size of PHOS modules number
116
117 if ( 0==(fEMCModules=new TClonesArray("AliPHOSCPVModule",fGeom->GetNModules())) ) {
118 Error("AliPHOSv1","Can not create array of EMC modules");
119 exit(1);
120 }
121 TClonesArray &lemcmodule = *fEMCModules;
122 for (Int_t i=0; i<fGeom->GetNModules(); i++) new(lemcmodule[i]) AliPHOSCPVModule();
123
fa412d9b 124 // Create array of CPV modules for the IHEP's version of CPV
125
ed4205d8 126 if ( strcmp(fGeom->GetName(),"IHEP") == 0 || strcmp(fGeom->GetName(),"MIXT") == 0 ) {
127 // Create array of CPV modules of the size of PHOS modules number
fa412d9b 128
ed4205d8 129 if ( 0==(fCPVModules=new TClonesArray("AliPHOSCPVModule",fGeom->GetNCPVModules())) ) {
fa412d9b 130 Error("AliPHOSv1","Can not create array of CPV modules");
131 exit(1);
132 }
133 TClonesArray &lcpvmodule = *fCPVModules;
ed4205d8 134 for (Int_t i=0; i<fGeom->GetNCPVModules(); i++) new(lcpvmodule[i]) AliPHOSCPVModule();
fa412d9b 135 }
136 else {
3d402178 137 // Create an empty array of AliPHOSCPVModule to satisfy
fa412d9b 138 // AliPHOSv1::Streamer when writing root file
139
3d402178 140 fCPVModules=new TClonesArray("AliPHOSCPVModule",0);
fa412d9b 141
142 }
5f20d3fb 143}
144
145//____________________________________________________________________________
146AliPHOSv1::AliPHOSv1(AliPHOSReconstructioner * Reconstructioner, const char *name, const char *title):
147 AliPHOSv0(name,title)
148{
149 // ctor : title is used to identify the layout
150 // GPS2 = 5 modules (EMC + PPSD)
151 // We use 2 arrays of hits :
152 //
153 // - fHits (the "normal" one), which retains the hits associated with
154 // the current primary particle being tracked
155 // (this array is reset after each primary has been tracked).
156 //
157 // - fTmpHits, which retains all the hits of the current event. It
158 // is used for the digitization part.
159
160 fPinElectronicNoise = 0.010 ;
161
162 // We do not want to save in TreeH the raw hits
163 //fHits = new TClonesArray("AliPHOSHit",100) ;
164
ed4205d8 165 fDigits = 0 ;
166 fHits= new TClonesArray("AliPHOSHit",1000) ;
b73f246d 167
ed4205d8 168 fNhits = 0 ;
5f20d3fb 169
170 fIshunt = 1 ; // All hits are associated with primary particles
171
172 // gets an instance of the geometry parameters class
173 fGeom = AliPHOSGeometry::GetInstance(title, "") ;
174
175 if (fGeom->IsInitialized() )
88bdfa12 176 cout << "AliPHOS" << Version() << " : PHOS geometry intialized for " << fGeom->GetName() << endl ;
5f20d3fb 177 else
fc879520 178 cout << "AliPHOS" << Version() << " : PHOS geometry initialization failed !" << endl ;
5f20d3fb 179
180 // Defining the PHOS Reconstructioner
181
182 fReconstructioner = Reconstructioner ;
183
7587f5a5 184}
b2a60966 185
7587f5a5 186//____________________________________________________________________________
bea63bea 187AliPHOSv1::~AliPHOSv1()
b2a60966 188{
bea63bea 189 // dtor
5f20d3fb 190
ed4205d8 191 if ( fHits) {
192 fHits->Delete() ;
193 delete fHits ;
194 fHits = 0 ;
8dfa469d 195 }
5f20d3fb 196
ed4205d8 197 if ( fDigits) {
198 fDigits->Delete() ;
199 delete fDigits ;
200 fDigits = 0 ;
201 }
5f20d3fb 202
ed4205d8 203 if ( fEmcRecPoints ) {
204 fEmcRecPoints->Delete() ;
205 delete fEmcRecPoints ;
206 fEmcRecPoints = 0 ;
207 }
208
209 if ( fPpsdRecPoints ) {
210 fPpsdRecPoints->Delete() ;
211 delete fPpsdRecPoints ;
212 fPpsdRecPoints = 0 ;
213 }
fad3e5b9 214
ed4205d8 215 if ( fTrackSegments ) {
216 fTrackSegments->Delete() ;
217 delete fTrackSegments ;
218 fTrackSegments = 0 ;
219 }
220
7587f5a5 221}
222
7587f5a5 223//____________________________________________________________________________
7eb9d12d 224void AliPHOSv1::AddHit(Int_t shunt, Int_t primary, Int_t tracknumber, Int_t Id, Float_t * hits, Int_t trackpid)
bea63bea 225{
226 // Add a hit to the hit list.
5f20d3fb 227 // A PHOS hit is the sum of all hits in a single crystal
228 // or in a single PPSD gas cell
bea63bea 229
5f20d3fb 230 Int_t hitCounter ;
bea63bea 231 AliPHOSHit *newHit ;
5f20d3fb 232 AliPHOSHit *curHit ;
233 Bool_t deja = kFALSE ;
bea63bea 234
7eb9d12d 235 newHit = new AliPHOSHit(shunt, primary, tracknumber, Id, hits, trackpid) ;
bea63bea 236
ed4205d8 237 for ( hitCounter = 0 ; hitCounter < fNhits && !deja ; hitCounter++ ) {
238 curHit = (AliPHOSHit*) (*fHits)[hitCounter] ;
239 if( *curHit == *newHit ) {
240 *curHit = *curHit + *newHit ;
241 deja = kTRUE ;
5f20d3fb 242 }
243 }
244
245 if ( !deja ) {
ed4205d8 246 new((*fHits)[fNhits]) AliPHOSHit(*newHit) ;
247 fNhits++ ;
5f20d3fb 248 }
249
bea63bea 250 delete newHit;
bea63bea 251}
252
5f20d3fb 253//___________________________________________________________________________
254Int_t AliPHOSv1::Digitize(Float_t Energy)
255{
256 // Applies the energy calibration
257
c3c187e5 258 Float_t fB = 10000000. ;
5f20d3fb 259 Float_t fA = 0. ;
260 Int_t chan = Int_t(fA + Energy*fB ) ;
261 return chan ;
262}
bea63bea 263
ed4205d8 264//____________________________________________________________________________
265void AliPHOSv1::Hit2Digit(Int_t ntracks){
266 //Collects all hits in the same active volume into digits
267
268 if(fDigits!= 0)
269 fDigits->Clear() ;
270 else
271 fDigits = new TClonesArray("AliPHOSDigit",1000) ;
272
273 // Branch address for digit tree
274 char branchname[20];
275 sprintf(branchname,"%s",GetName());
276 gAlice->TreeD()->Branch(branchname,&fDigits,fBufferSize);
bea63bea 277
ed4205d8 278 gAlice->TreeD()->GetEvent(0);
279
5f20d3fb 280
bea63bea 281 Int_t i ;
282 Int_t relid[4];
283 Int_t j ;
fa412d9b 284 AliPHOSHit * hit ;
bea63bea 285 AliPHOSDigit * newdigit ;
286 AliPHOSDigit * curdigit ;
287 Bool_t deja = kFALSE ;
288
ed4205d8 289 Int_t itrack ;
290 for (itrack=0; itrack<ntracks; itrack++){
291
292 //=========== Get the Hits Tree for the Primary track itrack
293 gAlice->ResetHits();
294 gAlice->TreeH()->GetEvent(itrack);
295
296 for ( i = 0 ; i < fHits->GetEntries() ; i++ ) {
297 hit = (AliPHOSHit*)fHits->At(i) ;
298
299 // Assign primary number only if contribution is significant
300 if( hit->GetEnergy() > fDigitThreshold)
301 newdigit = new AliPHOSDigit( hit->GetPrimary(), hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
302 else
303 newdigit = new AliPHOSDigit( -1 , hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
304 deja =kFALSE ;
305
306
307 for ( j = 0 ; j < fNdigits ; j++) {
308 curdigit = (AliPHOSDigit*) fDigits->At(j) ;
309 if ( *curdigit == *newdigit) {
310 *curdigit = *curdigit + *newdigit ;
311 deja = kTRUE ;
312 }
bea63bea 313 }
bea63bea 314
ed4205d8 315 if ( !deja ) {
316 new((*fDigits)[fNdigits]) AliPHOSDigit(* newdigit) ;
317 fNdigits++ ;
318 }
319
320 delete newdigit ;
321 }
322
323 } // loop over tracks
324
325 // Noise induced by the PIN diode of the PbWO crystals
326
bea63bea 327 Float_t energyandnoise ;
328 for ( i = 0 ; i < fNdigits ; i++ ) {
329 newdigit = (AliPHOSDigit * ) fDigits->At(i) ;
ed4205d8 330
bea63bea 331 fGeom->AbsToRelNumbering(newdigit->GetId(), relid) ;
ed4205d8 332
bea63bea 333 if (relid[1]==0){ // Digits belong to EMC (PbW0_4 crystals)
334 energyandnoise = newdigit->GetAmp() + Digitize(gRandom->Gaus(0., fPinElectronicNoise)) ;
ed4205d8 335
bea63bea 336 if (energyandnoise < 0 )
337 energyandnoise = 0 ;
ed4205d8 338
bea63bea 339 if ( newdigit->GetAmp() < fDigitThreshold ) // if threshold not surpassed, remove digit from list
340 fDigits->RemoveAt(i) ;
341 }
342 }
ed4205d8 343
bea63bea 344 fDigits->Compress() ;
ed4205d8 345
346 fNdigits = fDigits->GetEntries() ;
347 fDigits->Expand(fNdigits) ;
bea63bea 348
bea63bea 349 for (i = 0 ; i < fNdigits ; i++) {
350 newdigit = (AliPHOSDigit *) fDigits->At(i) ;
ed4205d8 351 newdigit->SetIndexInList(i) ;
bea63bea 352 }
fa412d9b 353
ed4205d8 354 gAlice->TreeD()->Fill() ;
5f20d3fb 355
ed4205d8 356 gAlice->TreeD()->Write(0,TObject::kOverwrite) ;
357
358}
2ab0c725 359
5f20d3fb 360//___________________________________________________________________________
2ab0c725 361void AliPHOSv1::MakeBranch(Option_t* opt, char *file)
5f20d3fb 362{
363 // Create new branche in the current Root Tree in the digit Tree
5f20d3fb 364 AliDetector::MakeBranch(opt) ;
365
ed4205d8 366 // Create new branches EMC<i> for hits in EMC modules
367
2ab0c725 368 for( Int_t i=0; i<fGeom->GetNModules(); i++ ) GetEMCModule(i).MakeBranch("EMC",i+1,file);
ed4205d8 369
fa412d9b 370 // Create new branches CPV<i> for hits in CPV modules for IHEP geometry
ed4205d8 371
372 if ( strcmp(fGeom->GetName(),"IHEP") == 0 || strcmp(fGeom->GetName(),"MIXT") == 0 ) {
2ab0c725 373 for( Int_t i=0; i<fGeom->GetNCPVModules(); i++ ) GetCPVModule(i).MakeBranch("CPV",i+1,file);
fa412d9b 374 }
ed4205d8 375
bea63bea 376}
377
5f20d3fb 378//_____________________________________________________________________________
379void AliPHOSv1::Reconstruction(AliPHOSReconstructioner * Reconstructioner)
380{
381 // 1. Reinitializes the existing RecPoint, TrackSegment, and RecParticles Lists and
382 // 2. Creates TreeR with a branch for each list
383 // 3. Steers the reconstruction processes
384 // 4. Saves the 3 lists in TreeR
385 // 5. Write the Tree to File
386
387 fReconstructioner = Reconstructioner ;
388
389 char branchname[10] ;
390
391 // 1.
392
393 // gAlice->MakeTree("R") ;
394 Int_t splitlevel = 0 ;
395
fc879520 396 fEmcRecPoints->Delete() ;
5f20d3fb 397
398 if ( fEmcRecPoints && gAlice->TreeR() ) {
399 sprintf(branchname,"%sEmcRP",GetName()) ;
5f20d3fb 400 gAlice->TreeR()->Branch(branchname, "TObjArray", &fEmcRecPoints, fBufferSize, splitlevel) ;
401 }
402
fa412d9b 403 fPpsdRecPoints->Delete() ;
5f20d3fb 404
5f20d3fb 405 if ( fPpsdRecPoints && gAlice->TreeR() ) {
406 sprintf(branchname,"%sPpsdRP",GetName()) ;
5f20d3fb 407 gAlice->TreeR()->Branch(branchname, "TObjArray", &fPpsdRecPoints, fBufferSize, splitlevel) ;
408 }
409
fc879520 410 fTrackSegments->Delete() ;
5f20d3fb 411
5f20d3fb 412 if ( fTrackSegments && gAlice->TreeR() ) {
413 sprintf(branchname,"%sTS",GetName()) ;
414 gAlice->TreeR()->Branch(branchname, &fTrackSegments, fBufferSize) ;
415 }
416
fa412d9b 417 fRecParticles->Delete() ;
fc879520 418
fad3e5b9 419 if ( fRecParticles && gAlice->TreeR() ) {
420 sprintf(branchname,"%sRP",GetName()) ;
421 gAlice->TreeR()->Branch(branchname, &fRecParticles, fBufferSize) ;
5f20d3fb 422 }
423
424 // 3.
b73f246d 425
fad3e5b9 426 fReconstructioner->Make(fDigits, fEmcRecPoints, fPpsdRecPoints, fTrackSegments, fRecParticles);
b73f246d 427
fad3e5b9 428 printf("Reconstruction: %d %d %d %d\n",
429 fEmcRecPoints->GetEntries(),fPpsdRecPoints->GetEntries(),
b73f246d 430 fTrackSegments->GetEntries(),fRecParticles->GetEntries());
5f20d3fb 431
432 // 4. Expand or Shrink the arrays to the proper size
433
434 Int_t size ;
435
436 size = fEmcRecPoints->GetEntries() ;
437 fEmcRecPoints->Expand(size) ;
ed4205d8 438
5f20d3fb 439 size = fPpsdRecPoints->GetEntries() ;
440 fPpsdRecPoints->Expand(size) ;
441
442 size = fTrackSegments->GetEntries() ;
443 fTrackSegments->Expand(size) ;
444
445 size = fRecParticles->GetEntries() ;
446 fRecParticles->Expand(size) ;
447
448 gAlice->TreeR()->Fill() ;
5f20d3fb 449 // 5.
450
4a2ca5e9 451 gAlice->TreeR()->Write(0,TObject::kOverwrite) ;
5f20d3fb 452
453 // Deleting reconstructed objects
454 ResetReconstruction();
455
456}
457
fa412d9b 458//____________________________________________________________________________
459void AliPHOSv1::ResetHits()
c4f224e7 460{
ed4205d8 461 // Reset hit tree for EMC and CPV
fa412d9b 462 // Yuri Kharlov, 28 September 2000
463
464 AliDetector::ResetHits();
ed4205d8 465 for (Int_t i=0; i<fGeom->GetNModules(); i++) ((AliPHOSCPVModule*)(*fEMCModules)[i]) -> Clear();
466 if ( strcmp(fGeom->GetName(),"IHEP") == 0 || strcmp(fGeom->GetName(),"MIXT") == 0 ) {
467 for (Int_t i=0; i<fGeom->GetNCPVModules(); i++) ((AliPHOSCPVModule*)(*fCPVModules)[i]) -> Clear();
5f20d3fb 468 }
ed4205d8 469
5f20d3fb 470}
471//____________________________________________________________________________
472void AliPHOSv1::ResetReconstruction()
473{
474 // Deleting reconstructed objects
475
b73f246d 476 if ( fEmcRecPoints ) fEmcRecPoints ->Delete();
5f20d3fb 477 if ( fPpsdRecPoints ) fPpsdRecPoints->Delete();
478 if ( fTrackSegments ) fTrackSegments->Delete();
b73f246d 479 if ( fRecParticles ) fRecParticles ->Delete();
5f20d3fb 480
481}
5f20d3fb 482
483//____________________________________________________________________________
484void AliPHOSv1::SetTreeAddress()
485{
486 // TBranch *branch;
487 AliPHOS::SetTreeAddress();
488
ed4205d8 489// //Branch address for TreeR: RecPpsdRecPoint
5f20d3fb 490// TTree *treeR = gAlice->TreeR();
491// if ( treeR && fPpsdRecPoints ) {
492// branch = treeR->GetBranch("PHOSPpsdRP");
493// if (branch) branch->SetAddress(&fPpsdRecPoints) ;
fa412d9b 494// }
495
ed4205d8 496 // Set branch address for the Hits Tree for hits in EMC modules
497 // Yuri Kharlov, 23 November 2000.
498
499 for( Int_t i=0; i<fGeom->GetNModules(); i++ ) GetEMCModule(i).SetTreeAddress("EMC",i+1);
500
fa412d9b 501 // Set branch address for the Hits Tree for hits in CPV modules for IHEP geometry
502 // Yuri Kharlov, 28 September 2000.
503
ed4205d8 504 if ( strcmp(fGeom->GetName(),"IHEP") == 0 || strcmp(fGeom->GetName(),"MIXT") == 0 ) {
505 for( Int_t i=0; i<fGeom->GetNCPVModules(); i++ ) GetCPVModule(i).SetTreeAddress("CPV",i+1);
fa412d9b 506 }
507
5f20d3fb 508}
509
510//____________________________________________________________________________
511
7587f5a5 512void AliPHOSv1::StepManager(void)
513{
b2a60966 514 // Accumulates hits as long as the track stays in a single crystal or PPSD gas Cell
b2a60966 515
ed4205d8 516// if (gMC->IsTrackEntering())
517// cout << "Track enters the volume " << gMC->CurrentVolName() << endl;
518// if (gMC->IsTrackExiting())
519// cout << "Track leaves the volume " << gMC->CurrentVolName() << endl;
520
7587f5a5 521 Int_t relid[4] ; // (box, layer, row, column) indices
fa412d9b 522 Int_t absid ; // absolute cell ID number
2168f43b 523 Float_t xyze[4]={0,0,0,0} ; // position wrt MRS and energy deposited
fa412d9b 524 TLorentzVector pos ; // Lorentz vector of the track current position
525 Int_t copy ;
ed4205d8 526 Int_t i ;
7587f5a5 527
bea63bea 528 Int_t tracknumber = gAlice->CurrentTrack() ;
fa412d9b 529 Int_t primary = gAlice->GetPrimary( gAlice->CurrentTrack() );
530 TString name = fGeom->GetName() ;
7eb9d12d 531 Int_t trackpid = gMC->TrackPid() ;
ed4205d8 532 if ( name == "GPS2" || name == "MIXT" ) { // ======> CPV is a GPS' PPSD
fa412d9b 533
b2a60966 534 if( gMC->CurrentVolID(copy) == gMC->VolId("GCEL") ) // We are inside a gas cell
7587f5a5 535 {
536 gMC->TrackPosition(pos) ;
537 xyze[0] = pos[0] ;
538 xyze[1] = pos[1] ;
539 xyze[2] = pos[2] ;
bea63bea 540 xyze[3] = gMC->Edep() ;
7587f5a5 541
542 if ( xyze[3] != 0 ) { // there is deposited energy
543 gMC->CurrentVolOffID(5, relid[0]) ; // get the PHOS Module number
fad3e5b9 544 if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(5),"PHO1") == 0 ){
ed4205d8 545 relid[0] += fGeom->GetNModules() - fGeom->GetNPPSDModules();
fad3e5b9 546 }
7587f5a5 547 gMC->CurrentVolOffID(3, relid[1]) ; // get the Micromegas Module number
ed4205d8 548 // 1-> fGeom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() upper
549 // > fGeom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() lower
7587f5a5 550 gMC->CurrentVolOffID(1, relid[2]) ; // get the row number of the cell
551 gMC->CurrentVolID(relid[3]) ; // get the column number
552
553 // get the absolute Id number
554
bea63bea 555 fGeom->RelToAbsNumbering(relid, absid) ;
7587f5a5 556
bea63bea 557 // add current hit to the hit list
7eb9d12d 558 AddHit(fIshunt, primary, tracknumber, absid, xyze, trackpid);
7587f5a5 559
560 } // there is deposited energy
fa412d9b 561 } // We are inside the gas of the CPV
562 } // GPS2 configuration
563
ed4205d8 564 if ( name == "IHEP" || name == "MIXT" ) { // ======> CPV is a IHEP's one
fa412d9b 565
566 // Yuri Kharlov, 28 September 2000
567
568 if( gMC->CurrentVolID(copy) == gMC->VolId("CPVQ") &&
569 gMC->IsTrackEntering() &&
570 gMC->TrackCharge() != 0) {
571
572 // Charged track has just entered to the CPV sensitive plane
573
a3dfe79c 574 AliPHOSv1 &phos = *(AliPHOSv1*)gAlice->GetModule("PHOS");
fa412d9b 575
a3dfe79c 576 Int_t moduleNumber;
577 gMC->CurrentVolOffID(3,moduleNumber);
578 moduleNumber--;
fa412d9b 579
580 // Current position of the hit in the CPV module ref. system
581
582 gMC -> TrackPosition(pos);
2168f43b 583 Float_t xyzm[3], xyzd[3], xyd[3]={0,0,0};
cd461ab8 584 Int_t i;
585 for (i=0; i<3; i++) xyzm[i] = pos[i];
fa412d9b 586 gMC -> Gmtod (xyzm, xyzd, 1); // transform coordinate from master to daughter system
587 xyd[0] = xyzd[0];
588 xyd[1] =-xyzd[2];
589
590 // Current momentum of the hit's track in the CPV module ref. system
591
592 TLorentzVector pmom;
593 gMC -> TrackMomentum(pmom);
594 Float_t pm[3], pd[3];
cd461ab8 595 for (i=0; i<3; i++) pm[i] = pmom[i];
fa412d9b 596 gMC -> Gmtod (pm, pd, 2); // transform 3-momentum from master to daughter system
597 pmom[0] = pd[0];
598 pmom[1] =-pd[1];
599 pmom[2] =-pd[2];
600
601 // Current particle type of the hit's track
602
603 Int_t ipart = gMC->TrackPid();
604
605 // Add the current particle in the list of the CPV hits.
606
ed4205d8 607 phos.GetCPVModule(moduleNumber).AddHit(fIshunt,primary,pmom,xyd,ipart);
fa412d9b 608
ed4205d8 609 if (fDebug == 1) {
a3dfe79c 610 printf("CPV hit added to module #%2d: p = (% .4f, % .4f, % .4f, % .4f) GeV,\n",
611 moduleNumber+1,pmom.Px(),pmom.Py(),pmom.Pz(),pmom.E());
ed4205d8 612 printf( " xy = (%8.4f, %8.4f) cm, ipart = %d, primary = %d\n",
613 xyd[0],xyd[1],ipart,primary);
a3dfe79c 614 }
fa412d9b 615
616 // Digitize the current CPV hit:
617
618 // 1. find pad response and
619
3d402178 620 TClonesArray *cpvDigits = new TClonesArray("AliPHOSCPVDigit",0); // array of digits for current hit
a3dfe79c 621 CPVDigitize(pmom,xyd,moduleNumber,cpvDigits);
fa412d9b 622
623 Float_t xmean = 0;
624 Float_t zmean = 0;
625 Float_t qsum = 0;
cd461ab8 626 Int_t idigit,ndigits;
fa412d9b 627
628 // 2. go through the current digit list and sum digits in pads
629
630 ndigits = cpvDigits->GetEntriesFast();
cd461ab8 631 for (idigit=0; idigit<ndigits-1; idigit++) {
3d402178 632 AliPHOSCPVDigit *cpvDigit1 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit);
fa412d9b 633 Float_t x1 = cpvDigit1->GetXpad() ;
634 Float_t z1 = cpvDigit1->GetYpad() ;
635 for (Int_t jdigit=idigit+1; jdigit<ndigits; jdigit++) {
3d402178 636 AliPHOSCPVDigit *cpvDigit2 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(jdigit);
fa412d9b 637 Float_t x2 = cpvDigit2->GetXpad() ;
638 Float_t z2 = cpvDigit2->GetYpad() ;
639 if (x1==x2 && z1==z2) {
640 Float_t qsum = cpvDigit1->GetQpad() + cpvDigit2->GetQpad() ;
641 cpvDigit2->SetQpad(qsum) ;
642 cpvDigits->RemoveAt(idigit) ;
643 }
644 }
645 }
646 cpvDigits->Compress() ;
647
648 // 3. add digits to temporary hit list fTmpHits
649
650 ndigits = cpvDigits->GetEntriesFast();
cd461ab8 651 for (idigit=0; idigit<ndigits; idigit++) {
3d402178 652 AliPHOSCPVDigit *cpvDigit = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit);
a3dfe79c 653 relid[0] = moduleNumber + 1 ; // CPV (or PHOS) module number
fa412d9b 654 relid[1] =-1 ; // means CPV
655 relid[2] = cpvDigit->GetXpad() ; // column number of a pad
656 relid[3] = cpvDigit->GetYpad() ; // row number of a pad
657
658 // get the absolute Id number
659 fGeom->RelToAbsNumbering(relid, absid) ;
660
661 // add current digit to the temporary hit list
662 xyze[0] = 0. ;
663 xyze[1] = 0. ;
664 xyze[2] = 0. ;
665 xyze[3] = cpvDigit->GetQpad() ; // amplitude in a pad
666 primary = -1; // No need in primary for CPV
7eb9d12d 667 AddHit(fIshunt, primary, tracknumber, absid, xyze, trackpid);
fa412d9b 668
669 if (cpvDigit->GetQpad() > 0.02) {
670 xmean += cpvDigit->GetQpad() * (cpvDigit->GetXpad() + 0.5);
671 zmean += cpvDigit->GetQpad() * (cpvDigit->GetYpad() + 0.5);
672 qsum += cpvDigit->GetQpad();
673 }
674 }
675 delete cpvDigits;
676 }
677 } // end of IHEP configuration
7587f5a5 678
fa412d9b 679 if(gMC->CurrentVolID(copy) == gMC->VolId("PXTL") ) { // We are inside a PBWO crystal
680 gMC->TrackPosition(pos) ;
681 xyze[0] = pos[0] ;
682 xyze[1] = pos[1] ;
683 xyze[2] = pos[2] ;
684 xyze[3] = gMC->Edep() ;
ed4205d8 685
686 // Track enters to the crystal from the top edge
687
688 if (gMC->IsTrackEntering()) {
689 Float_t posloc[3];
690 gMC -> Gmtod (xyze, posloc, 1);
691 if (posloc[1] > fGeom->GetCrystalSize(1)/2-0.01) {
692 Int_t row,cel;
76eb93f7 693 Float_t xyd[3]={0,0,0};
ed4205d8 694 AliPHOSv1 &phos = *(AliPHOSv1*)gAlice->GetModule("PHOS");
695
696 Int_t moduleNumber;
697 gMC->CurrentVolOffID(10,moduleNumber);
fad3e5b9 698
ed4205d8 699 if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(10),"PHO1") == 0 )
700 moduleNumber += fGeom->GetNModules() - fGeom->GetNPPSDModules();
701 moduleNumber--;
702
703 gMC->CurrentVolOffID(4, row) ; // get the row number inside the module
704 gMC->CurrentVolOffID(3, cel) ; // get the cell number inside the module
705 xyd[0] = -(posloc[2] + (cel-0.5-fGeom->GetNZ() /2) *
706 (fGeom->GetCrystalSize(2) + 2 * fGeom->GetGapBetweenCrystals()));
707 xyd[1] = posloc[0] + (row-0.5-fGeom->GetNPhi()/2) *
708 (fGeom->GetCrystalSize(0) + 2 * fGeom->GetGapBetweenCrystals());
709
710 // Current momentum of the hit's track in the CPV module ref. system
711
712 TLorentzVector pmom;
713 gMC -> TrackMomentum(pmom);
714 Float_t pm[3], pd[3];
715 for (i=0; i<3; i++) pm[i] = pmom[i];
716 gMC -> Gmtod (pm, pd, 2); // transform 3-momentum from master to daughter system
717 pmom[0] = pd[0];
718 pmom[1] =-pd[1];
719 pmom[2] =-pd[2];
720
721 // Current particle type of the hit's track
722
723 Int_t ipart = gMC->TrackPid();
724
725 // Add the current particle in the list of the EMC hits.
726
727 phos.GetEMCModule(moduleNumber).AddHit(fIshunt,primary,pmom,xyd,ipart);
728
729 if (fDebug == 1) {
730 printf("EMC hit added to module #%2d: p = (% .4f, % .4f, % .4f, % .4f) GeV,\n",
731 moduleNumber+1,pmom.Px(),pmom.Py(),pmom.Pz(),pmom.E());
732 printf( " xy = (%8.4f, %8.4f) cm, ipart = %d, primary = %d\n",
733 xyd[0],xyd[1],ipart,primary);
734 }
735 }
736 }
737
738 // Track is inside the crystal and deposits some energy
739
fa412d9b 740 if ( xyze[3] != 0 ) {
741 gMC->CurrentVolOffID(10, relid[0]) ; // get the PHOS module number ;
fad3e5b9 742 if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(10),"PHO1") == 0 )
743 relid[0] += fGeom->GetNModules() - fGeom->GetNPPSDModules();
744
fa412d9b 745 relid[1] = 0 ; // means PBW04
746 gMC->CurrentVolOffID(4, relid[2]) ; // get the row number inside the module
747 gMC->CurrentVolOffID(3, relid[3]) ; // get the cell number inside the module
748
749 // get the absolute Id number
fa412d9b 750 fGeom->RelToAbsNumbering(relid, absid) ;
ed4205d8 751
fa412d9b 752 // add current hit to the hit list
7eb9d12d 753 AddHit(fIshunt, primary,tracknumber, absid, xyze, trackpid);
ed4205d8 754
fa412d9b 755 } // there is deposited energy
756 } // we are inside a PHOS Xtal
757}
758
759//____________________________________________________________________________
760void AliPHOSv1::CPVDigitize (TLorentzVector p, Float_t *zxhit, Int_t moduleNumber, TClonesArray *cpvDigits)
761{
762 // ------------------------------------------------------------------------
763 // Digitize one CPV hit:
764 // On input take exact 4-momentum p and position zxhit of the hit,
765 // find the pad response around this hit and
766 // put the amplitudes in the pads into array digits
767 //
768 // Author: Yuri Kharlov (after Serguei Sadovsky)
769 // 2 October 2000
770 // ------------------------------------------------------------------------
771
a3dfe79c 772 const Float_t kCelWr = fGeom->GetPadSizePhi()/2; // Distance between wires (2 wires above 1 pad)
773 const Float_t kDetR = 0.1; // Relative energy fluctuation in track for 100 e-
774 const Float_t kdEdx = 4.0; // Average energy loss in CPV;
775 const Int_t kNgamz = 5; // Ionization size in Z
776 const Int_t kNgamx = 9; // Ionization size in Phi
777 const Float_t kNoise = 0.03; // charge noise in one pad
fa412d9b 778
779 Float_t rnor1,rnor2;
780
781 // Just a reminder on axes notation in the CPV module:
782 // axis Z goes along the beam
783 // axis X goes across the beam in the module plane
784 // axis Y is a normal to the module plane showing from the IP
785
786 Float_t hitX = zxhit[0];
787 Float_t hitZ =-zxhit[1];
788 Float_t pX = p.Px();
789 Float_t pZ =-p.Pz();
790 Float_t pNorm = p.Py();
a3dfe79c 791 Float_t eloss = kdEdx;
3d402178 792
7eb9d12d 793// cout << "CPVDigitize: YVK : "<<hitX<<" "<<hitZ<<" | "<<pX<<" "<<pZ<<" "<<pNorm<<endl;
794
fa412d9b 795 Float_t dZY = pZ/pNorm * fGeom->GetCPVGasThickness();
796 Float_t dXY = pX/pNorm * fGeom->GetCPVGasThickness();
797 gRandom->Rannor(rnor1,rnor2);
a3dfe79c 798 eloss *= (1 + kDetR*rnor1) *
799 TMath::Sqrt((1 + ( pow(dZY,2) + pow(dXY,2) ) / pow(fGeom->GetCPVGasThickness(),2)));
fa412d9b 800 Float_t zhit1 = hitZ + fGeom->GetCPVActiveSize(1)/2 - dZY/2;
801 Float_t xhit1 = hitX + fGeom->GetCPVActiveSize(0)/2 - dXY/2;
802 Float_t zhit2 = zhit1 + dZY;
803 Float_t xhit2 = xhit1 + dXY;
804
a3dfe79c 805 Int_t iwht1 = (Int_t) (xhit1 / kCelWr); // wire (x) coordinate "in"
806 Int_t iwht2 = (Int_t) (xhit2 / kCelWr); // wire (x) coordinate "out"
fa412d9b 807
808 Int_t nIter;
809 Float_t zxe[3][5];
810 if (iwht1==iwht2) { // incline 1-wire hit
811 nIter = 2;
812 zxe[0][0] = (zhit1 + zhit2 - dZY*0.57735) / 2;
a3dfe79c 813 zxe[1][0] = (iwht1 + 0.5) * kCelWr;
814 zxe[2][0] = eloss/2;
fa412d9b 815 zxe[0][1] = (zhit1 + zhit2 + dZY*0.57735) / 2;
a3dfe79c 816 zxe[1][1] = (iwht1 + 0.5) * kCelWr;
817 zxe[2][1] = eloss/2;
fa412d9b 818 }
819 else if (TMath::Abs(iwht1-iwht2) != 1) { // incline 3-wire hit
820 nIter = 3;
821 Int_t iwht3 = (iwht1 + iwht2) / 2;
a3dfe79c 822 Float_t xwht1 = (iwht1 + 0.5) * kCelWr; // wire 1
823 Float_t xwht2 = (iwht2 + 0.5) * kCelWr; // wire 2
824 Float_t xwht3 = (iwht3 + 0.5) * kCelWr; // wire 3
fa412d9b 825 Float_t xwr13 = (xwht1 + xwht3) / 2; // center 13
826 Float_t xwr23 = (xwht2 + xwht3) / 2; // center 23
827 Float_t dxw1 = xhit1 - xwr13;
828 Float_t dxw2 = xhit2 - xwr23;
a3dfe79c 829 Float_t egm1 = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
830 Float_t egm2 = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
831 Float_t egm3 = kCelWr / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
fa412d9b 832 zxe[0][0] = (dXY*(xwr13-xwht1)/dXY + zhit1 + zhit1) / 2;
833 zxe[1][0] = xwht1;
a3dfe79c 834 zxe[2][0] = eloss * egm1;
fa412d9b 835 zxe[0][1] = (dXY*(xwr23-xwht1)/dXY + zhit1 + zhit2) / 2;
836 zxe[1][1] = xwht2;
a3dfe79c 837 zxe[2][1] = eloss * egm2;
fa412d9b 838 zxe[0][2] = dXY*(xwht3-xwht1)/dXY + zhit1;
839 zxe[1][2] = xwht3;
a3dfe79c 840 zxe[2][2] = eloss * egm3;
fa412d9b 841 }
842 else { // incline 2-wire hit
843 nIter = 2;
a3dfe79c 844 Float_t xwht1 = (iwht1 + 0.5) * kCelWr;
845 Float_t xwht2 = (iwht2 + 0.5) * kCelWr;
fa412d9b 846 Float_t xwr12 = (xwht1 + xwht2) / 2;
847 Float_t dxw1 = xhit1 - xwr12;
848 Float_t dxw2 = xhit2 - xwr12;
849 Float_t egm1 = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) );
850 Float_t egm2 = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) );
851 zxe[0][0] = (zhit1 + zhit2 - dZY*egm1) / 2;
852 zxe[1][0] = xwht1;
a3dfe79c 853 zxe[2][0] = eloss * egm1;
fa412d9b 854 zxe[0][1] = (zhit1 + zhit2 + dZY*egm2) / 2;
855 zxe[1][1] = xwht2;
a3dfe79c 856 zxe[2][1] = eloss * egm2;
fa412d9b 857 }
bea63bea 858
fa412d9b 859 // Finite size of ionization region
860
ed4205d8 861 Int_t nCellZ = fGeom->GetNumberOfCPVPadsZ();
862 Int_t nCellX = fGeom->GetNumberOfCPVPadsPhi();
a3dfe79c 863 Int_t nz3 = (kNgamz+1)/2;
864 Int_t nx3 = (kNgamx+1)/2;
865 cpvDigits->Expand(nIter*kNgamx*kNgamz);
fa412d9b 866 TClonesArray &ldigits = *(TClonesArray *)cpvDigits;
867
868 for (Int_t iter=0; iter<nIter; iter++) {
869
870 Float_t zhit = zxe[0][iter];
871 Float_t xhit = zxe[1][iter];
872 Float_t qhit = zxe[2][iter];
873 Float_t zcell = zhit / fGeom->GetPadSizeZ();
874 Float_t xcell = xhit / fGeom->GetPadSizePhi();
875 if ( zcell<=0 || xcell<=0 ||
876 zcell>=nCellZ || xcell>=nCellX) return;
877 Int_t izcell = (Int_t) zcell;
878 Int_t ixcell = (Int_t) xcell;
879 Float_t zc = zcell - izcell - 0.5;
880 Float_t xc = xcell - ixcell - 0.5;
a3dfe79c 881 for (Int_t iz=1; iz<=kNgamz; iz++) {
fa412d9b 882 Int_t kzg = izcell + iz - nz3;
883 if (kzg<=0 || kzg>nCellZ) continue;
884 Float_t zg = (Float_t)(iz-nz3) - zc;
a3dfe79c 885 for (Int_t ix=1; ix<=kNgamx; ix++) {
fa412d9b 886 Int_t kxg = ixcell + ix - nx3;
887 if (kxg<=0 || kxg>nCellX) continue;
888 Float_t xg = (Float_t)(ix-nx3) - xc;
889
890 // Now calculate pad response
891 Float_t qpad = CPVPadResponseFunction(qhit,zg,xg);
a3dfe79c 892 qpad += kNoise*rnor2;
fa412d9b 893 if (qpad<0) continue;
894
895 // Fill the array with pad response ID and amplitude
3d402178 896 new(ldigits[cpvDigits->GetEntriesFast()]) AliPHOSCPVDigit(kxg,kzg,qpad);
fa412d9b 897 }
fa412d9b 898 }
fa412d9b 899 }
900}
901
902//____________________________________________________________________________
903Float_t AliPHOSv1::CPVPadResponseFunction(Float_t qhit, Float_t zhit, Float_t xhit) {
904 // ------------------------------------------------------------------------
905 // Calculate the amplitude in one CPV pad using the
906 // cumulative pad response function
907 // Author: Yuri Kharlov (after Serguei Sadovski)
908 // 3 October 2000
909 // ------------------------------------------------------------------------
910
911 Double_t dz = fGeom->GetPadSizeZ() / 2;
912 Double_t dx = fGeom->GetPadSizePhi() / 2;
913 Double_t z = zhit * fGeom->GetPadSizeZ();
914 Double_t x = xhit * fGeom->GetPadSizePhi();
915 Double_t amplitude = qhit *
916 (CPVCumulPadResponse(z+dz,x+dx) - CPVCumulPadResponse(z+dz,x-dx) -
917 CPVCumulPadResponse(z-dz,x+dx) + CPVCumulPadResponse(z-dz,x-dx));
918 return (Float_t)amplitude;
7587f5a5 919}
920
fa412d9b 921//____________________________________________________________________________
922Double_t AliPHOSv1::CPVCumulPadResponse(Double_t x, Double_t y) {
923 // ------------------------------------------------------------------------
924 // Cumulative pad response function
925 // It includes several terms from the CF decomposition in electrostatics
926 // Note: this cumulative function is wrong since omits some terms
927 // but the cell amplitude obtained with it is correct because
928 // these omitting terms cancel
929 // Author: Yuri Kharlov (after Serguei Sadovski)
930 // 3 October 2000
931 // ------------------------------------------------------------------------
932
a3dfe79c 933 const Double_t kA=1.0;
934 const Double_t kB=0.7;
fa412d9b 935
936 Double_t r2 = x*x + y*y;
937 Double_t xy = x*y;
938 Double_t cumulPRF = 0;
939 for (Int_t i=0; i<=4; i++) {
a3dfe79c 940 Double_t b1 = (2*i + 1) * kB;
fa412d9b 941 cumulPRF += TMath::Power(-1,i) * TMath::ATan( xy / (b1*TMath::Sqrt(b1*b1 + r2)) );
942 }
a3dfe79c 943 cumulPRF *= kA/(2*TMath::Pi());
fa412d9b 944 return cumulPRF;
945}
7eb9d12d 946