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():
66 fReconstructioner = 0;
67 fTrackSegmentMaker = 0;
71 //____________________________________________________________________________
72 AliPHOSv1::AliPHOSv1(const char *name, const char *title):
75 // ctor : title is used to identify the layout
76 // GPS2 = 5 modules (EMC + PPSD)
77 // IHEP = 5 modules (EMC + CPV )
78 // MIXT = 4 modules (EMC + CPV ) and 1 module (EMC + PPSD)
81 // - fHits (the "normal" one), which retains the hits associated with
82 // the current primary particle being tracked
83 // (this array is reset after each primary has been tracked).
86 fPinElectronicNoise = 0.010 ;
87 fDigitThreshold = 0.01 ; // 1 GeV
89 fDigitizeB = 10000000. ;
92 // We do not want to save in TreeH the raw hits
93 // But save the cumulated hits instead (need to create the branch myself)
94 // It is put in the Digit Tree because the TreeH is filled after each primary
95 // and the TreeD at the end of the event (branch is set in FinishEvent() ).
97 fHits= new TClonesArray("AliPHOSHit",1000) ;
101 fReconstructioner = 0;
102 fTrackSegmentMaker = 0;
104 fIshunt = 1 ; // All hits are associated with primary particles
108 //____________________________________________________________________________
109 AliPHOSv1::AliPHOSv1(AliPHOSReconstructioner * Reconstructioner, const char *name, const char *title):
110 AliPHOSv0(name,title)
112 // ctor : title is used to identify the layout
113 // GPS2 = 5 modules (EMC + PPSD)
115 fPinElectronicNoise = 0.010 ;
117 // We do not want to save in TreeH the raw hits
120 fHits= new TClonesArray("AliPHOSHit",1000) ;
124 fIshunt = 1 ; // All hits are associated with primary particles
126 // gets an instance of the geometry parameters class
127 fGeom = AliPHOSGeometry::GetInstance(title, "") ;
129 if (fGeom->IsInitialized() )
130 cout << "AliPHOS" << Version() << " : PHOS geometry intialized for " << fGeom->GetName() << endl ;
132 cout << "AliPHOS" << Version() << " : PHOS geometry initialization failed !" << endl ;
134 // Defining the PHOS Reconstructioner
136 fReconstructioner = Reconstructioner ;
140 //____________________________________________________________________________
141 AliPHOSv1::~AliPHOSv1()
163 if ( fEmcRecPoints ) {
164 fEmcRecPoints->Delete() ;
165 delete fEmcRecPoints ;
169 if ( fPpsdRecPoints ) {
170 fPpsdRecPoints->Delete() ;
171 delete fPpsdRecPoints ;
175 if ( fTrackSegments ) {
176 fTrackSegments->Delete() ;
177 delete fTrackSegments ;
183 //____________________________________________________________________________
184 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)
186 // Add a hit to the hit list.
187 // A PHOS hit is the sum of all hits in a single crystal
188 // or in a single PPSD gas cell
193 Bool_t deja = kFALSE ;
195 newHit = new AliPHOSHit(shunt, primary, tracknumber, Id, hits, trackpid, p, lpos) ;
197 for ( hitCounter = fNhits-1 ; hitCounter >= 0 && !deja ; hitCounter-- ) {
198 curHit = (AliPHOSHit*) (*fHits)[hitCounter] ;
199 if( *curHit == *newHit ) {
200 *curHit = *curHit + *newHit ;
206 new((*fHits)[fNhits]) AliPHOSHit(*newHit) ;
213 //____________________________________________________________________________
214 void AliPHOSv1::Hits2SDigits(){
215 //Collects all hits in the same active volume into digit
220 AliPHOSDigit * newdigit ;
221 AliPHOSDigit * curdigit ;
222 Bool_t deja = kFALSE ;
226 for (itrack=0; itrack<gAlice->GetNtrack(); itrack++){
228 //=========== Get the Hits Tree for the Primary track itrack
230 gAlice->TreeH()->GetEvent(itrack);
233 for ( i = 0 ; i < fHits->GetEntries() ; i++ ) {
234 hit = (AliPHOSHit*)fHits->At(i) ;
236 // Assign primary number only if contribution is significant
237 if( hit->GetEnergy() > fDigitThreshold)
238 newdigit = new AliPHOSDigit( hit->GetPrimary(), hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
240 newdigit = new AliPHOSDigit( -1 , hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
243 for ( j = 0 ; j < fnSdigits ; j++) {
244 curdigit = (AliPHOSDigit*) fSDigits->At(j) ;
245 if ( *curdigit == *newdigit) {
246 *curdigit = *curdigit + *newdigit ;
252 new((*fSDigits)[fnSdigits]) AliPHOSDigit(* newdigit) ;
259 } // loop over tracks
263 fnSdigits = fSDigits->GetEntries() ;
264 fSDigits->Expand(fnSdigits) ;
266 for (i = 0 ; i < fnSdigits ; i++) {
267 AliPHOSDigit * digit = (AliPHOSDigit *) fSDigits->At(i) ;
268 digit->SetIndexInList(i) ;
271 gAlice->TreeS()->Fill() ;
272 gAlice->TreeS()->Write(0,TObject::kOverwrite) ;
276 //____________________________________________________________________________
277 void AliPHOSv1::SDigits2Digits(){
278 //Adds noise to the summable digits and removes everething below thresholds
279 //Note, that sDigits should be SORTED in accordance with abs ID.
282 gAlice->TreeS()->GetEvent(0) ;
284 // First calculate noise induced by the PIN diode of the PbWO crystals
285 Int_t iCurSDigit = 0 ;
287 //we assume, that there is al least one EMC digit...
288 if(fSDigits->GetEntries() == 0) {
289 cout << "PHOS::SDigits2Digits> No SDigits !!! Do not produce Digits " << endl ;
293 Int_t idCurSDigit = ((AliPHOSDigit *)fSDigits->At(0))->GetId() ;
296 for(absID = 1; absID < fGeom->GetNModules()*fGeom->GetNPhi()*fGeom->GetNZ(); absID++){
297 Float_t noise = gRandom->Gaus(0., fPinElectronicNoise) ;
298 if(absID < idCurSDigit ){
299 if(noise >fDigitThreshold ){
300 new((*fDigits)[fNdigits]) AliPHOSDigit( -1,absID,Digitize(noise) ) ;
304 else{ //add noise and may be remove the true hit
305 Float_t signal = noise + Calibrate(((AliPHOSDigit *)fSDigits->At(iCurSDigit))->GetAmp()) ;
306 if( signal >fDigitThreshold ){
307 AliPHOSDigit * digit = (AliPHOSDigit*) fSDigits->At(iCurSDigit) ;
308 new((*fDigits)[fNdigits]) AliPHOSDigit( *digit ) ;
309 ((AliPHOSDigit *)fDigits->At(fNdigits))->SetAmp(Digitize(signal));
313 if(iCurSDigit < fSDigits->GetEntries()-1){
315 idCurSDigit = ((AliPHOSDigit*)fSDigits->At(iCurSDigit))->GetId() ;
318 idCurSDigit = 10000000; //no real hits left
323 //remove PPSD/CPV digits below thresholds
325 for(idigit = iCurSDigit; idigit < fSDigits->GetEntries() ; idigit++){ //loop over CPV/PPSD digits
327 AliPHOSDigit * digit = (AliPHOSDigit *) fSDigits->At(idigit) ;
328 Float_t ene = Calibrate(digit->GetAmp()) ;
331 fGeom->AbsToRelNumbering(digit->GetId(), relid) ;
332 if ( relid[0] > fGeom->GetNCPVModules() ){ //ppsd
333 if ( ( (relid[1] > 0) && (ene > fPpsdEnergyThreshold)) || //PPSD digit
334 ( (relid[1] < 0) && (ene > fCpvEnergyThreshold ) ) ) //CPV digit
335 new((*fDigits)[fNdigits]) AliPHOSDigit( *digit ) ;
340 fDigits->Compress() ;
342 fNdigits = fDigits->GetEntries() ;
343 fDigits->Expand(fNdigits) ;
346 for (i = 0 ; i < fNdigits ; i++) {
347 AliPHOSDigit * digit = (AliPHOSDigit *) fDigits->At(i) ;
348 digit->SetIndexInList(i) ;
351 gAlice->TreeD()->Fill() ;
353 gAlice->TreeD()->Write(0,TObject::kOverwrite) ;
357 //___________________________________________________________________________
358 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());
375 gAlice->MakeBranchInTree(gAlice->TreeS(),branchname,&fSDigits,fBufferSize,file);
378 cH = strstr(opt,"D");
379 //Create a branch for Digits
382 sprintf(branchname,"%s",GetName());
387 gAlice->MakeBranchInTree(gAlice->TreeD(),branchname,&fDigits,fBufferSize,file);
390 cH = strstr(opt,"R");
391 //Create a branch for Reconstruction
395 Int_t splitlevel = 0 ;
398 fEmcRecPoints->Delete() ;
400 if ( fEmcRecPoints && gAlice->TreeR() ) {
401 sprintf(branchname,"%sEmcRP",GetName()) ;
402 gAlice->MakeBranchInTree(gAlice->TreeR(),branchname,"TObjArray",&fEmcRecPoints, fBufferSize, splitlevel,file);
406 fPpsdRecPoints->Delete() ;
408 if ( fPpsdRecPoints && gAlice->TreeR() ) {
409 sprintf(branchname,"%sPpsdRP",GetName()) ;
410 gAlice->MakeBranchInTree(gAlice->TreeR(),branchname,"TObjArray",&fPpsdRecPoints, fBufferSize, splitlevel,file);
414 fTrackSegments->Clear() ;
416 if ( fTrackSegments && gAlice->TreeR() ) {
417 sprintf(branchname,"%sTS",GetName()) ;
418 gAlice->MakeBranchInTree(gAlice->TreeR(),branchname,&fTrackSegments,fBufferSize,file);
422 fRecParticles->Clear() ;
424 if ( fRecParticles && gAlice->TreeR() ) {
425 sprintf(branchname,"%sRP",GetName()) ;
426 gAlice->MakeBranchInTree(gAlice->TreeR(),branchname,&fRecParticles,fBufferSize,file);
433 //_____________________________________________________________________________
434 void AliPHOSv1::Reconstruction(AliPHOSReconstructioner * Reconstructioner)
436 // 1. Reinitializes the existing RecPoint, TrackSegment, and RecParticles Lists and
437 // 2. Creates TreeR with a branch for each list
438 // 3. Steers the reconstruction processes
439 // 4. Saves the 3 lists in TreeR
440 // 5. Write the Tree to File
442 fReconstructioner = Reconstructioner ;
446 // gAlice->MakeTree("R") ;
452 fReconstructioner->Make(fDigits, fEmcRecPoints, fPpsdRecPoints, fTrackSegments, fRecParticles);
454 printf("Reconstruction: %d %d %d %d\n",
455 fEmcRecPoints->GetEntries(),fPpsdRecPoints->GetEntries(),
456 fTrackSegments->GetEntries(),fRecParticles->GetEntries());
458 // 4. Expand or Shrink the arrays to the proper size
462 size = fEmcRecPoints->GetEntries() ;
463 fEmcRecPoints->Expand(size) ;
465 size = fPpsdRecPoints->GetEntries() ;
466 fPpsdRecPoints->Expand(size) ;
468 size = fTrackSegments->GetEntries() ;
469 fTrackSegments->Expand(size) ;
471 size = fRecParticles->GetEntries() ;
472 fRecParticles->Expand(size) ;
474 gAlice->TreeR()->Fill() ;
477 gAlice->TreeR()->Write(0,TObject::kOverwrite) ;
479 // Deleting reconstructed objects
480 ResetReconstruction();
484 //____________________________________________________________________________
485 void AliPHOSv1::ResetReconstruction()
487 // Deleting reconstructed objects
489 if ( fEmcRecPoints ) fEmcRecPoints ->Delete();
490 if ( fPpsdRecPoints ) fPpsdRecPoints->Delete();
491 if ( fTrackSegments ) fTrackSegments->Delete();
492 if ( fRecParticles ) fRecParticles ->Delete();
496 //____________________________________________________________________________
498 void AliPHOSv1::StepManager(void)
500 // Accumulates hits as long as the track stays in a single crystal or PPSD gas Cell
502 Int_t relid[4] ; // (box, layer, row, column) indices
503 Int_t absid ; // absolute cell ID number
504 Float_t xyze[4]={0,0,0,0} ; // position wrt MRS and energy deposited
505 TLorentzVector pos ; // Lorentz vector of the track current position
506 TLorentzVector pmom ; //momentum of the particle initiated hit
507 Float_t xyd[3]={0,0,0} ; //local posiiton of the entering
508 Bool_t entered = kFALSE ;
511 Int_t tracknumber = gAlice->CurrentTrack() ;
512 Int_t primary = gAlice->GetPrimary( gAlice->CurrentTrack() );
513 TString name = fGeom->GetName() ;
516 if( gMC->IsTrackEntering() ){ // create hit with position and momentum of new particle,
517 // but may be without energy deposition
519 // Current position of the hit in the local ref. system
520 gMC -> TrackPosition(pos);
521 Float_t xyzm[3], xyzd[3] ;
523 for (i=0; i<3; i++) xyzm[i] = pos[i];
524 gMC -> Gmtod (xyzm, xyzd, 1); // transform coordinate from master to daughter system
530 // Current momentum of the hit's track in the local ref. system
531 gMC -> TrackMomentum(pmom);
532 Float_t pm[3], pd[3];
533 for (i=0; i<3; i++) pm[i] = pmom[i];
534 gMC -> Gmtod (pm, pd, 2); // transform 3-momentum from master to daughter system
539 trackpid = gMC->TrackPid();
540 entered = kTRUE ; // Mark to create hit even withou energy deposition
545 if ( name == "GPS2" || name == "MIXT" ) { // ======> CPV is a GPS' PPSD
547 if( gMC->CurrentVolID(copy) == gMC->VolId("GCEL") ) // We are inside a gas cell
549 gMC->TrackPosition(pos) ;
553 xyze[3] = gMC->Edep() ;
555 if ( (xyze[3] != 0) || entered ) { // there is deposited energy or new particle entering PPSD
556 gMC->CurrentVolOffID(5, relid[0]) ; // get the PHOS Module number
557 if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(5),"PHO1") == 0 ){
558 relid[0] += fGeom->GetNModules() - fGeom->GetNPPSDModules();
560 gMC->CurrentVolOffID(3, relid[1]) ; // get the Micromegas Module number
561 // 1-> fGeom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() upper
562 // > fGeom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() lower
563 gMC->CurrentVolOffID(1, relid[2]) ; // get the row number of the cell
564 gMC->CurrentVolID(relid[3]) ; // get the column number
566 // get the absolute Id number
568 fGeom->RelToAbsNumbering(relid, absid) ;
570 // add current hit to the hit list
571 AddHit(fIshunt, primary, tracknumber, absid, xyze, trackpid, pmom, xyd);
574 } // there is deposited energy
575 } // We are inside the gas of the CPV
576 } // GPS2 configuration
578 if ( name == "IHEP" || name == "MIXT" ) { // ======> CPV is a IHEP's one
580 // Yuri Kharlov, 28 September 2000
582 if( gMC->CurrentVolID(copy) == gMC->VolId("CPVQ") &&
584 gMC->TrackCharge() != 0) {
586 // Digitize the current CPV hit:
588 // 1. find pad response and
591 gMC->CurrentVolOffID(3,moduleNumber);
595 TClonesArray *cpvDigits = new TClonesArray("AliPHOSCPVDigit",0); // array of digits for current hit
596 CPVDigitize(pmom,xyd,moduleNumber,cpvDigits);
601 Int_t idigit,ndigits;
603 // 2. go through the current digit list and sum digits in pads
605 ndigits = cpvDigits->GetEntriesFast();
606 for (idigit=0; idigit<ndigits-1; idigit++) {
607 AliPHOSCPVDigit *cpvDigit1 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit);
608 Float_t x1 = cpvDigit1->GetXpad() ;
609 Float_t z1 = cpvDigit1->GetYpad() ;
610 for (Int_t jdigit=idigit+1; jdigit<ndigits; jdigit++) {
611 AliPHOSCPVDigit *cpvDigit2 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(jdigit);
612 Float_t x2 = cpvDigit2->GetXpad() ;
613 Float_t z2 = cpvDigit2->GetYpad() ;
614 if (x1==x2 && z1==z2) {
615 Float_t qsum = cpvDigit1->GetQpad() + cpvDigit2->GetQpad() ;
616 cpvDigit2->SetQpad(qsum) ;
617 cpvDigits->RemoveAt(idigit) ;
621 cpvDigits->Compress() ;
623 // 3. add digits to temporary hit list fTmpHits
625 ndigits = cpvDigits->GetEntriesFast();
626 for (idigit=0; idigit<ndigits; idigit++) {
627 AliPHOSCPVDigit *cpvDigit = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit);
628 relid[0] = moduleNumber + 1 ; // CPV (or PHOS) module number
629 relid[1] =-1 ; // means CPV
630 relid[2] = cpvDigit->GetXpad() ; // column number of a pad
631 relid[3] = cpvDigit->GetYpad() ; // row number of a pad
633 // get the absolute Id number
634 fGeom->RelToAbsNumbering(relid, absid) ;
636 // add current digit to the temporary hit list
640 xyze[3] = cpvDigit->GetQpad() ; // amplitude in a pad
641 primary = -1; // No need in primary for CPV
642 AddHit(fIshunt, primary, tracknumber, absid, xyze, trackpid, pmom, xyd);
644 if (cpvDigit->GetQpad() > 0.02) {
645 xmean += cpvDigit->GetQpad() * (cpvDigit->GetXpad() + 0.5);
646 zmean += cpvDigit->GetQpad() * (cpvDigit->GetYpad() + 0.5);
647 qsum += cpvDigit->GetQpad();
652 } // end of IHEP configuration
655 if(gMC->CurrentVolID(copy) == gMC->VolId("PXTL") ) { // We are inside a PBWO crystal
656 gMC->TrackPosition(pos) ;
660 xyze[3] = gMC->Edep() ;
663 if ( (xyze[3] != 0) || entered ) { // Track is inside the crystal and deposits some energy or just entered
665 gMC->CurrentVolOffID(10, relid[0]) ; // get the PHOS module number ;
667 if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(10),"PHO1") == 0 )
668 relid[0] += fGeom->GetNModules() - fGeom->GetNPPSDModules();
670 relid[1] = 0 ; // means PBW04
671 gMC->CurrentVolOffID(4, relid[2]) ; // get the row number inside the module
672 gMC->CurrentVolOffID(3, relid[3]) ; // get the cell number inside the module
674 // get the absolute Id number
675 fGeom->RelToAbsNumbering(relid, absid) ;
677 // add current hit to the hit list
678 AddHit(fIshunt, primary,tracknumber, absid, xyze, trackpid,pmom, xyd);
681 } // there is deposited energy
682 } // we are inside a PHOS Xtal
687 //____________________________________________________________________________
688 void AliPHOSv1::CPVDigitize (TLorentzVector p, Float_t *zxhit, Int_t moduleNumber, TClonesArray *cpvDigits)
690 // ------------------------------------------------------------------------
691 // Digitize one CPV hit:
692 // On input take exact 4-momentum p and position zxhit of the hit,
693 // find the pad response around this hit and
694 // put the amplitudes in the pads into array digits
696 // Author: Yuri Kharlov (after Serguei Sadovsky)
698 // ------------------------------------------------------------------------
700 const Float_t kCelWr = fGeom->GetPadSizePhi()/2; // Distance between wires (2 wires above 1 pad)
701 const Float_t kDetR = 0.1; // Relative energy fluctuation in track for 100 e-
702 const Float_t kdEdx = 4.0; // Average energy loss in CPV;
703 const Int_t kNgamz = 5; // Ionization size in Z
704 const Int_t kNgamx = 9; // Ionization size in Phi
705 const Float_t kNoise = 0.03; // charge noise in one pad
709 // Just a reminder on axes notation in the CPV module:
710 // axis Z goes along the beam
711 // axis X goes across the beam in the module plane
712 // axis Y is a normal to the module plane showing from the IP
714 Float_t hitX = zxhit[0];
715 Float_t hitZ =-zxhit[1];
718 Float_t pNorm = p.Py();
719 Float_t eloss = kdEdx;
721 // cout << "CPVDigitize: YVK : "<<hitX<<" "<<hitZ<<" | "<<pX<<" "<<pZ<<" "<<pNorm<<endl;
723 Float_t dZY = pZ/pNorm * fGeom->GetCPVGasThickness();
724 Float_t dXY = pX/pNorm * fGeom->GetCPVGasThickness();
725 gRandom->Rannor(rnor1,rnor2);
726 eloss *= (1 + kDetR*rnor1) *
727 TMath::Sqrt((1 + ( pow(dZY,2) + pow(dXY,2) ) / pow(fGeom->GetCPVGasThickness(),2)));
728 Float_t zhit1 = hitZ + fGeom->GetCPVActiveSize(1)/2 - dZY/2;
729 Float_t xhit1 = hitX + fGeom->GetCPVActiveSize(0)/2 - dXY/2;
730 Float_t zhit2 = zhit1 + dZY;
731 Float_t xhit2 = xhit1 + dXY;
733 Int_t iwht1 = (Int_t) (xhit1 / kCelWr); // wire (x) coordinate "in"
734 Int_t iwht2 = (Int_t) (xhit2 / kCelWr); // wire (x) coordinate "out"
738 if (iwht1==iwht2) { // incline 1-wire hit
740 zxe[0][0] = (zhit1 + zhit2 - dZY*0.57735) / 2;
741 zxe[1][0] = (iwht1 + 0.5) * kCelWr;
743 zxe[0][1] = (zhit1 + zhit2 + dZY*0.57735) / 2;
744 zxe[1][1] = (iwht1 + 0.5) * kCelWr;
747 else if (TMath::Abs(iwht1-iwht2) != 1) { // incline 3-wire hit
749 Int_t iwht3 = (iwht1 + iwht2) / 2;
750 Float_t xwht1 = (iwht1 + 0.5) * kCelWr; // wire 1
751 Float_t xwht2 = (iwht2 + 0.5) * kCelWr; // wire 2
752 Float_t xwht3 = (iwht3 + 0.5) * kCelWr; // wire 3
753 Float_t xwr13 = (xwht1 + xwht3) / 2; // center 13
754 Float_t xwr23 = (xwht2 + xwht3) / 2; // center 23
755 Float_t dxw1 = xhit1 - xwr13;
756 Float_t dxw2 = xhit2 - xwr23;
757 Float_t egm1 = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
758 Float_t egm2 = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
759 Float_t egm3 = kCelWr / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
760 zxe[0][0] = (dXY*(xwr13-xwht1)/dXY + zhit1 + zhit1) / 2;
762 zxe[2][0] = eloss * egm1;
763 zxe[0][1] = (dXY*(xwr23-xwht1)/dXY + zhit1 + zhit2) / 2;
765 zxe[2][1] = eloss * egm2;
766 zxe[0][2] = dXY*(xwht3-xwht1)/dXY + zhit1;
768 zxe[2][2] = eloss * egm3;
770 else { // incline 2-wire hit
772 Float_t xwht1 = (iwht1 + 0.5) * kCelWr;
773 Float_t xwht2 = (iwht2 + 0.5) * kCelWr;
774 Float_t xwr12 = (xwht1 + xwht2) / 2;
775 Float_t dxw1 = xhit1 - xwr12;
776 Float_t dxw2 = xhit2 - xwr12;
777 Float_t egm1 = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) );
778 Float_t egm2 = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) );
779 zxe[0][0] = (zhit1 + zhit2 - dZY*egm1) / 2;
781 zxe[2][0] = eloss * egm1;
782 zxe[0][1] = (zhit1 + zhit2 + dZY*egm2) / 2;
784 zxe[2][1] = eloss * egm2;
787 // Finite size of ionization region
789 Int_t nCellZ = fGeom->GetNumberOfCPVPadsZ();
790 Int_t nCellX = fGeom->GetNumberOfCPVPadsPhi();
791 Int_t nz3 = (kNgamz+1)/2;
792 Int_t nx3 = (kNgamx+1)/2;
793 cpvDigits->Expand(nIter*kNgamx*kNgamz);
794 TClonesArray &ldigits = *(TClonesArray *)cpvDigits;
796 for (Int_t iter=0; iter<nIter; iter++) {
798 Float_t zhit = zxe[0][iter];
799 Float_t xhit = zxe[1][iter];
800 Float_t qhit = zxe[2][iter];
801 Float_t zcell = zhit / fGeom->GetPadSizeZ();
802 Float_t xcell = xhit / fGeom->GetPadSizePhi();
803 if ( zcell<=0 || xcell<=0 ||
804 zcell>=nCellZ || xcell>=nCellX) return;
805 Int_t izcell = (Int_t) zcell;
806 Int_t ixcell = (Int_t) xcell;
807 Float_t zc = zcell - izcell - 0.5;
808 Float_t xc = xcell - ixcell - 0.5;
809 for (Int_t iz=1; iz<=kNgamz; iz++) {
810 Int_t kzg = izcell + iz - nz3;
811 if (kzg<=0 || kzg>nCellZ) continue;
812 Float_t zg = (Float_t)(iz-nz3) - zc;
813 for (Int_t ix=1; ix<=kNgamx; ix++) {
814 Int_t kxg = ixcell + ix - nx3;
815 if (kxg<=0 || kxg>nCellX) continue;
816 Float_t xg = (Float_t)(ix-nx3) - xc;
818 // Now calculate pad response
819 Float_t qpad = CPVPadResponseFunction(qhit,zg,xg);
820 qpad += kNoise*rnor2;
821 if (qpad<0) continue;
823 // Fill the array with pad response ID and amplitude
824 new(ldigits[cpvDigits->GetEntriesFast()]) AliPHOSCPVDigit(kxg,kzg,qpad);
830 //____________________________________________________________________________
831 Float_t AliPHOSv1::CPVPadResponseFunction(Float_t qhit, Float_t zhit, Float_t xhit) {
832 // ------------------------------------------------------------------------
833 // Calculate the amplitude in one CPV pad using the
834 // cumulative pad response function
835 // Author: Yuri Kharlov (after Serguei Sadovski)
837 // ------------------------------------------------------------------------
839 Double_t dz = fGeom->GetPadSizeZ() / 2;
840 Double_t dx = fGeom->GetPadSizePhi() / 2;
841 Double_t z = zhit * fGeom->GetPadSizeZ();
842 Double_t x = xhit * fGeom->GetPadSizePhi();
843 Double_t amplitude = qhit *
844 (CPVCumulPadResponse(z+dz,x+dx) - CPVCumulPadResponse(z+dz,x-dx) -
845 CPVCumulPadResponse(z-dz,x+dx) + CPVCumulPadResponse(z-dz,x-dx));
846 return (Float_t)amplitude;
849 //____________________________________________________________________________
850 Double_t AliPHOSv1::CPVCumulPadResponse(Double_t x, Double_t y) {
851 // ------------------------------------------------------------------------
852 // Cumulative pad response function
853 // It includes several terms from the CF decomposition in electrostatics
854 // Note: this cumulative function is wrong since omits some terms
855 // but the cell amplitude obtained with it is correct because
856 // these omitting terms cancel
857 // Author: Yuri Kharlov (after Serguei Sadovski)
859 // ------------------------------------------------------------------------
861 const Double_t kA=1.0;
862 const Double_t kB=0.7;
864 Double_t r2 = x*x + y*y;
866 Double_t cumulPRF = 0;
867 for (Int_t i=0; i<=4; i++) {
868 Double_t b1 = (2*i + 1) * kB;
869 cumulPRF += TMath::Power(-1,i) * TMath::ATan( xy / (b1*TMath::Sqrt(b1*b1 + r2)) );
871 cumulPRF *= kA/(2*TMath::Pi());