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