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 fReconstructioner = 0;
66 fTrackSegmentMaker = 0;
70 //____________________________________________________________________________
71 AliPHOSv1::AliPHOSv1(const char *name, const char *title):
74 // ctor : title is used to identify the layout
75 // GPS2 = 5 modules (EMC + PPSD)
76 // IHEP = 5 modules (EMC + CPV )
77 // MIXT = 4 modules (EMC + CPV ) and 1 module (EMC + PPSD)
80 // - fHits (the "normal" one), which retains the hits associated with
81 // the current primary particle being tracked
82 // (this array is reset after each primary has been tracked).
85 fPinElectronicNoise = 0.010 ;
86 fDigitThreshold = 0.01 ; // 1 GeV
88 fDigitizeB = 10000000. ;
91 // We do not want to save in TreeH the raw hits
92 // But save the cumulated hits instead (need to create the branch myself)
93 // It is put in the Digit Tree because the TreeH is filled after each primary
94 // and the TreeD at the end of the event (branch is set in FinishEvent() ).
96 fHits= new TClonesArray("AliPHOSHit",1000) ;
100 fReconstructioner = 0;
101 fTrackSegmentMaker = 0;
103 fIshunt = 1 ; // All hits are associated with primary particles
107 //____________________________________________________________________________
108 AliPHOSv1::AliPHOSv1(AliPHOSReconstructioner * Reconstructioner, const char *name, const char *title):
109 AliPHOSv0(name,title)
111 // ctor : title is used to identify the layout
112 // GPS2 = 5 modules (EMC + PPSD)
114 fPinElectronicNoise = 0.010 ;
116 // We do not want to save in TreeH the raw hits
119 fHits= new TClonesArray("AliPHOSHit",1000) ;
123 fIshunt = 1 ; // All hits are associated with primary particles
125 // gets an instance of the geometry parameters class
126 fGeom = AliPHOSGeometry::GetInstance(title, "") ;
128 if (fGeom->IsInitialized() )
129 cout << "AliPHOS" << Version() << " : PHOS geometry intialized for " << fGeom->GetName() << endl ;
131 cout << "AliPHOS" << Version() << " : PHOS geometry initialization failed !" << endl ;
133 // Defining the PHOS Reconstructioner
135 fReconstructioner = Reconstructioner ;
139 //____________________________________________________________________________
140 AliPHOSv1::~AliPHOSv1()
162 if ( fEmcRecPoints ) {
163 fEmcRecPoints->Delete() ;
164 delete fEmcRecPoints ;
168 if ( fPpsdRecPoints ) {
169 fPpsdRecPoints->Delete() ;
170 delete fPpsdRecPoints ;
174 if ( fTrackSegments ) {
175 fTrackSegments->Delete() ;
176 delete fTrackSegments ;
182 //____________________________________________________________________________
183 void AliPHOSv1::AddHit(Int_t shunt, Int_t primary, Int_t tracknumber, Int_t Id, Float_t * hits, Int_t trackpid, TLorentzVector p, Float_t * lpos)
185 // Add a hit to the hit list.
186 // A PHOS hit is the sum of all hits in a single crystal
187 // or in a single PPSD gas cell
192 Bool_t deja = kFALSE ;
194 newHit = new AliPHOSHit(shunt, primary, tracknumber, Id, hits, trackpid, p, lpos) ;
196 for ( hitCounter = fNhits-1 ; hitCounter >= 0 && !deja ; hitCounter-- ) {
197 curHit = (AliPHOSHit*) (*fHits)[hitCounter] ;
198 if( *curHit == *newHit ) {
199 *curHit = *curHit + *newHit ;
205 new((*fHits)[fNhits]) AliPHOSHit(*newHit) ;
212 //____________________________________________________________________________
213 void AliPHOSv1::Hits2SDigits(){
214 //Collects all hits in the same active volume into digit
219 AliPHOSDigit * newdigit ;
220 AliPHOSDigit * curdigit ;
221 Bool_t deja = kFALSE ;
225 for (itrack=0; itrack<gAlice->GetNtrack(); itrack++){
227 //=========== Get the Hits Tree for the Primary track itrack
229 gAlice->TreeH()->GetEvent(itrack);
232 for ( i = 0 ; i < fHits->GetEntries() ; i++ ) {
233 hit = (AliPHOSHit*)fHits->At(i) ;
235 // Assign primary number only if contribution is significant
236 if( hit->GetEnergy() > fDigitThreshold)
237 newdigit = new AliPHOSDigit( hit->GetPrimary(), hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
239 newdigit = new AliPHOSDigit( -1 , hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
242 for ( j = 0 ; j < fnSdigits ; j++) {
243 curdigit = (AliPHOSDigit*) fSDigits->At(j) ;
244 if ( *curdigit == *newdigit) {
245 *curdigit = *curdigit + *newdigit ;
251 new((*fSDigits)[fnSdigits]) AliPHOSDigit(* newdigit) ;
258 } // loop over tracks
262 fnSdigits = fSDigits->GetEntries() ;
263 fSDigits->Expand(fnSdigits) ;
265 for (i = 0 ; i < fnSdigits ; i++) {
266 AliPHOSDigit * digit = (AliPHOSDigit *) fSDigits->At(i) ;
267 digit->SetIndexInList(i) ;
270 gAlice->TreeS()->Fill() ;
271 gAlice->TreeS()->Write(0,TObject::kOverwrite) ;
275 //____________________________________________________________________________
276 void AliPHOSv1::SDigits2Digits(){
277 //Adds noise to the summable digits and removes everething below thresholds
278 //Note, that sDigits should be SORTED in accordance with abs ID.
281 gAlice->TreeS()->GetEvent(0) ;
283 // First calculate noise induced by the PIN diode of the PbWO crystals
284 Int_t iCurSDigit = 0 ;
286 //we assume, that there is al least one EMC digit...
287 if(fSDigits->GetEntries() == 0) {
288 cout << "PHOS::SDigits2Digits> No SDigits !!! Do not produce Digits " << endl ;
292 Int_t idCurSDigit = ((AliPHOSDigit *)fSDigits->At(0))->GetId() ;
295 for(absID = 1; absID < fGeom->GetNModules()*fGeom->GetNPhi()*fGeom->GetNZ(); absID++){
296 Float_t noise = gRandom->Gaus(0., fPinElectronicNoise) ;
297 if(absID < idCurSDigit ){
298 if(noise >fDigitThreshold ){
299 new((*fDigits)[fNdigits]) AliPHOSDigit( -1,absID,Digitize(noise) ) ;
303 else{ //add noise and may be remove the true hit
304 Float_t signal = noise + Calibrate(((AliPHOSDigit *)fSDigits->At(iCurSDigit))->GetAmp()) ;
305 if( signal >fDigitThreshold ){
306 AliPHOSDigit * digit = (AliPHOSDigit*) fSDigits->At(iCurSDigit) ;
307 new((*fDigits)[fNdigits]) AliPHOSDigit( *digit ) ;
308 ((AliPHOSDigit *)fDigits->At(fNdigits))->SetAmp(Digitize(signal));
312 if(iCurSDigit < fSDigits->GetEntries()-1){
314 idCurSDigit = ((AliPHOSDigit*)fSDigits->At(iCurSDigit))->GetId() ;
317 idCurSDigit = 10000000; //no real hits left
322 //remove PPSD/CPV digits below thresholds
324 for(idigit = iCurSDigit; idigit < fSDigits->GetEntries() ; idigit++){ //loop over CPV/PPSD digits
326 AliPHOSDigit * digit = (AliPHOSDigit *) fSDigits->At(idigit) ;
327 Float_t ene = Calibrate(digit->GetAmp()) ;
330 fGeom->AbsToRelNumbering(digit->GetId(), relid) ;
331 if ( relid[0] > fGeom->GetNCPVModules() ){ //ppsd
332 if ( ( (relid[1] > 0) && (ene > fPpsdEnergyThreshold)) || //PPSD digit
333 ( (relid[1] < 0) && (ene > fCpvEnergyThreshold ) ) ) //CPV digit
334 new((*fDigits)[fNdigits]) AliPHOSDigit( *digit ) ;
339 fDigits->Compress() ;
341 fNdigits = fDigits->GetEntries() ;
342 fDigits->Expand(fNdigits) ;
345 for (i = 0 ; i < fNdigits ; i++) {
346 AliPHOSDigit * digit = (AliPHOSDigit *) fDigits->At(i) ;
347 digit->SetIndexInList(i) ;
350 gAlice->TreeD()->Fill() ;
352 gAlice->TreeD()->Write(0,TObject::kOverwrite) ;
356 //___________________________________________________________________________
357 void AliPHOSv1::MakeBranch(Option_t* opt, char *file)
362 // Create new branche in the current Root Tree in the digit Tree
363 AliDetector::MakeBranch(opt) ;
366 cH = strstr(opt,"S");
367 //Create a branch for SDigits
370 sprintf(branchname,"%s",GetName());
374 fSDigits = new TClonesArray("AliPHOSDigit",1000);
377 gAlice->MakeBranchInTree(gAlice->TreeS(),branchname,&fSDigits,fBufferSize,file);
380 cH = strstr(opt,"D");
381 //Create a branch for Digits
384 sprintf(branchname,"%s",GetName());
389 fDigits = new TClonesArray("AliPHOSDigit",1000);
392 gAlice->MakeBranchInTree(gAlice->TreeD(),branchname,&fDigits,fBufferSize,file);
395 cH = strstr(opt,"R");
396 //Create a branch for Reconstruction
400 Int_t splitlevel = 0 ;
403 fEmcRecPoints->Delete() ;
405 fEmcRecPoints = new TObjArray(100) ;
407 if ( fEmcRecPoints && gAlice->TreeR() ) {
408 sprintf(branchname,"%sEmcRP",GetName()) ;
409 gAlice->TreeR()->Branch(branchname, "TObjArray", &fEmcRecPoints, fBufferSize, splitlevel) ;
413 fPpsdRecPoints->Delete() ;
415 fPpsdRecPoints = new TObjArray(100) ;
417 if ( fPpsdRecPoints && gAlice->TreeR() ) {
418 sprintf(branchname,"%sPpsdRP",GetName()) ;
419 gAlice->TreeR()->Branch(branchname, "TObjArray", &fPpsdRecPoints, fBufferSize, splitlevel) ;
423 fTrackSegments->Delete() ;
425 fTrackSegments = new TClonesArray("AliPHOSTrackSegment",100) ;
427 if ( fTrackSegments && gAlice->TreeR() ) {
428 sprintf(branchname,"%sTS",GetName()) ;
429 gAlice->TreeR()->Branch(branchname, &fTrackSegments, fBufferSize) ;
433 fRecParticles->Delete() ;
435 fRecParticles = new TClonesArray("AliPHOSRecParticle",100) ;
437 if ( fRecParticles && gAlice->TreeR() ) {
438 sprintf(branchname,"%sRP",GetName()) ;
439 gAlice->TreeR()->Branch(branchname, &fRecParticles, fBufferSize) ;
446 //_____________________________________________________________________________
447 void AliPHOSv1::Reconstruction(AliPHOSReconstructioner * Reconstructioner)
449 // 1. Reinitializes the existing RecPoint, TrackSegment, and RecParticles Lists and
450 // 2. Creates TreeR with a branch for each list
451 // 3. Steers the reconstruction processes
452 // 4. Saves the 3 lists in TreeR
453 // 5. Write the Tree to File
455 fReconstructioner = Reconstructioner ;
459 // gAlice->MakeTree("R") ;
465 fReconstructioner->Make(fDigits, fEmcRecPoints, fPpsdRecPoints, fTrackSegments, fRecParticles);
467 printf("Reconstruction: %d %d %d %d\n",
468 fEmcRecPoints->GetEntries(),fPpsdRecPoints->GetEntries(),
469 fTrackSegments->GetEntries(),fRecParticles->GetEntries());
471 // 4. Expand or Shrink the arrays to the proper size
475 size = fEmcRecPoints->GetEntries() ;
476 fEmcRecPoints->Expand(size) ;
478 size = fPpsdRecPoints->GetEntries() ;
479 fPpsdRecPoints->Expand(size) ;
481 size = fTrackSegments->GetEntries() ;
482 fTrackSegments->Expand(size) ;
484 size = fRecParticles->GetEntries() ;
485 fRecParticles->Expand(size) ;
487 gAlice->TreeR()->Fill() ;
490 gAlice->TreeR()->Write(0,TObject::kOverwrite) ;
492 // Deleting reconstructed objects
493 ResetReconstruction();
497 //____________________________________________________________________________
498 void AliPHOSv1::ResetReconstruction()
500 // Deleting reconstructed objects
502 if ( fEmcRecPoints ) fEmcRecPoints ->Delete();
503 if ( fPpsdRecPoints ) fPpsdRecPoints->Delete();
504 if ( fTrackSegments ) fTrackSegments->Delete();
505 if ( fRecParticles ) fRecParticles ->Delete();
509 //____________________________________________________________________________
511 void AliPHOSv1::StepManager(void)
513 // Accumulates hits as long as the track stays in a single crystal or PPSD gas Cell
515 Int_t relid[4] ; // (box, layer, row, column) indices
516 Int_t absid ; // absolute cell ID number
517 Float_t xyze[4]={0,0,0,0} ; // position wrt MRS and energy deposited
518 TLorentzVector pos ; // Lorentz vector of the track current position
519 TLorentzVector pmom ; //momentum of the particle initiated hit
520 Float_t xyd[3]={0,0,0} ; //local posiiton of the entering
521 Bool_t entered = kFALSE ;
524 Int_t tracknumber = gAlice->CurrentTrack() ;
525 Int_t primary = gAlice->GetPrimary( gAlice->CurrentTrack() );
526 TString name = fGeom->GetName() ;
529 if( gMC->IsTrackEntering() ){ // create hit with position and momentum of new particle,
530 // but may be without energy deposition
532 // Current position of the hit in the local ref. system
533 gMC -> TrackPosition(pos);
534 Float_t xyzm[3], xyzd[3] ;
536 for (i=0; i<3; i++) xyzm[i] = pos[i];
537 gMC -> Gmtod (xyzm, xyzd, 1); // transform coordinate from master to daughter system
543 // Current momentum of the hit's track in the local ref. system
544 gMC -> TrackMomentum(pmom);
545 Float_t pm[3], pd[3];
546 for (i=0; i<3; i++) pm[i] = pmom[i];
547 gMC -> Gmtod (pm, pd, 2); // transform 3-momentum from master to daughter system
552 trackpid = gMC->TrackPid();
553 entered = kTRUE ; // Mark to create hit even withou energy deposition
558 if ( name == "GPS2" || name == "MIXT" ) { // ======> CPV is a GPS' PPSD
560 if( gMC->CurrentVolID(copy) == gMC->VolId("GCEL") ) // We are inside a gas cell
562 gMC->TrackPosition(pos) ;
566 xyze[3] = gMC->Edep() ;
568 if ( (xyze[3] != 0) || entered ) { // there is deposited energy or new particle entering PPSD
569 gMC->CurrentVolOffID(5, relid[0]) ; // get the PHOS Module number
570 if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(5),"PHO1") == 0 ){
571 relid[0] += fGeom->GetNModules() - fGeom->GetNPPSDModules();
573 gMC->CurrentVolOffID(3, relid[1]) ; // get the Micromegas Module number
574 // 1-> fGeom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() upper
575 // > fGeom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() lower
576 gMC->CurrentVolOffID(1, relid[2]) ; // get the row number of the cell
577 gMC->CurrentVolID(relid[3]) ; // get the column number
579 // get the absolute Id number
581 fGeom->RelToAbsNumbering(relid, absid) ;
583 // add current hit to the hit list
584 AddHit(fIshunt, primary, tracknumber, absid, xyze, trackpid, pmom, xyd);
587 } // there is deposited energy
588 } // We are inside the gas of the CPV
589 } // GPS2 configuration
591 if ( name == "IHEP" || name == "MIXT" ) { // ======> CPV is a IHEP's one
593 // Yuri Kharlov, 28 September 2000
595 if( gMC->CurrentVolID(copy) == gMC->VolId("CPVQ") &&
597 gMC->TrackCharge() != 0) {
599 // Digitize the current CPV hit:
601 // 1. find pad response and
604 gMC->CurrentVolOffID(3,moduleNumber);
608 TClonesArray *cpvDigits = new TClonesArray("AliPHOSCPVDigit",0); // array of digits for current hit
609 CPVDigitize(pmom,xyd,moduleNumber,cpvDigits);
614 Int_t idigit,ndigits;
616 // 2. go through the current digit list and sum digits in pads
618 ndigits = cpvDigits->GetEntriesFast();
619 for (idigit=0; idigit<ndigits-1; idigit++) {
620 AliPHOSCPVDigit *cpvDigit1 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit);
621 Float_t x1 = cpvDigit1->GetXpad() ;
622 Float_t z1 = cpvDigit1->GetYpad() ;
623 for (Int_t jdigit=idigit+1; jdigit<ndigits; jdigit++) {
624 AliPHOSCPVDigit *cpvDigit2 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(jdigit);
625 Float_t x2 = cpvDigit2->GetXpad() ;
626 Float_t z2 = cpvDigit2->GetYpad() ;
627 if (x1==x2 && z1==z2) {
628 Float_t qsum = cpvDigit1->GetQpad() + cpvDigit2->GetQpad() ;
629 cpvDigit2->SetQpad(qsum) ;
630 cpvDigits->RemoveAt(idigit) ;
634 cpvDigits->Compress() ;
636 // 3. add digits to temporary hit list fTmpHits
638 ndigits = cpvDigits->GetEntriesFast();
639 for (idigit=0; idigit<ndigits; idigit++) {
640 AliPHOSCPVDigit *cpvDigit = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit);
641 relid[0] = moduleNumber + 1 ; // CPV (or PHOS) module number
642 relid[1] =-1 ; // means CPV
643 relid[2] = cpvDigit->GetXpad() ; // column number of a pad
644 relid[3] = cpvDigit->GetYpad() ; // row number of a pad
646 // get the absolute Id number
647 fGeom->RelToAbsNumbering(relid, absid) ;
649 // add current digit to the temporary hit list
653 xyze[3] = cpvDigit->GetQpad() ; // amplitude in a pad
654 primary = -1; // No need in primary for CPV
655 AddHit(fIshunt, primary, tracknumber, absid, xyze, trackpid, pmom, xyd);
657 if (cpvDigit->GetQpad() > 0.02) {
658 xmean += cpvDigit->GetQpad() * (cpvDigit->GetXpad() + 0.5);
659 zmean += cpvDigit->GetQpad() * (cpvDigit->GetYpad() + 0.5);
660 qsum += cpvDigit->GetQpad();
665 } // end of IHEP configuration
668 if(gMC->CurrentVolID(copy) == gMC->VolId("PXTL") ) { // We are inside a PBWO crystal
669 gMC->TrackPosition(pos) ;
673 xyze[3] = gMC->Edep() ;
676 if ( (xyze[3] != 0) || entered ) { // Track is inside the crystal and deposits some energy or just entered
678 gMC->CurrentVolOffID(10, relid[0]) ; // get the PHOS module number ;
680 if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(10),"PHO1") == 0 )
681 relid[0] += fGeom->GetNModules() - fGeom->GetNPPSDModules();
683 relid[1] = 0 ; // means PBW04
684 gMC->CurrentVolOffID(4, relid[2]) ; // get the row number inside the module
685 gMC->CurrentVolOffID(3, relid[3]) ; // get the cell number inside the module
687 // get the absolute Id number
688 fGeom->RelToAbsNumbering(relid, absid) ;
690 // add current hit to the hit list
691 AddHit(fIshunt, primary,tracknumber, absid, xyze, trackpid,pmom, xyd);
694 } // there is deposited energy
695 } // we are inside a PHOS Xtal
700 //____________________________________________________________________________
701 void AliPHOSv1::CPVDigitize (TLorentzVector p, Float_t *zxhit, Int_t moduleNumber, TClonesArray *cpvDigits)
703 // ------------------------------------------------------------------------
704 // Digitize one CPV hit:
705 // On input take exact 4-momentum p and position zxhit of the hit,
706 // find the pad response around this hit and
707 // put the amplitudes in the pads into array digits
709 // Author: Yuri Kharlov (after Serguei Sadovsky)
711 // ------------------------------------------------------------------------
713 const Float_t kCelWr = fGeom->GetPadSizePhi()/2; // Distance between wires (2 wires above 1 pad)
714 const Float_t kDetR = 0.1; // Relative energy fluctuation in track for 100 e-
715 const Float_t kdEdx = 4.0; // Average energy loss in CPV;
716 const Int_t kNgamz = 5; // Ionization size in Z
717 const Int_t kNgamx = 9; // Ionization size in Phi
718 const Float_t kNoise = 0.03; // charge noise in one pad
722 // Just a reminder on axes notation in the CPV module:
723 // axis Z goes along the beam
724 // axis X goes across the beam in the module plane
725 // axis Y is a normal to the module plane showing from the IP
727 Float_t hitX = zxhit[0];
728 Float_t hitZ =-zxhit[1];
731 Float_t pNorm = p.Py();
732 Float_t eloss = kdEdx;
734 // cout << "CPVDigitize: YVK : "<<hitX<<" "<<hitZ<<" | "<<pX<<" "<<pZ<<" "<<pNorm<<endl;
736 Float_t dZY = pZ/pNorm * fGeom->GetCPVGasThickness();
737 Float_t dXY = pX/pNorm * fGeom->GetCPVGasThickness();
738 gRandom->Rannor(rnor1,rnor2);
739 eloss *= (1 + kDetR*rnor1) *
740 TMath::Sqrt((1 + ( pow(dZY,2) + pow(dXY,2) ) / pow(fGeom->GetCPVGasThickness(),2)));
741 Float_t zhit1 = hitZ + fGeom->GetCPVActiveSize(1)/2 - dZY/2;
742 Float_t xhit1 = hitX + fGeom->GetCPVActiveSize(0)/2 - dXY/2;
743 Float_t zhit2 = zhit1 + dZY;
744 Float_t xhit2 = xhit1 + dXY;
746 Int_t iwht1 = (Int_t) (xhit1 / kCelWr); // wire (x) coordinate "in"
747 Int_t iwht2 = (Int_t) (xhit2 / kCelWr); // wire (x) coordinate "out"
751 if (iwht1==iwht2) { // incline 1-wire hit
753 zxe[0][0] = (zhit1 + zhit2 - dZY*0.57735) / 2;
754 zxe[1][0] = (iwht1 + 0.5) * kCelWr;
756 zxe[0][1] = (zhit1 + zhit2 + dZY*0.57735) / 2;
757 zxe[1][1] = (iwht1 + 0.5) * kCelWr;
760 else if (TMath::Abs(iwht1-iwht2) != 1) { // incline 3-wire hit
762 Int_t iwht3 = (iwht1 + iwht2) / 2;
763 Float_t xwht1 = (iwht1 + 0.5) * kCelWr; // wire 1
764 Float_t xwht2 = (iwht2 + 0.5) * kCelWr; // wire 2
765 Float_t xwht3 = (iwht3 + 0.5) * kCelWr; // wire 3
766 Float_t xwr13 = (xwht1 + xwht3) / 2; // center 13
767 Float_t xwr23 = (xwht2 + xwht3) / 2; // center 23
768 Float_t dxw1 = xhit1 - xwr13;
769 Float_t dxw2 = xhit2 - xwr23;
770 Float_t egm1 = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
771 Float_t egm2 = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
772 Float_t egm3 = kCelWr / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
773 zxe[0][0] = (dXY*(xwr13-xwht1)/dXY + zhit1 + zhit1) / 2;
775 zxe[2][0] = eloss * egm1;
776 zxe[0][1] = (dXY*(xwr23-xwht1)/dXY + zhit1 + zhit2) / 2;
778 zxe[2][1] = eloss * egm2;
779 zxe[0][2] = dXY*(xwht3-xwht1)/dXY + zhit1;
781 zxe[2][2] = eloss * egm3;
783 else { // incline 2-wire hit
785 Float_t xwht1 = (iwht1 + 0.5) * kCelWr;
786 Float_t xwht2 = (iwht2 + 0.5) * kCelWr;
787 Float_t xwr12 = (xwht1 + xwht2) / 2;
788 Float_t dxw1 = xhit1 - xwr12;
789 Float_t dxw2 = xhit2 - xwr12;
790 Float_t egm1 = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) );
791 Float_t egm2 = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) );
792 zxe[0][0] = (zhit1 + zhit2 - dZY*egm1) / 2;
794 zxe[2][0] = eloss * egm1;
795 zxe[0][1] = (zhit1 + zhit2 + dZY*egm2) / 2;
797 zxe[2][1] = eloss * egm2;
800 // Finite size of ionization region
802 Int_t nCellZ = fGeom->GetNumberOfCPVPadsZ();
803 Int_t nCellX = fGeom->GetNumberOfCPVPadsPhi();
804 Int_t nz3 = (kNgamz+1)/2;
805 Int_t nx3 = (kNgamx+1)/2;
806 cpvDigits->Expand(nIter*kNgamx*kNgamz);
807 TClonesArray &ldigits = *(TClonesArray *)cpvDigits;
809 for (Int_t iter=0; iter<nIter; iter++) {
811 Float_t zhit = zxe[0][iter];
812 Float_t xhit = zxe[1][iter];
813 Float_t qhit = zxe[2][iter];
814 Float_t zcell = zhit / fGeom->GetPadSizeZ();
815 Float_t xcell = xhit / fGeom->GetPadSizePhi();
816 if ( zcell<=0 || xcell<=0 ||
817 zcell>=nCellZ || xcell>=nCellX) return;
818 Int_t izcell = (Int_t) zcell;
819 Int_t ixcell = (Int_t) xcell;
820 Float_t zc = zcell - izcell - 0.5;
821 Float_t xc = xcell - ixcell - 0.5;
822 for (Int_t iz=1; iz<=kNgamz; iz++) {
823 Int_t kzg = izcell + iz - nz3;
824 if (kzg<=0 || kzg>nCellZ) continue;
825 Float_t zg = (Float_t)(iz-nz3) - zc;
826 for (Int_t ix=1; ix<=kNgamx; ix++) {
827 Int_t kxg = ixcell + ix - nx3;
828 if (kxg<=0 || kxg>nCellX) continue;
829 Float_t xg = (Float_t)(ix-nx3) - xc;
831 // Now calculate pad response
832 Float_t qpad = CPVPadResponseFunction(qhit,zg,xg);
833 qpad += kNoise*rnor2;
834 if (qpad<0) continue;
836 // Fill the array with pad response ID and amplitude
837 new(ldigits[cpvDigits->GetEntriesFast()]) AliPHOSCPVDigit(kxg,kzg,qpad);
843 //____________________________________________________________________________
844 Float_t AliPHOSv1::CPVPadResponseFunction(Float_t qhit, Float_t zhit, Float_t xhit) {
845 // ------------------------------------------------------------------------
846 // Calculate the amplitude in one CPV pad using the
847 // cumulative pad response function
848 // Author: Yuri Kharlov (after Serguei Sadovski)
850 // ------------------------------------------------------------------------
852 Double_t dz = fGeom->GetPadSizeZ() / 2;
853 Double_t dx = fGeom->GetPadSizePhi() / 2;
854 Double_t z = zhit * fGeom->GetPadSizeZ();
855 Double_t x = xhit * fGeom->GetPadSizePhi();
856 Double_t amplitude = qhit *
857 (CPVCumulPadResponse(z+dz,x+dx) - CPVCumulPadResponse(z+dz,x-dx) -
858 CPVCumulPadResponse(z-dz,x+dx) + CPVCumulPadResponse(z-dz,x-dx));
859 return (Float_t)amplitude;
862 //____________________________________________________________________________
863 Double_t AliPHOSv1::CPVCumulPadResponse(Double_t x, Double_t y) {
864 // ------------------------------------------------------------------------
865 // Cumulative pad response function
866 // It includes several terms from the CF decomposition in electrostatics
867 // Note: this cumulative function is wrong since omits some terms
868 // but the cell amplitude obtained with it is correct because
869 // these omitting terms cancel
870 // Author: Yuri Kharlov (after Serguei Sadovski)
872 // ------------------------------------------------------------------------
874 const Double_t kA=1.0;
875 const Double_t kB=0.7;
877 Double_t r2 = x*x + y*y;
879 Double_t cumulPRF = 0;
880 for (Int_t i=0; i<=4; i++) {
881 Double_t b1 = (2*i + 1) * kB;
882 cumulPRF += TMath::Power(-1,i) * TMath::ATan( xy / (b1*TMath::Sqrt(b1*b1 + r2)) );
884 cumulPRF *= kA/(2*TMath::Pi());