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;
73 //____________________________________________________________________________
74 AliPHOSv1::AliPHOSv1(const char *name, const char *title):
77 // ctor : title is used to identify the layout
78 // GPS2 = 5 modules (EMC + PPSD)
79 // IHEP = 5 modules (EMC + CPV )
80 // MIXT = 4 modules (EMC + CPV ) and 1 module (EMC + PPSD)
83 // - fHits (the "normal" one), which retains the hits associated with
84 // the current primary particle being tracked
85 // (this array is reset after each primary has been tracked).
88 fPinElectronicNoise = 0.010 ;
89 fDigitThreshold = 0.01 ; // 1 GeV
91 fDigitizeB = 10000000. ;
94 // We do not want to save in TreeH the raw hits
95 // But save the cumulated hits instead (need to create the branch myself)
96 // It is put in the Digit Tree because the TreeH is filled after each primary
97 // and the TreeD at the end of the event (branch is set in FinishEvent() ).
99 fHits= new TClonesArray("AliPHOSHit",1000) ;
103 fReconstructioner = 0;
104 fTrackSegmentMaker = 0;
106 fIshunt = 1 ; // All hits are associated with primary particles
110 //____________________________________________________________________________
111 AliPHOSv1::AliPHOSv1(AliPHOSReconstructioner * Reconstructioner, const char *name, const char *title):
112 AliPHOSv0(name,title)
114 // ctor : title is used to identify the layout
115 // GPS2 = 5 modules (EMC + PPSD)
117 fPinElectronicNoise = 0.010 ;
119 // We do not want to save in TreeH the raw hits
122 fHits= new TClonesArray("AliPHOSHit",1000) ;
126 fIshunt = 1 ; // All hits are associated with primary particles
128 // gets an instance of the geometry parameters class
129 fGeom = AliPHOSGeometry::GetInstance(title, "") ;
131 if (fGeom->IsInitialized() )
132 cout << "AliPHOS" << Version() << " : PHOS geometry intialized for " << fGeom->GetName() << endl ;
134 cout << "AliPHOS" << Version() << " : PHOS geometry initialization failed !" << endl ;
136 // Defining the PHOS Reconstructioner
138 fReconstructioner = Reconstructioner ;
142 //____________________________________________________________________________
143 AliPHOSv1::~AliPHOSv1()
165 if ( fEmcRecPoints ) {
166 fEmcRecPoints->Delete() ;
167 delete fEmcRecPoints ;
171 if ( fPpsdRecPoints ) {
172 fPpsdRecPoints->Delete() ;
173 delete fPpsdRecPoints ;
177 if ( fTrackSegments ) {
178 fTrackSegments->Delete() ;
179 delete fTrackSegments ;
185 //____________________________________________________________________________
186 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)
188 // Add a hit to the hit list.
189 // A PHOS hit is the sum of all hits in a single crystal
190 // or in a single PPSD gas cell
195 Bool_t deja = kFALSE ;
197 newHit = new AliPHOSHit(shunt, primary, tracknumber, Id, hits, trackpid, p, lpos) ;
199 for ( hitCounter = fNhits-1 ; hitCounter >= 0 && !deja ; hitCounter-- ) {
200 curHit = (AliPHOSHit*) (*fHits)[hitCounter] ;
201 if( *curHit == *newHit ) {
202 *curHit = *curHit + *newHit ;
208 new((*fHits)[fNhits]) AliPHOSHit(*newHit) ;
215 //____________________________________________________________________________
216 void AliPHOSv1::Hits2SDigits(){
217 //Collects all hits in the same active volume into digit
222 AliPHOSDigit * newdigit ;
223 AliPHOSDigit * curdigit ;
224 Bool_t deja = kFALSE ;
228 for (itrack=0; itrack<gAlice->GetNtrack(); itrack++){
230 //=========== Get the Hits Tree for the Primary track itrack
232 gAlice->TreeH()->GetEvent(itrack);
235 for ( i = 0 ; i < fHits->GetEntries() ; i++ ) {
236 hit = (AliPHOSHit*)fHits->At(i) ;
238 // Assign primary number only if contribution is significant
239 if( hit->GetEnergy() > fDigitThreshold)
240 newdigit = new AliPHOSDigit( hit->GetPrimary(), hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
242 newdigit = new AliPHOSDigit( -1 , hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
245 for ( j = 0 ; j < fnSdigits ; j++) {
246 curdigit = (AliPHOSDigit*) fSDigits->At(j) ;
247 if ( *curdigit == *newdigit) {
248 *curdigit = *curdigit + *newdigit ;
254 new((*fSDigits)[fnSdigits]) AliPHOSDigit(* newdigit) ;
261 } // loop over tracks
265 fnSdigits = fSDigits->GetEntries() ;
266 fSDigits->Expand(fnSdigits) ;
268 for (i = 0 ; i < fnSdigits ; i++) {
269 AliPHOSDigit * digit = (AliPHOSDigit *) fSDigits->At(i) ;
270 digit->SetIndexInList(i) ;
273 gAlice->TreeS()->Fill() ;
274 gAlice->TreeS()->Write(0,TObject::kOverwrite) ;
278 //____________________________________________________________________________
279 void AliPHOSv1::SDigits2Digits(){
280 //Adds noise to the summable digits and removes everething below thresholds
281 //Note, that sDigits should be SORTED in accordance with abs ID.
284 gAlice->TreeS()->GetEvent(0) ;
286 // First calculate noise induced by the PIN diode of the PbWO crystals
287 Int_t iCurSDigit = 0 ;
289 //we assume, that there is al least one EMC digit...
290 if(fSDigits->GetEntries() == 0) {
291 cout << "PHOS::SDigits2Digits> No SDigits !!! Do not produce Digits " << endl ;
295 Int_t idCurSDigit = ((AliPHOSDigit *)fSDigits->At(0))->GetId() ;
298 for(absID = 1; absID < fGeom->GetNModules()*fGeom->GetNPhi()*fGeom->GetNZ(); absID++){
299 Float_t noise = gRandom->Gaus(0., fPinElectronicNoise) ;
300 if(absID < idCurSDigit ){
301 if(noise >fDigitThreshold ){
302 new((*fDigits)[fNdigits]) AliPHOSDigit( -1,absID,Digitize(noise) ) ;
306 else{ //add noise and may be remove the true hit
307 Float_t signal = noise + Calibrate(((AliPHOSDigit *)fSDigits->At(iCurSDigit))->GetAmp()) ;
308 if( signal >fDigitThreshold ){
309 AliPHOSDigit * digit = (AliPHOSDigit*) fSDigits->At(iCurSDigit) ;
310 new((*fDigits)[fNdigits]) AliPHOSDigit( *digit ) ;
311 ((AliPHOSDigit *)fDigits->At(fNdigits))->SetAmp(Digitize(signal));
315 if(iCurSDigit < fSDigits->GetEntries()-1){
317 idCurSDigit = ((AliPHOSDigit*)fSDigits->At(iCurSDigit))->GetId() ;
320 idCurSDigit = 10000000; //no real hits left
325 //remove PPSD/CPV digits below thresholds
327 for(idigit = iCurSDigit; idigit < fSDigits->GetEntries() ; idigit++){ //loop over CPV/PPSD digits
329 AliPHOSDigit * digit = (AliPHOSDigit *) fSDigits->At(idigit) ;
330 Float_t ene = Calibrate(digit->GetAmp()) ;
333 fGeom->AbsToRelNumbering(digit->GetId(), relid) ;
334 if ( relid[0] > fGeom->GetNCPVModules() ){ //ppsd
335 if ( ( (relid[1] > 0) && (ene > fPpsdEnergyThreshold)) || //PPSD digit
336 ( (relid[1] < 0) && (ene > fCpvEnergyThreshold ) ) ) //CPV digit
337 new((*fDigits)[fNdigits]) AliPHOSDigit( *digit ) ;
342 fDigits->Compress() ;
344 fNdigits = fDigits->GetEntries() ;
345 fDigits->Expand(fNdigits) ;
348 for (i = 0 ; i < fNdigits ; i++) {
349 AliPHOSDigit * digit = (AliPHOSDigit *) fDigits->At(i) ;
350 digit->SetIndexInList(i) ;
353 gAlice->TreeD()->Fill() ;
355 gAlice->TreeD()->Write(0,TObject::kOverwrite) ;
359 //___________________________________________________________________________
360 void AliPHOSv1::MakeBranch(Option_t* opt, char *file)
365 // Create new branche in the current Root Tree in the digit Tree
366 AliDetector::MakeBranch(opt) ;
369 cH = strstr(opt,"S");
370 //Create a branch for SDigits
373 sprintf(branchname,"%s",GetName());
377 fSDigits = new TClonesArray("AliPHOSDigit",1000);
380 gAlice->MakeBranchInTree(gAlice->TreeS(),branchname,&fSDigits,fBufferSize,file);
383 cH = strstr(opt,"D");
384 //Create a branch for Digits
387 sprintf(branchname,"%s",GetName());
392 fDigits = new TClonesArray("AliPHOSDigit",1000);
395 gAlice->MakeBranchInTree(gAlice->TreeD(),branchname,&fDigits,fBufferSize,file);
398 cH = strstr(opt,"R");
399 //Create a branch for Reconstruction
403 Int_t splitlevel = 0 ;
406 fEmcRecPoints->Delete() ;
408 fEmcRecPoints = new TObjArray(100) ;
410 if ( fEmcRecPoints && gAlice->TreeR() ) {
411 sprintf(branchname,"%sEmcRP",GetName()) ;
412 gAlice->TreeR()->Branch(branchname, "TObjArray", &fEmcRecPoints, fBufferSize, splitlevel) ;
416 fPpsdRecPoints->Delete() ;
418 fPpsdRecPoints = new TObjArray(100) ;
420 if ( fPpsdRecPoints && gAlice->TreeR() ) {
421 sprintf(branchname,"%sPpsdRP",GetName()) ;
422 gAlice->TreeR()->Branch(branchname, "TObjArray", &fPpsdRecPoints, fBufferSize, splitlevel) ;
426 fTrackSegments->Delete() ;
428 fTrackSegments = new TClonesArray("AliPHOSTrackSegment",100) ;
430 if ( fTrackSegments && gAlice->TreeR() ) {
431 sprintf(branchname,"%sTS",GetName()) ;
432 gAlice->TreeR()->Branch(branchname, &fTrackSegments, fBufferSize) ;
436 fRecParticles->Delete() ;
438 fRecParticles = new TClonesArray("AliPHOSRecParticle",100) ;
440 if ( fRecParticles && gAlice->TreeR() ) {
441 sprintf(branchname,"%sRP",GetName()) ;
442 gAlice->TreeR()->Branch(branchname, &fRecParticles, fBufferSize) ;
449 //_____________________________________________________________________________
450 void AliPHOSv1::Reconstruction(AliPHOSReconstructioner * Reconstructioner)
452 // 1. Reinitializes the existing RecPoint, TrackSegment, and RecParticles Lists and
453 // 2. Creates TreeR with a branch for each list
454 // 3. Steers the reconstruction processes
455 // 4. Saves the 3 lists in TreeR
456 // 5. Write the Tree to File
458 fReconstructioner = Reconstructioner ;
462 // gAlice->MakeTree("R") ;
468 fReconstructioner->Make(fDigits, fEmcRecPoints, fPpsdRecPoints, fTrackSegments, fRecParticles);
470 printf("Reconstruction: %d %d %d %d\n",
471 fEmcRecPoints->GetEntries(),fPpsdRecPoints->GetEntries(),
472 fTrackSegments->GetEntries(),fRecParticles->GetEntries());
474 // 4. Expand or Shrink the arrays to the proper size
478 size = fEmcRecPoints->GetEntries() ;
479 fEmcRecPoints->Expand(size) ;
481 size = fPpsdRecPoints->GetEntries() ;
482 fPpsdRecPoints->Expand(size) ;
484 size = fTrackSegments->GetEntries() ;
485 fTrackSegments->Expand(size) ;
487 size = fRecParticles->GetEntries() ;
488 fRecParticles->Expand(size) ;
490 gAlice->TreeR()->Fill() ;
493 gAlice->TreeR()->Write(0,TObject::kOverwrite) ;
495 // Deleting reconstructed objects
496 ResetReconstruction();
500 //____________________________________________________________________________
501 void AliPHOSv1::ResetReconstruction()
503 // Deleting reconstructed objects
505 if ( fEmcRecPoints ) fEmcRecPoints ->Delete();
506 if ( fPpsdRecPoints ) fPpsdRecPoints->Delete();
507 if ( fTrackSegments ) fTrackSegments->Delete();
508 if ( fRecParticles ) fRecParticles ->Delete();
512 //____________________________________________________________________________
514 void AliPHOSv1::StepManager(void)
516 // Accumulates hits as long as the track stays in a single crystal or PPSD gas Cell
518 Int_t relid[4] ; // (box, layer, row, column) indices
519 Int_t absid ; // absolute cell ID number
520 Float_t xyze[4]={0,0,0,0} ; // position wrt MRS and energy deposited
521 TLorentzVector pos ; // Lorentz vector of the track current position
522 TLorentzVector pmom = 0 ; //momentum of the particle initiated hit
523 Float_t xyd[3]={0,0,0} ; //local posiiton of the entering
524 Bool_t entered = kFALSE ;
527 Int_t tracknumber = gAlice->CurrentTrack() ;
528 Int_t primary = gAlice->GetPrimary( gAlice->CurrentTrack() );
529 TString name = fGeom->GetName() ;
532 if( gMC->IsTrackEntering() ){ // create hit with position and momentum of new particle,
533 // but may be without energy deposition
535 // Current position of the hit in the local ref. system
536 gMC -> TrackPosition(pos);
537 Float_t xyzm[3], xyzd[3] ;
539 for (i=0; i<3; i++) xyzm[i] = pos[i];
540 gMC -> Gmtod (xyzm, xyzd, 1); // transform coordinate from master to daughter system
546 // Current momentum of the hit's track in the local ref. system
547 gMC -> TrackMomentum(pmom);
548 Float_t pm[3], pd[3];
549 for (i=0; i<3; i++) pm[i] = pmom[i];
550 gMC -> Gmtod (pm, pd, 2); // transform 3-momentum from master to daughter system
555 trackpid = gMC->TrackPid();
556 entered = kTRUE ; // Mark to create hit even withou energy deposition
561 if ( name == "GPS2" || name == "MIXT" ) { // ======> CPV is a GPS' PPSD
563 if( gMC->CurrentVolID(copy) == gMC->VolId("GCEL") ) // We are inside a gas cell
565 gMC->TrackPosition(pos) ;
569 xyze[3] = gMC->Edep() ;
571 if ( (xyze[3] != 0) || entered ) { // there is deposited energy or new particle entering PPSD
572 gMC->CurrentVolOffID(5, relid[0]) ; // get the PHOS Module number
573 if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(5),"PHO1") == 0 ){
574 relid[0] += fGeom->GetNModules() - fGeom->GetNPPSDModules();
576 gMC->CurrentVolOffID(3, relid[1]) ; // get the Micromegas Module number
577 // 1-> fGeom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() upper
578 // > fGeom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() lower
579 gMC->CurrentVolOffID(1, relid[2]) ; // get the row number of the cell
580 gMC->CurrentVolID(relid[3]) ; // get the column number
582 // get the absolute Id number
584 fGeom->RelToAbsNumbering(relid, absid) ;
586 // add current hit to the hit list
587 AddHit(fIshunt, primary, tracknumber, absid, xyze, trackpid, pmom, xyd);
590 } // there is deposited energy
591 } // We are inside the gas of the CPV
592 } // GPS2 configuration
594 if ( name == "IHEP" || name == "MIXT" ) { // ======> CPV is a IHEP's one
596 // Yuri Kharlov, 28 September 2000
598 if( gMC->CurrentVolID(copy) == gMC->VolId("CPVQ") &&
600 gMC->TrackCharge() != 0) {
602 // Digitize the current CPV hit:
604 // 1. find pad response and
607 gMC->CurrentVolOffID(3,moduleNumber);
611 TClonesArray *cpvDigits = new TClonesArray("AliPHOSCPVDigit",0); // array of digits for current hit
612 CPVDigitize(pmom,xyd,moduleNumber,cpvDigits);
617 Int_t idigit,ndigits;
619 // 2. go through the current digit list and sum digits in pads
621 ndigits = cpvDigits->GetEntriesFast();
622 for (idigit=0; idigit<ndigits-1; idigit++) {
623 AliPHOSCPVDigit *cpvDigit1 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit);
624 Float_t x1 = cpvDigit1->GetXpad() ;
625 Float_t z1 = cpvDigit1->GetYpad() ;
626 for (Int_t jdigit=idigit+1; jdigit<ndigits; jdigit++) {
627 AliPHOSCPVDigit *cpvDigit2 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(jdigit);
628 Float_t x2 = cpvDigit2->GetXpad() ;
629 Float_t z2 = cpvDigit2->GetYpad() ;
630 if (x1==x2 && z1==z2) {
631 Float_t qsum = cpvDigit1->GetQpad() + cpvDigit2->GetQpad() ;
632 cpvDigit2->SetQpad(qsum) ;
633 cpvDigits->RemoveAt(idigit) ;
637 cpvDigits->Compress() ;
639 // 3. add digits to temporary hit list fTmpHits
641 ndigits = cpvDigits->GetEntriesFast();
642 for (idigit=0; idigit<ndigits; idigit++) {
643 AliPHOSCPVDigit *cpvDigit = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit);
644 relid[0] = moduleNumber + 1 ; // CPV (or PHOS) module number
645 relid[1] =-1 ; // means CPV
646 relid[2] = cpvDigit->GetXpad() ; // column number of a pad
647 relid[3] = cpvDigit->GetYpad() ; // row number of a pad
649 // get the absolute Id number
650 fGeom->RelToAbsNumbering(relid, absid) ;
652 // add current digit to the temporary hit list
656 xyze[3] = cpvDigit->GetQpad() ; // amplitude in a pad
657 primary = -1; // No need in primary for CPV
658 AddHit(fIshunt, primary, tracknumber, absid, xyze, trackpid, pmom, xyd);
660 if (cpvDigit->GetQpad() > 0.02) {
661 xmean += cpvDigit->GetQpad() * (cpvDigit->GetXpad() + 0.5);
662 zmean += cpvDigit->GetQpad() * (cpvDigit->GetYpad() + 0.5);
663 qsum += cpvDigit->GetQpad();
668 } // end of IHEP configuration
671 if(gMC->CurrentVolID(copy) == gMC->VolId("PXTL") ) { // We are inside a PBWO crystal
672 gMC->TrackPosition(pos) ;
676 xyze[3] = gMC->Edep() ;
679 if ( (xyze[3] != 0) || entered ) { // Track is inside the crystal and deposits some energy or just entered
681 gMC->CurrentVolOffID(10, relid[0]) ; // get the PHOS module number ;
683 if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(10),"PHO1") == 0 )
684 relid[0] += fGeom->GetNModules() - fGeom->GetNPPSDModules();
686 relid[1] = 0 ; // means PBW04
687 gMC->CurrentVolOffID(4, relid[2]) ; // get the row number inside the module
688 gMC->CurrentVolOffID(3, relid[3]) ; // get the cell number inside the module
690 // get the absolute Id number
691 fGeom->RelToAbsNumbering(relid, absid) ;
693 // add current hit to the hit list
694 AddHit(fIshunt, primary,tracknumber, absid, xyze, trackpid,pmom, xyd);
697 } // there is deposited energy
698 } // we are inside a PHOS Xtal
703 //____________________________________________________________________________
704 void AliPHOSv1::CPVDigitize (TLorentzVector p, Float_t *zxhit, Int_t moduleNumber, TClonesArray *cpvDigits)
706 // ------------------------------------------------------------------------
707 // Digitize one CPV hit:
708 // On input take exact 4-momentum p and position zxhit of the hit,
709 // find the pad response around this hit and
710 // put the amplitudes in the pads into array digits
712 // Author: Yuri Kharlov (after Serguei Sadovsky)
714 // ------------------------------------------------------------------------
716 const Float_t kCelWr = fGeom->GetPadSizePhi()/2; // Distance between wires (2 wires above 1 pad)
717 const Float_t kDetR = 0.1; // Relative energy fluctuation in track for 100 e-
718 const Float_t kdEdx = 4.0; // Average energy loss in CPV;
719 const Int_t kNgamz = 5; // Ionization size in Z
720 const Int_t kNgamx = 9; // Ionization size in Phi
721 const Float_t kNoise = 0.03; // charge noise in one pad
725 // Just a reminder on axes notation in the CPV module:
726 // axis Z goes along the beam
727 // axis X goes across the beam in the module plane
728 // axis Y is a normal to the module plane showing from the IP
730 Float_t hitX = zxhit[0];
731 Float_t hitZ =-zxhit[1];
734 Float_t pNorm = p.Py();
735 Float_t eloss = kdEdx;
737 // cout << "CPVDigitize: YVK : "<<hitX<<" "<<hitZ<<" | "<<pX<<" "<<pZ<<" "<<pNorm<<endl;
739 Float_t dZY = pZ/pNorm * fGeom->GetCPVGasThickness();
740 Float_t dXY = pX/pNorm * fGeom->GetCPVGasThickness();
741 gRandom->Rannor(rnor1,rnor2);
742 eloss *= (1 + kDetR*rnor1) *
743 TMath::Sqrt((1 + ( pow(dZY,2) + pow(dXY,2) ) / pow(fGeom->GetCPVGasThickness(),2)));
744 Float_t zhit1 = hitZ + fGeom->GetCPVActiveSize(1)/2 - dZY/2;
745 Float_t xhit1 = hitX + fGeom->GetCPVActiveSize(0)/2 - dXY/2;
746 Float_t zhit2 = zhit1 + dZY;
747 Float_t xhit2 = xhit1 + dXY;
749 Int_t iwht1 = (Int_t) (xhit1 / kCelWr); // wire (x) coordinate "in"
750 Int_t iwht2 = (Int_t) (xhit2 / kCelWr); // wire (x) coordinate "out"
754 if (iwht1==iwht2) { // incline 1-wire hit
756 zxe[0][0] = (zhit1 + zhit2 - dZY*0.57735) / 2;
757 zxe[1][0] = (iwht1 + 0.5) * kCelWr;
759 zxe[0][1] = (zhit1 + zhit2 + dZY*0.57735) / 2;
760 zxe[1][1] = (iwht1 + 0.5) * kCelWr;
763 else if (TMath::Abs(iwht1-iwht2) != 1) { // incline 3-wire hit
765 Int_t iwht3 = (iwht1 + iwht2) / 2;
766 Float_t xwht1 = (iwht1 + 0.5) * kCelWr; // wire 1
767 Float_t xwht2 = (iwht2 + 0.5) * kCelWr; // wire 2
768 Float_t xwht3 = (iwht3 + 0.5) * kCelWr; // wire 3
769 Float_t xwr13 = (xwht1 + xwht3) / 2; // center 13
770 Float_t xwr23 = (xwht2 + xwht3) / 2; // center 23
771 Float_t dxw1 = xhit1 - xwr13;
772 Float_t dxw2 = xhit2 - xwr23;
773 Float_t egm1 = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
774 Float_t egm2 = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
775 Float_t egm3 = kCelWr / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
776 zxe[0][0] = (dXY*(xwr13-xwht1)/dXY + zhit1 + zhit1) / 2;
778 zxe[2][0] = eloss * egm1;
779 zxe[0][1] = (dXY*(xwr23-xwht1)/dXY + zhit1 + zhit2) / 2;
781 zxe[2][1] = eloss * egm2;
782 zxe[0][2] = dXY*(xwht3-xwht1)/dXY + zhit1;
784 zxe[2][2] = eloss * egm3;
786 else { // incline 2-wire hit
788 Float_t xwht1 = (iwht1 + 0.5) * kCelWr;
789 Float_t xwht2 = (iwht2 + 0.5) * kCelWr;
790 Float_t xwr12 = (xwht1 + xwht2) / 2;
791 Float_t dxw1 = xhit1 - xwr12;
792 Float_t dxw2 = xhit2 - xwr12;
793 Float_t egm1 = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) );
794 Float_t egm2 = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) );
795 zxe[0][0] = (zhit1 + zhit2 - dZY*egm1) / 2;
797 zxe[2][0] = eloss * egm1;
798 zxe[0][1] = (zhit1 + zhit2 + dZY*egm2) / 2;
800 zxe[2][1] = eloss * egm2;
803 // Finite size of ionization region
805 Int_t nCellZ = fGeom->GetNumberOfCPVPadsZ();
806 Int_t nCellX = fGeom->GetNumberOfCPVPadsPhi();
807 Int_t nz3 = (kNgamz+1)/2;
808 Int_t nx3 = (kNgamx+1)/2;
809 cpvDigits->Expand(nIter*kNgamx*kNgamz);
810 TClonesArray &ldigits = *(TClonesArray *)cpvDigits;
812 for (Int_t iter=0; iter<nIter; iter++) {
814 Float_t zhit = zxe[0][iter];
815 Float_t xhit = zxe[1][iter];
816 Float_t qhit = zxe[2][iter];
817 Float_t zcell = zhit / fGeom->GetPadSizeZ();
818 Float_t xcell = xhit / fGeom->GetPadSizePhi();
819 if ( zcell<=0 || xcell<=0 ||
820 zcell>=nCellZ || xcell>=nCellX) return;
821 Int_t izcell = (Int_t) zcell;
822 Int_t ixcell = (Int_t) xcell;
823 Float_t zc = zcell - izcell - 0.5;
824 Float_t xc = xcell - ixcell - 0.5;
825 for (Int_t iz=1; iz<=kNgamz; iz++) {
826 Int_t kzg = izcell + iz - nz3;
827 if (kzg<=0 || kzg>nCellZ) continue;
828 Float_t zg = (Float_t)(iz-nz3) - zc;
829 for (Int_t ix=1; ix<=kNgamx; ix++) {
830 Int_t kxg = ixcell + ix - nx3;
831 if (kxg<=0 || kxg>nCellX) continue;
832 Float_t xg = (Float_t)(ix-nx3) - xc;
834 // Now calculate pad response
835 Float_t qpad = CPVPadResponseFunction(qhit,zg,xg);
836 qpad += kNoise*rnor2;
837 if (qpad<0) continue;
839 // Fill the array with pad response ID and amplitude
840 new(ldigits[cpvDigits->GetEntriesFast()]) AliPHOSCPVDigit(kxg,kzg,qpad);
846 //____________________________________________________________________________
847 Float_t AliPHOSv1::CPVPadResponseFunction(Float_t qhit, Float_t zhit, Float_t xhit) {
848 // ------------------------------------------------------------------------
849 // Calculate the amplitude in one CPV pad using the
850 // cumulative pad response function
851 // Author: Yuri Kharlov (after Serguei Sadovski)
853 // ------------------------------------------------------------------------
855 Double_t dz = fGeom->GetPadSizeZ() / 2;
856 Double_t dx = fGeom->GetPadSizePhi() / 2;
857 Double_t z = zhit * fGeom->GetPadSizeZ();
858 Double_t x = xhit * fGeom->GetPadSizePhi();
859 Double_t amplitude = qhit *
860 (CPVCumulPadResponse(z+dz,x+dx) - CPVCumulPadResponse(z+dz,x-dx) -
861 CPVCumulPadResponse(z-dz,x+dx) + CPVCumulPadResponse(z-dz,x-dx));
862 return (Float_t)amplitude;
865 //____________________________________________________________________________
866 Double_t AliPHOSv1::CPVCumulPadResponse(Double_t x, Double_t y) {
867 // ------------------------------------------------------------------------
868 // Cumulative pad response function
869 // It includes several terms from the CF decomposition in electrostatics
870 // Note: this cumulative function is wrong since omits some terms
871 // but the cell amplitude obtained with it is correct because
872 // these omitting terms cancel
873 // Author: Yuri Kharlov (after Serguei Sadovski)
875 // ------------------------------------------------------------------------
877 const Double_t kA=1.0;
878 const Double_t kB=0.7;
880 Double_t r2 = x*x + y*y;
882 Double_t cumulPRF = 0;
883 for (Int_t i=0; i<=4; i++) {
884 Double_t b1 = (2*i + 1) * kB;
885 cumulPRF += TMath::Power(-1,i) * TMath::ATan( xy / (b1*TMath::Sqrt(b1*b1 + r2)) );
887 cumulPRF *= kA/(2*TMath::Pi());