1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //_________________________________________________________________________
19 // Implementation version v1 of PHOS Manager class
21 // Layout EMC + PPSD has name GPS2:
22 // Produces cumulated hits
24 // Layout EMC + CPV has name IHEP:
25 // Produces hits for CPV, cumulated hits
27 // Layout EMC + CPV + PPSD has name GPS:
28 // Produces hits for CPV, cumulated hits
30 //*-- Author: Yves Schutz (SUBATECH)
33 // --- ROOT system ---
41 // --- Standard library ---
46 #include <strstream.h>
48 // --- AliRoot header files ---
50 #include "AliPHOSv1.h"
51 #include "AliPHOSHit.h"
52 #include "AliPHOSDigit.h"
53 #include "AliPHOSReconstructioner.h"
60 //____________________________________________________________________________
61 AliPHOSv1::AliPHOSv1()
65 // Create an empty array of AliPHOSCPVModule to satisfy
66 // AliPHOSv1::Streamer when reading root file
68 fReconstructioner = 0;
69 fTrackSegmentMaker = 0;
71 if ( 0==(fEMCModules=new TClonesArray("AliPHOSCPVModule",0)) ) {
72 Error("AliPHOSv1","Can not create array of EMC modules");
76 if ( 0==(fCPVModules=new TClonesArray("AliPHOSCPVModule",0)) ) {
77 Error("AliPHOSv1","Can not create array of CPV modules");
83 //____________________________________________________________________________
84 AliPHOSv1::AliPHOSv1(const char *name, const char *title):
87 // ctor : title is used to identify the layout
88 // GPS2 = 5 modules (EMC + PPSD)
89 // IHEP = 5 modules (EMC + CPV )
90 // MIXT = 4 modules (EMC + CPV ) and 1 module (EMC + PPSD)
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).
98 fPinElectronicNoise = 0.010 ;
99 fDigitThreshold = 0.1 ; // 1 GeV
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
104 // and the TreeD at the end of the event (branch is set in FinishEvent() ).
106 fHits= new TClonesArray("AliPHOSHit",1000) ;
110 fReconstructioner = 0;
111 fTrackSegmentMaker = 0;
113 fIshunt = 1 ; // All hits are associated with primary particles
115 // Create array of EMC modules of the size of PHOS modules number
117 if ( 0==(fEMCModules=new TClonesArray("AliPHOSCPVModule",fGeom->GetNModules())) ) {
118 Error("AliPHOSv1","Can not create array of EMC modules");
121 TClonesArray &lemcmodule = *fEMCModules;
122 for (Int_t i=0; i<fGeom->GetNModules(); i++) new(lemcmodule[i]) AliPHOSCPVModule();
124 // Create array of CPV modules for the IHEP's version of CPV
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
129 if ( 0==(fCPVModules=new TClonesArray("AliPHOSCPVModule",fGeom->GetNCPVModules())) ) {
130 Error("AliPHOSv1","Can not create array of CPV modules");
133 TClonesArray &lcpvmodule = *fCPVModules;
134 for (Int_t i=0; i<fGeom->GetNCPVModules(); i++) new(lcpvmodule[i]) AliPHOSCPVModule();
137 // Create an empty array of AliPHOSCPVModule to satisfy
138 // AliPHOSv1::Streamer when writing root file
140 fCPVModules=new TClonesArray("AliPHOSCPVModule",0);
145 //____________________________________________________________________________
146 AliPHOSv1::AliPHOSv1(AliPHOSReconstructioner * Reconstructioner, const char *name, const char *title):
147 AliPHOSv0(name,title)
149 // ctor : title is used to identify the layout
150 // GPS2 = 5 modules (EMC + PPSD)
151 // We use 2 arrays of hits :
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).
157 // - fTmpHits, which retains all the hits of the current event. It
158 // is used for the digitization part.
160 fPinElectronicNoise = 0.010 ;
162 // We do not want to save in TreeH the raw hits
163 //fHits = new TClonesArray("AliPHOSHit",100) ;
166 fHits= new TClonesArray("AliPHOSHit",1000) ;
170 fIshunt = 1 ; // All hits are associated with primary particles
172 // gets an instance of the geometry parameters class
173 fGeom = AliPHOSGeometry::GetInstance(title, "") ;
175 if (fGeom->IsInitialized() )
176 cout << "AliPHOS" << Version() << " : PHOS geometry intialized for " << fGeom->GetName() << endl ;
178 cout << "AliPHOS" << Version() << " : PHOS geometry initialization failed !" << endl ;
180 // Defining the PHOS Reconstructioner
182 fReconstructioner = Reconstructioner ;
186 //____________________________________________________________________________
187 AliPHOSv1::~AliPHOSv1()
203 if ( fEmcRecPoints ) {
204 fEmcRecPoints->Delete() ;
205 delete fEmcRecPoints ;
209 if ( fPpsdRecPoints ) {
210 fPpsdRecPoints->Delete() ;
211 delete fPpsdRecPoints ;
215 if ( fTrackSegments ) {
216 fTrackSegments->Delete() ;
217 delete fTrackSegments ;
223 //____________________________________________________________________________
224 void AliPHOSv1::AddHit(Int_t shunt, Int_t primary, Int_t tracknumber, Int_t Id, Float_t * hits, Int_t trackpid)
226 // Add a hit to the hit list.
227 // A PHOS hit is the sum of all hits in a single crystal
228 // or in a single PPSD gas cell
233 Bool_t deja = kFALSE ;
235 newHit = new AliPHOSHit(shunt, primary, tracknumber, Id, hits, trackpid) ;
237 for ( hitCounter = 0 ; hitCounter < fNhits && !deja ; hitCounter++ ) {
238 curHit = (AliPHOSHit*) (*fHits)[hitCounter] ;
239 if( *curHit == *newHit ) {
240 *curHit = *curHit + *newHit ;
246 new((*fHits)[fNhits]) AliPHOSHit(*newHit) ;
253 //___________________________________________________________________________
254 Int_t AliPHOSv1::Digitize(Float_t Energy)
256 // Applies the energy calibration
258 Float_t fB = 100000000. ;
260 Int_t chan = Int_t(fA + Energy*fB ) ;
264 //____________________________________________________________________________
265 void AliPHOSv1::Hit2Digit(Int_t ntracks){
266 //Collects all hits in the same active volume into digits
271 fDigits = new TClonesArray("AliPHOSDigit",1000) ;
273 // Branch address for digit tree
275 sprintf(branchname,"%s",GetName());
276 gAlice->TreeD()->Branch(branchname,&fDigits,fBufferSize);
278 gAlice->TreeD()->GetEvent(0);
285 AliPHOSDigit * newdigit ;
286 AliPHOSDigit * curdigit ;
287 Bool_t deja = kFALSE ;
290 for (itrack=0; itrack<ntracks; itrack++){
292 //=========== Get the Hits Tree for the Primary track itrack
294 gAlice->TreeH()->GetEvent(itrack);
296 for ( i = 0 ; i < fHits->GetEntries() ; i++ ) {
297 hit = (AliPHOSHit*)fHits->At(i) ;
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() ) ) ;
303 newdigit = new AliPHOSDigit( -1 , hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
307 for ( j = 0 ; j < fNdigits ; j++) {
308 curdigit = (AliPHOSDigit*) fDigits->At(j) ;
309 if ( *curdigit == *newdigit) {
310 *curdigit = *curdigit + *newdigit ;
316 new((*fDigits)[fNdigits]) AliPHOSDigit(* newdigit) ;
323 } // loop over tracks
325 // Noise induced by the PIN diode of the PbWO crystals
327 Float_t energyandnoise ;
328 for ( i = 0 ; i < fNdigits ; i++ ) {
329 newdigit = (AliPHOSDigit * ) fDigits->At(i) ;
331 fGeom->AbsToRelNumbering(newdigit->GetId(), relid) ;
333 if (relid[1]==0){ // Digits belong to EMC (PbW0_4 crystals)
334 energyandnoise = newdigit->GetAmp() + Digitize(gRandom->Gaus(0., fPinElectronicNoise)) ;
336 if (energyandnoise < 0 )
339 if ( newdigit->GetAmp() < fDigitThreshold ) // if threshold not surpassed, remove digit from list
340 fDigits->RemoveAt(i) ;
344 fDigits->Compress() ;
346 fNdigits = fDigits->GetEntries() ;
347 fDigits->Expand(fNdigits) ;
349 for (i = 0 ; i < fNdigits ; i++) {
350 newdigit = (AliPHOSDigit *) fDigits->At(i) ;
351 newdigit->SetIndexInList(i) ;
354 gAlice->TreeD()->Fill() ;
356 gAlice->TreeD()->Write(0,TObject::kOverwrite) ;
359 //___________________________________________________________________________
360 void AliPHOSv1::MakeBranch(Option_t* opt)
362 // Create new branche in the current Root Tree in the digit Tree
363 AliDetector::MakeBranch(opt) ;
365 // Create new branches EMC<i> for hits in EMC modules
367 for( Int_t i=0; i<fGeom->GetNModules(); i++ ) GetEMCModule(i).MakeBranch("EMC",i+1);
369 // Create new branches CPV<i> for hits in CPV modules for IHEP geometry
371 if ( strcmp(fGeom->GetName(),"IHEP") == 0 || strcmp(fGeom->GetName(),"MIXT") == 0 ) {
372 for( Int_t i=0; i<fGeom->GetNCPVModules(); i++ ) GetCPVModule(i).MakeBranch("CPV",i+1);
377 //_____________________________________________________________________________
378 void AliPHOSv1::Reconstruction(AliPHOSReconstructioner * Reconstructioner)
380 // 1. Reinitializes the existing RecPoint, TrackSegment, and RecParticles Lists and
381 // 2. Creates TreeR with a branch for each list
382 // 3. Steers the reconstruction processes
383 // 4. Saves the 3 lists in TreeR
384 // 5. Write the Tree to File
386 fReconstructioner = Reconstructioner ;
388 char branchname[10] ;
392 // gAlice->MakeTree("R") ;
393 Int_t splitlevel = 0 ;
395 fEmcRecPoints->Delete() ;
397 if ( fEmcRecPoints && gAlice->TreeR() ) {
398 sprintf(branchname,"%sEmcRP",GetName()) ;
399 gAlice->TreeR()->Branch(branchname, "TObjArray", &fEmcRecPoints, fBufferSize, splitlevel) ;
402 fPpsdRecPoints->Delete() ;
404 if ( fPpsdRecPoints && gAlice->TreeR() ) {
405 sprintf(branchname,"%sPpsdRP",GetName()) ;
406 gAlice->TreeR()->Branch(branchname, "TObjArray", &fPpsdRecPoints, fBufferSize, splitlevel) ;
409 fTrackSegments->Delete() ;
411 if ( fTrackSegments && gAlice->TreeR() ) {
412 sprintf(branchname,"%sTS",GetName()) ;
413 gAlice->TreeR()->Branch(branchname, &fTrackSegments, fBufferSize) ;
416 fRecParticles->Delete() ;
418 if ( fRecParticles && gAlice->TreeR() ) {
419 sprintf(branchname,"%sRP",GetName()) ;
420 gAlice->TreeR()->Branch(branchname, &fRecParticles, fBufferSize) ;
425 fReconstructioner->Make(fDigits, fEmcRecPoints, fPpsdRecPoints, fTrackSegments, fRecParticles);
427 printf("Reconstruction: %d %d %d %d\n",
428 fEmcRecPoints->GetEntries(),fPpsdRecPoints->GetEntries(),
429 fTrackSegments->GetEntries(),fRecParticles->GetEntries());
431 // 4. Expand or Shrink the arrays to the proper size
435 size = fEmcRecPoints->GetEntries() ;
436 fEmcRecPoints->Expand(size) ;
438 size = fPpsdRecPoints->GetEntries() ;
439 fPpsdRecPoints->Expand(size) ;
441 size = fTrackSegments->GetEntries() ;
442 fTrackSegments->Expand(size) ;
444 size = fRecParticles->GetEntries() ;
445 fRecParticles->Expand(size) ;
447 gAlice->TreeR()->Fill() ;
450 gAlice->TreeR()->Write(0,TObject::kOverwrite) ;
452 // Deleting reconstructed objects
453 ResetReconstruction();
457 //____________________________________________________________________________
458 void AliPHOSv1::ResetHits()
460 // Reset hit tree for EMC and CPV
461 // Yuri Kharlov, 28 September 2000
463 AliDetector::ResetHits();
464 for (Int_t i=0; i<fGeom->GetNModules(); i++) ((AliPHOSCPVModule*)(*fEMCModules)[i]) -> Clear();
465 if ( strcmp(fGeom->GetName(),"IHEP") == 0 || strcmp(fGeom->GetName(),"MIXT") == 0 ) {
466 for (Int_t i=0; i<fGeom->GetNCPVModules(); i++) ((AliPHOSCPVModule*)(*fCPVModules)[i]) -> Clear();
470 //____________________________________________________________________________
471 void AliPHOSv1::ResetReconstruction()
473 // Deleting reconstructed objects
475 if ( fEmcRecPoints ) fEmcRecPoints ->Delete();
476 if ( fPpsdRecPoints ) fPpsdRecPoints->Delete();
477 if ( fTrackSegments ) fTrackSegments->Delete();
478 if ( fRecParticles ) fRecParticles ->Delete();
482 //____________________________________________________________________________
483 void AliPHOSv1::SetTreeAddress()
486 AliPHOS::SetTreeAddress();
488 // //Branch address for TreeR: RecPpsdRecPoint
489 // TTree *treeR = gAlice->TreeR();
490 // if ( treeR && fPpsdRecPoints ) {
491 // branch = treeR->GetBranch("PHOSPpsdRP");
492 // if (branch) branch->SetAddress(&fPpsdRecPoints) ;
495 // Set branch address for the Hits Tree for hits in EMC modules
496 // Yuri Kharlov, 23 November 2000.
498 for( Int_t i=0; i<fGeom->GetNModules(); i++ ) GetEMCModule(i).SetTreeAddress("EMC",i+1);
500 // Set branch address for the Hits Tree for hits in CPV modules for IHEP geometry
501 // Yuri Kharlov, 28 September 2000.
503 if ( strcmp(fGeom->GetName(),"IHEP") == 0 || strcmp(fGeom->GetName(),"MIXT") == 0 ) {
504 for( Int_t i=0; i<fGeom->GetNCPVModules(); i++ ) GetCPVModule(i).SetTreeAddress("CPV",i+1);
509 //____________________________________________________________________________
511 void AliPHOSv1::StepManager(void)
513 // Accumulates hits as long as the track stays in a single crystal or PPSD gas Cell
515 // if (gMC->IsTrackEntering())
516 // cout << "Track enters the volume " << gMC->CurrentVolName() << endl;
517 // if (gMC->IsTrackExiting())
518 // cout << "Track leaves the volume " << gMC->CurrentVolName() << endl;
520 Int_t relid[4] ; // (box, layer, row, column) indices
521 Int_t absid ; // absolute cell ID number
522 Float_t xyze[4]={0,0,0,0} ; // position wrt MRS and energy deposited
523 TLorentzVector pos ; // Lorentz vector of the track current position
527 Int_t tracknumber = gAlice->CurrentTrack() ;
528 Int_t primary = gAlice->GetPrimary( gAlice->CurrentTrack() );
529 TString name = fGeom->GetName() ;
530 Int_t trackpid = gMC->TrackPid() ;
531 if ( name == "GPS2" || name == "MIXT" ) { // ======> CPV is a GPS' PPSD
533 if( gMC->CurrentVolID(copy) == gMC->VolId("GCEL") ) // We are inside a gas cell
535 gMC->TrackPosition(pos) ;
539 xyze[3] = gMC->Edep() ;
541 if ( xyze[3] != 0 ) { // there is deposited energy
542 gMC->CurrentVolOffID(5, relid[0]) ; // get the PHOS Module number
543 if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(5),"PHO1") == 0 ){
544 relid[0] += fGeom->GetNModules() - fGeom->GetNPPSDModules();
546 gMC->CurrentVolOffID(3, relid[1]) ; // get the Micromegas Module number
547 // 1-> fGeom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() upper
548 // > fGeom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() lower
549 gMC->CurrentVolOffID(1, relid[2]) ; // get the row number of the cell
550 gMC->CurrentVolID(relid[3]) ; // get the column number
552 // get the absolute Id number
554 fGeom->RelToAbsNumbering(relid, absid) ;
556 // add current hit to the hit list
557 AddHit(fIshunt, primary, tracknumber, absid, xyze, trackpid);
559 } // there is deposited energy
560 } // We are inside the gas of the CPV
561 } // GPS2 configuration
563 if ( name == "IHEP" || name == "MIXT" ) { // ======> CPV is a IHEP's one
565 // Yuri Kharlov, 28 September 2000
567 if( gMC->CurrentVolID(copy) == gMC->VolId("CPVQ") &&
568 gMC->IsTrackEntering() &&
569 gMC->TrackCharge() != 0) {
571 // Charged track has just entered to the CPV sensitive plane
573 AliPHOSv1 &phos = *(AliPHOSv1*)gAlice->GetModule("PHOS");
576 gMC->CurrentVolOffID(3,moduleNumber);
579 // Current position of the hit in the CPV module ref. system
581 gMC -> TrackPosition(pos);
582 Float_t xyzm[3], xyzd[3], xyd[3]={0,0,0};
584 for (i=0; i<3; i++) xyzm[i] = pos[i];
585 gMC -> Gmtod (xyzm, xyzd, 1); // transform coordinate from master to daughter system
589 // Current momentum of the hit's track in the CPV module ref. system
592 gMC -> TrackMomentum(pmom);
593 Float_t pm[3], pd[3];
594 for (i=0; i<3; i++) pm[i] = pmom[i];
595 gMC -> Gmtod (pm, pd, 2); // transform 3-momentum from master to daughter system
600 // Current particle type of the hit's track
602 Int_t ipart = gMC->TrackPid();
604 // Add the current particle in the list of the CPV hits.
606 phos.GetCPVModule(moduleNumber).AddHit(fIshunt,primary,pmom,xyd,ipart);
609 printf("CPV hit added to module #%2d: p = (% .4f, % .4f, % .4f, % .4f) GeV,\n",
610 moduleNumber+1,pmom.Px(),pmom.Py(),pmom.Pz(),pmom.E());
611 printf( " xy = (%8.4f, %8.4f) cm, ipart = %d, primary = %d\n",
612 xyd[0],xyd[1],ipart,primary);
615 // Digitize the current CPV hit:
617 // 1. find pad response and
619 TClonesArray *cpvDigits = new TClonesArray("AliPHOSCPVDigit",0); // array of digits for current hit
620 CPVDigitize(pmom,xyd,moduleNumber,cpvDigits);
625 Int_t idigit,ndigits;
627 // 2. go through the current digit list and sum digits in pads
629 ndigits = cpvDigits->GetEntriesFast();
630 for (idigit=0; idigit<ndigits-1; idigit++) {
631 AliPHOSCPVDigit *cpvDigit1 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit);
632 Float_t x1 = cpvDigit1->GetXpad() ;
633 Float_t z1 = cpvDigit1->GetYpad() ;
634 for (Int_t jdigit=idigit+1; jdigit<ndigits; jdigit++) {
635 AliPHOSCPVDigit *cpvDigit2 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(jdigit);
636 Float_t x2 = cpvDigit2->GetXpad() ;
637 Float_t z2 = cpvDigit2->GetYpad() ;
638 if (x1==x2 && z1==z2) {
639 Float_t qsum = cpvDigit1->GetQpad() + cpvDigit2->GetQpad() ;
640 cpvDigit2->SetQpad(qsum) ;
641 cpvDigits->RemoveAt(idigit) ;
645 cpvDigits->Compress() ;
647 // 3. add digits to temporary hit list fTmpHits
649 ndigits = cpvDigits->GetEntriesFast();
650 for (idigit=0; idigit<ndigits; idigit++) {
651 AliPHOSCPVDigit *cpvDigit = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit);
652 relid[0] = moduleNumber + 1 ; // CPV (or PHOS) module number
653 relid[1] =-1 ; // means CPV
654 relid[2] = cpvDigit->GetXpad() ; // column number of a pad
655 relid[3] = cpvDigit->GetYpad() ; // row number of a pad
657 // get the absolute Id number
658 fGeom->RelToAbsNumbering(relid, absid) ;
660 // add current digit to the temporary hit list
664 xyze[3] = cpvDigit->GetQpad() ; // amplitude in a pad
665 primary = -1; // No need in primary for CPV
666 AddHit(fIshunt, primary, tracknumber, absid, xyze, trackpid);
668 if (cpvDigit->GetQpad() > 0.02) {
669 xmean += cpvDigit->GetQpad() * (cpvDigit->GetXpad() + 0.5);
670 zmean += cpvDigit->GetQpad() * (cpvDigit->GetYpad() + 0.5);
671 qsum += cpvDigit->GetQpad();
676 } // end of IHEP configuration
678 if(gMC->CurrentVolID(copy) == gMC->VolId("PXTL") ) { // We are inside a PBWO crystal
679 gMC->TrackPosition(pos) ;
683 xyze[3] = gMC->Edep() ;
685 // Track enters to the crystal from the top edge
687 if (gMC->IsTrackEntering()) {
689 gMC -> Gmtod (xyze, posloc, 1);
690 if (posloc[1] > fGeom->GetCrystalSize(1)/2-0.01) {
692 Float_t xyd[3]={0,0,0};
693 AliPHOSv1 &phos = *(AliPHOSv1*)gAlice->GetModule("PHOS");
696 gMC->CurrentVolOffID(10,moduleNumber);
698 if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(10),"PHO1") == 0 )
699 moduleNumber += fGeom->GetNModules() - fGeom->GetNPPSDModules();
702 gMC->CurrentVolOffID(4, row) ; // get the row number inside the module
703 gMC->CurrentVolOffID(3, cel) ; // get the cell number inside the module
704 xyd[0] = -(posloc[2] + (cel-0.5-fGeom->GetNZ() /2) *
705 (fGeom->GetCrystalSize(2) + 2 * fGeom->GetGapBetweenCrystals()));
706 xyd[1] = posloc[0] + (row-0.5-fGeom->GetNPhi()/2) *
707 (fGeom->GetCrystalSize(0) + 2 * fGeom->GetGapBetweenCrystals());
709 // Current momentum of the hit's track in the CPV module ref. system
712 gMC -> TrackMomentum(pmom);
713 Float_t pm[3], pd[3];
714 for (i=0; i<3; i++) pm[i] = pmom[i];
715 gMC -> Gmtod (pm, pd, 2); // transform 3-momentum from master to daughter system
720 // Current particle type of the hit's track
722 Int_t ipart = gMC->TrackPid();
724 // Add the current particle in the list of the EMC hits.
726 phos.GetEMCModule(moduleNumber).AddHit(fIshunt,primary,pmom,xyd,ipart);
729 printf("EMC hit added to module #%2d: p = (% .4f, % .4f, % .4f, % .4f) GeV,\n",
730 moduleNumber+1,pmom.Px(),pmom.Py(),pmom.Pz(),pmom.E());
731 printf( " xy = (%8.4f, %8.4f) cm, ipart = %d, primary = %d\n",
732 xyd[0],xyd[1],ipart,primary);
737 // Track is inside the crystal and deposits some energy
739 if ( xyze[3] != 0 ) {
740 gMC->CurrentVolOffID(10, relid[0]) ; // get the PHOS module number ;
741 if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(10),"PHO1") == 0 )
742 relid[0] += fGeom->GetNModules() - fGeom->GetNPPSDModules();
744 relid[1] = 0 ; // means PBW04
745 gMC->CurrentVolOffID(4, relid[2]) ; // get the row number inside the module
746 gMC->CurrentVolOffID(3, relid[3]) ; // get the cell number inside the module
748 // get the absolute Id number
749 fGeom->RelToAbsNumbering(relid, absid) ;
751 // add current hit to the hit list
752 AddHit(fIshunt, primary,tracknumber, absid, xyze, trackpid);
754 } // there is deposited energy
755 } // we are inside a PHOS Xtal
758 //____________________________________________________________________________
759 void AliPHOSv1::CPVDigitize (TLorentzVector p, Float_t *zxhit, Int_t moduleNumber, TClonesArray *cpvDigits)
761 // ------------------------------------------------------------------------
762 // Digitize one CPV hit:
763 // On input take exact 4-momentum p and position zxhit of the hit,
764 // find the pad response around this hit and
765 // put the amplitudes in the pads into array digits
767 // Author: Yuri Kharlov (after Serguei Sadovsky)
769 // ------------------------------------------------------------------------
771 const Float_t kCelWr = fGeom->GetPadSizePhi()/2; // Distance between wires (2 wires above 1 pad)
772 const Float_t kDetR = 0.1; // Relative energy fluctuation in track for 100 e-
773 const Float_t kdEdx = 4.0; // Average energy loss in CPV;
774 const Int_t kNgamz = 5; // Ionization size in Z
775 const Int_t kNgamx = 9; // Ionization size in Phi
776 const Float_t kNoise = 0.03; // charge noise in one pad
780 // Just a reminder on axes notation in the CPV module:
781 // axis Z goes along the beam
782 // axis X goes across the beam in the module plane
783 // axis Y is a normal to the module plane showing from the IP
785 Float_t hitX = zxhit[0];
786 Float_t hitZ =-zxhit[1];
789 Float_t pNorm = p.Py();
790 Float_t eloss = kdEdx;
792 // cout << "CPVDigitize: YVK : "<<hitX<<" "<<hitZ<<" | "<<pX<<" "<<pZ<<" "<<pNorm<<endl;
794 Float_t dZY = pZ/pNorm * fGeom->GetCPVGasThickness();
795 Float_t dXY = pX/pNorm * fGeom->GetCPVGasThickness();
796 gRandom->Rannor(rnor1,rnor2);
797 eloss *= (1 + kDetR*rnor1) *
798 TMath::Sqrt((1 + ( pow(dZY,2) + pow(dXY,2) ) / pow(fGeom->GetCPVGasThickness(),2)));
799 Float_t zhit1 = hitZ + fGeom->GetCPVActiveSize(1)/2 - dZY/2;
800 Float_t xhit1 = hitX + fGeom->GetCPVActiveSize(0)/2 - dXY/2;
801 Float_t zhit2 = zhit1 + dZY;
802 Float_t xhit2 = xhit1 + dXY;
804 Int_t iwht1 = (Int_t) (xhit1 / kCelWr); // wire (x) coordinate "in"
805 Int_t iwht2 = (Int_t) (xhit2 / kCelWr); // wire (x) coordinate "out"
809 if (iwht1==iwht2) { // incline 1-wire hit
811 zxe[0][0] = (zhit1 + zhit2 - dZY*0.57735) / 2;
812 zxe[1][0] = (iwht1 + 0.5) * kCelWr;
814 zxe[0][1] = (zhit1 + zhit2 + dZY*0.57735) / 2;
815 zxe[1][1] = (iwht1 + 0.5) * kCelWr;
818 else if (TMath::Abs(iwht1-iwht2) != 1) { // incline 3-wire hit
820 Int_t iwht3 = (iwht1 + iwht2) / 2;
821 Float_t xwht1 = (iwht1 + 0.5) * kCelWr; // wire 1
822 Float_t xwht2 = (iwht2 + 0.5) * kCelWr; // wire 2
823 Float_t xwht3 = (iwht3 + 0.5) * kCelWr; // wire 3
824 Float_t xwr13 = (xwht1 + xwht3) / 2; // center 13
825 Float_t xwr23 = (xwht2 + xwht3) / 2; // center 23
826 Float_t dxw1 = xhit1 - xwr13;
827 Float_t dxw2 = xhit2 - xwr23;
828 Float_t egm1 = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
829 Float_t egm2 = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
830 Float_t egm3 = kCelWr / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
831 zxe[0][0] = (dXY*(xwr13-xwht1)/dXY + zhit1 + zhit1) / 2;
833 zxe[2][0] = eloss * egm1;
834 zxe[0][1] = (dXY*(xwr23-xwht1)/dXY + zhit1 + zhit2) / 2;
836 zxe[2][1] = eloss * egm2;
837 zxe[0][2] = dXY*(xwht3-xwht1)/dXY + zhit1;
839 zxe[2][2] = eloss * egm3;
841 else { // incline 2-wire hit
843 Float_t xwht1 = (iwht1 + 0.5) * kCelWr;
844 Float_t xwht2 = (iwht2 + 0.5) * kCelWr;
845 Float_t xwr12 = (xwht1 + xwht2) / 2;
846 Float_t dxw1 = xhit1 - xwr12;
847 Float_t dxw2 = xhit2 - xwr12;
848 Float_t egm1 = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) );
849 Float_t egm2 = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) );
850 zxe[0][0] = (zhit1 + zhit2 - dZY*egm1) / 2;
852 zxe[2][0] = eloss * egm1;
853 zxe[0][1] = (zhit1 + zhit2 + dZY*egm2) / 2;
855 zxe[2][1] = eloss * egm2;
858 // Finite size of ionization region
860 Int_t nCellZ = fGeom->GetNumberOfCPVPadsZ();
861 Int_t nCellX = fGeom->GetNumberOfCPVPadsPhi();
862 Int_t nz3 = (kNgamz+1)/2;
863 Int_t nx3 = (kNgamx+1)/2;
864 cpvDigits->Expand(nIter*kNgamx*kNgamz);
865 TClonesArray &ldigits = *(TClonesArray *)cpvDigits;
867 for (Int_t iter=0; iter<nIter; iter++) {
869 Float_t zhit = zxe[0][iter];
870 Float_t xhit = zxe[1][iter];
871 Float_t qhit = zxe[2][iter];
872 Float_t zcell = zhit / fGeom->GetPadSizeZ();
873 Float_t xcell = xhit / fGeom->GetPadSizePhi();
874 if ( zcell<=0 || xcell<=0 ||
875 zcell>=nCellZ || xcell>=nCellX) return;
876 Int_t izcell = (Int_t) zcell;
877 Int_t ixcell = (Int_t) xcell;
878 Float_t zc = zcell - izcell - 0.5;
879 Float_t xc = xcell - ixcell - 0.5;
880 for (Int_t iz=1; iz<=kNgamz; iz++) {
881 Int_t kzg = izcell + iz - nz3;
882 if (kzg<=0 || kzg>nCellZ) continue;
883 Float_t zg = (Float_t)(iz-nz3) - zc;
884 for (Int_t ix=1; ix<=kNgamx; ix++) {
885 Int_t kxg = ixcell + ix - nx3;
886 if (kxg<=0 || kxg>nCellX) continue;
887 Float_t xg = (Float_t)(ix-nx3) - xc;
889 // Now calculate pad response
890 Float_t qpad = CPVPadResponseFunction(qhit,zg,xg);
891 qpad += kNoise*rnor2;
892 if (qpad<0) continue;
894 // Fill the array with pad response ID and amplitude
895 new(ldigits[cpvDigits->GetEntriesFast()]) AliPHOSCPVDigit(kxg,kzg,qpad);
901 //____________________________________________________________________________
902 Float_t AliPHOSv1::CPVPadResponseFunction(Float_t qhit, Float_t zhit, Float_t xhit) {
903 // ------------------------------------------------------------------------
904 // Calculate the amplitude in one CPV pad using the
905 // cumulative pad response function
906 // Author: Yuri Kharlov (after Serguei Sadovski)
908 // ------------------------------------------------------------------------
910 Double_t dz = fGeom->GetPadSizeZ() / 2;
911 Double_t dx = fGeom->GetPadSizePhi() / 2;
912 Double_t z = zhit * fGeom->GetPadSizeZ();
913 Double_t x = xhit * fGeom->GetPadSizePhi();
914 Double_t amplitude = qhit *
915 (CPVCumulPadResponse(z+dz,x+dx) - CPVCumulPadResponse(z+dz,x-dx) -
916 CPVCumulPadResponse(z-dz,x+dx) + CPVCumulPadResponse(z-dz,x-dx));
917 return (Float_t)amplitude;
920 //____________________________________________________________________________
921 Double_t AliPHOSv1::CPVCumulPadResponse(Double_t x, Double_t y) {
922 // ------------------------------------------------------------------------
923 // Cumulative pad response function
924 // It includes several terms from the CF decomposition in electrostatics
925 // Note: this cumulative function is wrong since omits some terms
926 // but the cell amplitude obtained with it is correct because
927 // these omitting terms cancel
928 // Author: Yuri Kharlov (after Serguei Sadovski)
930 // ------------------------------------------------------------------------
932 const Double_t kA=1.0;
933 const Double_t kB=0.7;
935 Double_t r2 = x*x + y*y;
937 Double_t cumulPRF = 0;
938 for (Int_t i=0; i<=4; i++) {
939 Double_t b1 = (2*i + 1) * kB;
940 cumulPRF += TMath::Power(-1,i) * TMath::ATan( xy / (b1*TMath::Sqrt(b1*b1 + r2)) );
942 cumulPRF *= kA/(2*TMath::Pi());