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 = 10000000. ;
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) ;
360 //___________________________________________________________________________
361 void AliPHOSv1::MakeBranch(Option_t* opt, char *file)
363 // Create new branche in the current Root Tree in the digit Tree
364 AliDetector::MakeBranch(opt) ;
366 // Create new branches EMC<i> for hits in EMC modules
368 for( Int_t i=0; i<fGeom->GetNModules(); i++ ) GetEMCModule(i).MakeBranch("EMC",i+1,file);
370 // Create new branches CPV<i> for hits in CPV modules for IHEP geometry
372 if ( strcmp(fGeom->GetName(),"IHEP") == 0 || strcmp(fGeom->GetName(),"MIXT") == 0 ) {
373 for( Int_t i=0; i<fGeom->GetNCPVModules(); i++ ) GetCPVModule(i).MakeBranch("CPV",i+1,file);
378 //_____________________________________________________________________________
379 void AliPHOSv1::Reconstruction(AliPHOSReconstructioner * Reconstructioner)
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
387 fReconstructioner = Reconstructioner ;
389 char branchname[10] ;
393 // gAlice->MakeTree("R") ;
394 Int_t splitlevel = 0 ;
396 fEmcRecPoints->Delete() ;
398 if ( fEmcRecPoints && gAlice->TreeR() ) {
399 sprintf(branchname,"%sEmcRP",GetName()) ;
400 gAlice->TreeR()->Branch(branchname, "TObjArray", &fEmcRecPoints, fBufferSize, splitlevel) ;
403 fPpsdRecPoints->Delete() ;
405 if ( fPpsdRecPoints && gAlice->TreeR() ) {
406 sprintf(branchname,"%sPpsdRP",GetName()) ;
407 gAlice->TreeR()->Branch(branchname, "TObjArray", &fPpsdRecPoints, fBufferSize, splitlevel) ;
410 fTrackSegments->Delete() ;
412 if ( fTrackSegments && gAlice->TreeR() ) {
413 sprintf(branchname,"%sTS",GetName()) ;
414 gAlice->TreeR()->Branch(branchname, &fTrackSegments, fBufferSize) ;
417 fRecParticles->Delete() ;
419 if ( fRecParticles && gAlice->TreeR() ) {
420 sprintf(branchname,"%sRP",GetName()) ;
421 gAlice->TreeR()->Branch(branchname, &fRecParticles, fBufferSize) ;
426 fReconstructioner->Make(fDigits, fEmcRecPoints, fPpsdRecPoints, fTrackSegments, fRecParticles);
428 printf("Reconstruction: %d %d %d %d\n",
429 fEmcRecPoints->GetEntries(),fPpsdRecPoints->GetEntries(),
430 fTrackSegments->GetEntries(),fRecParticles->GetEntries());
432 // 4. Expand or Shrink the arrays to the proper size
436 size = fEmcRecPoints->GetEntries() ;
437 fEmcRecPoints->Expand(size) ;
439 size = fPpsdRecPoints->GetEntries() ;
440 fPpsdRecPoints->Expand(size) ;
442 size = fTrackSegments->GetEntries() ;
443 fTrackSegments->Expand(size) ;
445 size = fRecParticles->GetEntries() ;
446 fRecParticles->Expand(size) ;
448 gAlice->TreeR()->Fill() ;
451 gAlice->TreeR()->Write(0,TObject::kOverwrite) ;
453 // Deleting reconstructed objects
454 ResetReconstruction();
458 //____________________________________________________________________________
459 void AliPHOSv1::ResetHits()
461 // Reset hit tree for EMC and CPV
462 // Yuri Kharlov, 28 September 2000
464 AliDetector::ResetHits();
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();
471 //____________________________________________________________________________
472 void AliPHOSv1::ResetReconstruction()
474 // Deleting reconstructed objects
476 if ( fEmcRecPoints ) fEmcRecPoints ->Delete();
477 if ( fPpsdRecPoints ) fPpsdRecPoints->Delete();
478 if ( fTrackSegments ) fTrackSegments->Delete();
479 if ( fRecParticles ) fRecParticles ->Delete();
483 //____________________________________________________________________________
484 void AliPHOSv1::SetTreeAddress()
487 AliPHOS::SetTreeAddress();
489 // //Branch address for TreeR: RecPpsdRecPoint
490 // TTree *treeR = gAlice->TreeR();
491 // if ( treeR && fPpsdRecPoints ) {
492 // branch = treeR->GetBranch("PHOSPpsdRP");
493 // if (branch) branch->SetAddress(&fPpsdRecPoints) ;
496 // Set branch address for the Hits Tree for hits in EMC modules
497 // Yuri Kharlov, 23 November 2000.
499 for( Int_t i=0; i<fGeom->GetNModules(); i++ ) GetEMCModule(i).SetTreeAddress("EMC",i+1);
501 // Set branch address for the Hits Tree for hits in CPV modules for IHEP geometry
502 // Yuri Kharlov, 28 September 2000.
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);
510 //____________________________________________________________________________
512 void AliPHOSv1::StepManager(void)
514 // Accumulates hits as long as the track stays in a single crystal or PPSD gas Cell
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;
521 Int_t relid[4] ; // (box, layer, row, column) indices
522 Int_t absid ; // absolute cell ID number
523 Float_t xyze[4]={0,0,0,0} ; // position wrt MRS and energy deposited
524 TLorentzVector pos ; // Lorentz vector of the track current position
528 Int_t tracknumber = gAlice->CurrentTrack() ;
529 Int_t primary = gAlice->GetPrimary( gAlice->CurrentTrack() );
530 TString name = fGeom->GetName() ;
531 Int_t trackpid = gMC->TrackPid() ;
532 if ( name == "GPS2" || name == "MIXT" ) { // ======> CPV is a GPS' PPSD
534 if( gMC->CurrentVolID(copy) == gMC->VolId("GCEL") ) // We are inside a gas cell
536 gMC->TrackPosition(pos) ;
540 xyze[3] = gMC->Edep() ;
542 if ( xyze[3] != 0 ) { // there is deposited energy
543 gMC->CurrentVolOffID(5, relid[0]) ; // get the PHOS Module number
544 if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(5),"PHO1") == 0 ){
545 relid[0] += fGeom->GetNModules() - fGeom->GetNPPSDModules();
547 gMC->CurrentVolOffID(3, relid[1]) ; // get the Micromegas Module number
548 // 1-> fGeom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() upper
549 // > fGeom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() lower
550 gMC->CurrentVolOffID(1, relid[2]) ; // get the row number of the cell
551 gMC->CurrentVolID(relid[3]) ; // get the column number
553 // get the absolute Id number
555 fGeom->RelToAbsNumbering(relid, absid) ;
557 // add current hit to the hit list
558 AddHit(fIshunt, primary, tracknumber, absid, xyze, trackpid);
560 } // there is deposited energy
561 } // We are inside the gas of the CPV
562 } // GPS2 configuration
564 if ( name == "IHEP" || name == "MIXT" ) { // ======> CPV is a IHEP's one
566 // Yuri Kharlov, 28 September 2000
568 if( gMC->CurrentVolID(copy) == gMC->VolId("CPVQ") &&
569 gMC->IsTrackEntering() &&
570 gMC->TrackCharge() != 0) {
572 // Charged track has just entered to the CPV sensitive plane
574 AliPHOSv1 &phos = *(AliPHOSv1*)gAlice->GetModule("PHOS");
577 gMC->CurrentVolOffID(3,moduleNumber);
580 // Current position of the hit in the CPV module ref. system
582 gMC -> TrackPosition(pos);
583 Float_t xyzm[3], xyzd[3], xyd[3]={0,0,0};
585 for (i=0; i<3; i++) xyzm[i] = pos[i];
586 gMC -> Gmtod (xyzm, xyzd, 1); // transform coordinate from master to daughter system
590 // Current momentum of the hit's track in the CPV module ref. system
593 gMC -> TrackMomentum(pmom);
594 Float_t pm[3], pd[3];
595 for (i=0; i<3; i++) pm[i] = pmom[i];
596 gMC -> Gmtod (pm, pd, 2); // transform 3-momentum from master to daughter system
601 // Current particle type of the hit's track
603 Int_t ipart = gMC->TrackPid();
605 // Add the current particle in the list of the CPV hits.
607 phos.GetCPVModule(moduleNumber).AddHit(fIshunt,primary,pmom,xyd,ipart);
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());
612 printf( " xy = (%8.4f, %8.4f) cm, ipart = %d, primary = %d\n",
613 xyd[0],xyd[1],ipart,primary);
616 // Digitize the current CPV hit:
618 // 1. find pad response and
620 TClonesArray *cpvDigits = new TClonesArray("AliPHOSCPVDigit",0); // array of digits for current hit
621 CPVDigitize(pmom,xyd,moduleNumber,cpvDigits);
626 Int_t idigit,ndigits;
628 // 2. go through the current digit list and sum digits in pads
630 ndigits = cpvDigits->GetEntriesFast();
631 for (idigit=0; idigit<ndigits-1; idigit++) {
632 AliPHOSCPVDigit *cpvDigit1 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit);
633 Float_t x1 = cpvDigit1->GetXpad() ;
634 Float_t z1 = cpvDigit1->GetYpad() ;
635 for (Int_t jdigit=idigit+1; jdigit<ndigits; jdigit++) {
636 AliPHOSCPVDigit *cpvDigit2 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(jdigit);
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) ;
646 cpvDigits->Compress() ;
648 // 3. add digits to temporary hit list fTmpHits
650 ndigits = cpvDigits->GetEntriesFast();
651 for (idigit=0; idigit<ndigits; idigit++) {
652 AliPHOSCPVDigit *cpvDigit = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit);
653 relid[0] = moduleNumber + 1 ; // CPV (or PHOS) module number
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
658 // get the absolute Id number
659 fGeom->RelToAbsNumbering(relid, absid) ;
661 // add current digit to the temporary hit list
665 xyze[3] = cpvDigit->GetQpad() ; // amplitude in a pad
666 primary = -1; // No need in primary for CPV
667 AddHit(fIshunt, primary, tracknumber, absid, xyze, trackpid);
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();
677 } // end of IHEP configuration
679 if(gMC->CurrentVolID(copy) == gMC->VolId("PXTL") ) { // We are inside a PBWO crystal
680 gMC->TrackPosition(pos) ;
684 xyze[3] = gMC->Edep() ;
686 // Track enters to the crystal from the top edge
688 if (gMC->IsTrackEntering()) {
690 gMC -> Gmtod (xyze, posloc, 1);
691 if (posloc[1] > fGeom->GetCrystalSize(1)/2-0.01) {
693 Float_t xyd[3]={0,0,0};
694 AliPHOSv1 &phos = *(AliPHOSv1*)gAlice->GetModule("PHOS");
697 gMC->CurrentVolOffID(10,moduleNumber);
699 if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(10),"PHO1") == 0 )
700 moduleNumber += fGeom->GetNModules() - fGeom->GetNPPSDModules();
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());
710 // Current momentum of the hit's track in the CPV module ref. system
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
721 // Current particle type of the hit's track
723 Int_t ipart = gMC->TrackPid();
725 // Add the current particle in the list of the EMC hits.
727 phos.GetEMCModule(moduleNumber).AddHit(fIshunt,primary,pmom,xyd,ipart);
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);
738 // Track is inside the crystal and deposits some energy
740 if ( xyze[3] != 0 ) {
741 gMC->CurrentVolOffID(10, relid[0]) ; // get the PHOS module number ;
742 if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(10),"PHO1") == 0 )
743 relid[0] += fGeom->GetNModules() - fGeom->GetNPPSDModules();
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
749 // get the absolute Id number
750 fGeom->RelToAbsNumbering(relid, absid) ;
752 // add current hit to the hit list
753 AddHit(fIshunt, primary,tracknumber, absid, xyze, trackpid);
755 } // there is deposited energy
756 } // we are inside a PHOS Xtal
759 //____________________________________________________________________________
760 void AliPHOSv1::CPVDigitize (TLorentzVector p, Float_t *zxhit, Int_t moduleNumber, TClonesArray *cpvDigits)
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
768 // Author: Yuri Kharlov (after Serguei Sadovsky)
770 // ------------------------------------------------------------------------
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
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
786 Float_t hitX = zxhit[0];
787 Float_t hitZ =-zxhit[1];
790 Float_t pNorm = p.Py();
791 Float_t eloss = kdEdx;
793 // cout << "CPVDigitize: YVK : "<<hitX<<" "<<hitZ<<" | "<<pX<<" "<<pZ<<" "<<pNorm<<endl;
795 Float_t dZY = pZ/pNorm * fGeom->GetCPVGasThickness();
796 Float_t dXY = pX/pNorm * fGeom->GetCPVGasThickness();
797 gRandom->Rannor(rnor1,rnor2);
798 eloss *= (1 + kDetR*rnor1) *
799 TMath::Sqrt((1 + ( pow(dZY,2) + pow(dXY,2) ) / pow(fGeom->GetCPVGasThickness(),2)));
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;
805 Int_t iwht1 = (Int_t) (xhit1 / kCelWr); // wire (x) coordinate "in"
806 Int_t iwht2 = (Int_t) (xhit2 / kCelWr); // wire (x) coordinate "out"
810 if (iwht1==iwht2) { // incline 1-wire hit
812 zxe[0][0] = (zhit1 + zhit2 - dZY*0.57735) / 2;
813 zxe[1][0] = (iwht1 + 0.5) * kCelWr;
815 zxe[0][1] = (zhit1 + zhit2 + dZY*0.57735) / 2;
816 zxe[1][1] = (iwht1 + 0.5) * kCelWr;
819 else if (TMath::Abs(iwht1-iwht2) != 1) { // incline 3-wire hit
821 Int_t iwht3 = (iwht1 + iwht2) / 2;
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
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;
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 );
832 zxe[0][0] = (dXY*(xwr13-xwht1)/dXY + zhit1 + zhit1) / 2;
834 zxe[2][0] = eloss * egm1;
835 zxe[0][1] = (dXY*(xwr23-xwht1)/dXY + zhit1 + zhit2) / 2;
837 zxe[2][1] = eloss * egm2;
838 zxe[0][2] = dXY*(xwht3-xwht1)/dXY + zhit1;
840 zxe[2][2] = eloss * egm3;
842 else { // incline 2-wire hit
844 Float_t xwht1 = (iwht1 + 0.5) * kCelWr;
845 Float_t xwht2 = (iwht2 + 0.5) * kCelWr;
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;
853 zxe[2][0] = eloss * egm1;
854 zxe[0][1] = (zhit1 + zhit2 + dZY*egm2) / 2;
856 zxe[2][1] = eloss * egm2;
859 // Finite size of ionization region
861 Int_t nCellZ = fGeom->GetNumberOfCPVPadsZ();
862 Int_t nCellX = fGeom->GetNumberOfCPVPadsPhi();
863 Int_t nz3 = (kNgamz+1)/2;
864 Int_t nx3 = (kNgamx+1)/2;
865 cpvDigits->Expand(nIter*kNgamx*kNgamz);
866 TClonesArray &ldigits = *(TClonesArray *)cpvDigits;
868 for (Int_t iter=0; iter<nIter; iter++) {
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;
881 for (Int_t iz=1; iz<=kNgamz; iz++) {
882 Int_t kzg = izcell + iz - nz3;
883 if (kzg<=0 || kzg>nCellZ) continue;
884 Float_t zg = (Float_t)(iz-nz3) - zc;
885 for (Int_t ix=1; ix<=kNgamx; ix++) {
886 Int_t kxg = ixcell + ix - nx3;
887 if (kxg<=0 || kxg>nCellX) continue;
888 Float_t xg = (Float_t)(ix-nx3) - xc;
890 // Now calculate pad response
891 Float_t qpad = CPVPadResponseFunction(qhit,zg,xg);
892 qpad += kNoise*rnor2;
893 if (qpad<0) continue;
895 // Fill the array with pad response ID and amplitude
896 new(ldigits[cpvDigits->GetEntriesFast()]) AliPHOSCPVDigit(kxg,kzg,qpad);
902 //____________________________________________________________________________
903 Float_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)
909 // ------------------------------------------------------------------------
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;
921 //____________________________________________________________________________
922 Double_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)
931 // ------------------------------------------------------------------------
933 const Double_t kA=1.0;
934 const Double_t kB=0.7;
936 Double_t r2 = x*x + y*y;
938 Double_t cumulPRF = 0;
939 for (Int_t i=0; i<=4; i++) {
940 Double_t b1 = (2*i + 1) * kB;
941 cumulPRF += TMath::Power(-1,i) * TMath::ATan( xy / (b1*TMath::Sqrt(b1*b1 + r2)) );
943 cumulPRF *= kA/(2*TMath::Pi());