]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSDebug.cxx
AliHMPIDDigitN no longer needed
[u/mrichter/AliRoot.git] / PHOS / AliPHOSDebug.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /* $Id$ */
17
18 //_________________________________________________________________________
19 // Implementation version v1 of PHOS Manager class 
20 //---
21 // Layout EMC + CPV  has name IHEP:
22 // Produces hits for CPV, cumulated hits
23 //---
24 //*-- Author: Yves Schutz (SUBATECH)
25
26
27 // --- ROOT system ---
28
29 #include "TBRIK.h"
30 #include "TNode.h"
31 #include "TRandom.h"
32 #include "TTree.h"
33
34
35 // --- Standard library ---
36
37 #include <stdio.h>
38 #include <string.h>
39 #include <stdlib.h>
40 #include <strstream.h>
41
42 // --- AliRoot header files ---
43
44 #include "AliPHOSv1.h"
45 #include "AliPHOSHit.h"
46 #include "AliPHOSDigit.h"
47 #include "AliPHOSReconstructor.h"
48 #include "AliRun.h"
49 #include "AliConst.h"
50
51 ClassImp(AliPHOSv1)
52
53 //____________________________________________________________________________
54 AliPHOSv1::AliPHOSv1():
55 AliPHOSv0()
56 {
57   // ctor
58  
59   fReconstructioner  = 0;
60   fTrackSegmentMaker = 0;
61
62 }
63
64 //____________________________________________________________________________
65 AliPHOSv1::AliPHOSv1(const char *name, const char *title):
66 AliPHOSv0(name,title) 
67 {
68   // ctor : title is used to identify the layout
69   //        IHEP = 5 modules (EMC + CPV )
70   //
71   // We store hits :
72   //   - fHits (the "normal" one), which retains the hits associated with
73   //     the current primary particle being tracked
74   //     (this array is reset after each primary has been tracked).
75   //
76
77   fPinElectronicNoise = 0.010 ;
78   fDigitThreshold      = 0.01 ;   // 1 GeV 
79   fDigitizeA= 0. ;             
80   fDigitizeB = 10000000. ;    
81
82
83   // We do not want to save in TreeH the raw hits
84   // But save the cumulated hits instead (need to create the branch myself)
85   // It is put in the Digit Tree because the TreeH is filled after each primary
86   // and the TreeD at the end of the event (branch is set in FinishEvent() ).
87   
88   fHits= new TClonesArray("AliPHOSHit",1000) ;
89
90   fNhits = 0 ;
91
92   fReconstructioner  = 0;
93   fTrackSegmentMaker = 0;
94
95   fIshunt     =  1 ; // All hits are associated with primary particles
96   
97 }
98
99 //____________________________________________________________________________
100 AliPHOSv1::AliPHOSv1(AliPHOSReconstructor * Reconstructioner, const char *name, const char *title):
101   AliPHOSv0(name,title)
102 {
103   // ctor : title is used to identify the layout
104
105   fPinElectronicNoise = 0.010 ;
106
107   // We do not want to save in TreeH the raw hits
108
109   fDigits = 0 ;
110   fHits= new TClonesArray("AliPHOSHit",1000) ;
111
112   fNhits = 0 ;
113
114   fIshunt     =  1 ; // All hits are associated with primary particles
115  
116   // gets an instance of the geometry parameters class  
117   AliPHOSGeometry::GetInstance(title, "") ; 
118
119   if (GetGeometry()->IsInitialized() ) 
120     Info("AliPHOSv1", "AliPHOS %d : PHOS geometry intialized for %s", Version(), GetGeometry()->GetName() );
121   else
122     Info("AliPHOSv1", "AliPHOS %d : PHOS geometry initialization failed !", Version() ) ;   
123
124   // Defining the PHOS Reconstructioner
125  
126  fReconstructioner = Reconstructioner ;
127
128 }
129
130 //____________________________________________________________________________
131 AliPHOSv1::~AliPHOSv1()
132 {
133   // dtor
134
135   if ( fHits) {
136     fHits->Delete() ; 
137     delete fHits ;
138     fHits = 0 ; 
139   }
140
141   if ( fSDigits) {
142     fSDigits->Delete() ; 
143     delete fSDigits ;
144     fSDigits = 0 ; 
145   }
146
147   if ( fDigits) {
148     fDigits->Delete() ; 
149     delete fDigits ;
150     fDigits = 0 ; 
151   }
152
153   if ( fEmcRecPoints ) {
154     fEmcRecPoints->Delete() ; 
155     delete fEmcRecPoints ; 
156     fEmcRecPoints = 0 ; 
157   }
158   
159   if ( fPpsdRecPoints ) { 
160     fPpsdRecPoints->Delete() ;
161     delete fPpsdRecPoints ;
162     fPpsdRecPoints = 0 ; 
163   }
164     
165   if ( fTrackSegments ) {
166     fTrackSegments->Delete() ; 
167     delete fTrackSegments ;
168     fTrackSegments = 0 ; 
169   }
170   
171 }
172
173 //____________________________________________________________________________
174 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)
175 {
176   // Add a hit to the hit list.
177   // A PHOS hit is the sum of all hits in a single crystal
178   //   or in a single PPSD gas cell
179
180   Int_t hitCounter ;
181   AliPHOSHit *newHit ;
182   AliPHOSHit *curHit ;
183   Bool_t deja = kFALSE ;
184
185   newHit = new AliPHOSHit(shunt, primary, tracknumber, Id, hits, trackpid, p, lpos) ;
186
187   for ( hitCounter = fNhits-1 ; hitCounter >= 0 && !deja ; hitCounter-- ) {
188     curHit = (AliPHOSHit*) (*fHits)[hitCounter] ;
189     if( *curHit == *newHit ) {
190       *curHit = *curHit + *newHit ;
191       deja = kTRUE ;
192     }
193   }
194          
195   if ( !deja ) {
196     new((*fHits)[fNhits]) AliPHOSHit(*newHit) ;
197     fNhits++ ;
198   }
199
200   delete newHit;
201 }
202
203 //____________________________________________________________________________
204 void AliPHOSv1::Hits2SDigits()
205 {
206   // Collects all hits in the same active volume into digit
207   // OBSOLETE replace by SDigitizer
208
209   Int_t i ;
210   Int_t j ; 
211   AliPHOSHit   * hit ;
212   AliPHOSDigit * newdigit ;
213   AliPHOSDigit * curdigit ;
214   Bool_t deja = kFALSE ; 
215   
216
217   Int_t itrack ;
218   for (itrack=0; itrack<gAlice->GetNtrack(); itrack++){
219         
220     //=========== Get the Hits Tree for the Primary track itrack
221     gAlice->ResetHits();
222     if (TreeH() == 0x0)
223      {
224        Error("Hits2SDigits","Can not find TreeH in the folder");
225        return;
226      }
227     TreeH()->GetEvent(itrack);
228       
229
230     for ( i = 0 ; i < fHits->GetEntries() ; i++ ) {
231       hit = (AliPHOSHit*)fHits->At(i) ;
232     
233       // Assign primary number only if contribution is significant
234       if( hit->GetEnergy() > fDigitThreshold)
235         newdigit = new AliPHOSDigit( hit->GetPrimary(), hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
236       else
237         newdigit = new AliPHOSDigit( -1               , hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
238       deja =kFALSE ;
239
240       for ( j = 0 ; j < fnSdigits ;  j++) { 
241         curdigit = (AliPHOSDigit*) fSDigits->At(j) ;
242         if ( *curdigit == *newdigit) {
243           *curdigit = *curdigit + *newdigit ; 
244           deja = kTRUE ; 
245         }
246       }
247
248       if ( !deja ) {
249         new((*fSDigits)[fnSdigits]) AliPHOSDigit(* newdigit) ;
250         fnSdigits++ ;  
251       }
252       
253       delete newdigit ;    
254     } 
255
256   } // loop over tracks
257
258   fSDigits->Sort() ;
259
260   fnSdigits = fSDigits->GetEntries() ;
261   fSDigits->Expand(fnSdigits) ;
262
263   for (i = 0 ; i < fnSdigits ; i++) { 
264     AliPHOSDigit * digit = (AliPHOSDigit *) fSDigits->At(i) ; 
265     digit->SetIndexInList(i) ;     
266   }
267
268   gAlice->TreeS()->Fill() ;
269   gAlice->TreeS()->Write(0,TObject::kOverwrite) ;
270  
271
272 }
273 //____________________________________________________________________________
274 void AliPHOSv1::SDigits2Digits()
275 {
276   // Adds noise to the summable digits and removes everething below thresholds
277   // Note, that sDigits should be SORTED in accordance with abs ID.
278   // OBSOLETE Replaced by Digitzer
279
280   gAlice->TreeS()->GetEvent(0) ;
281
282   // First calculate noise induced by the PIN diode of the PbWO crystals
283   Int_t iCurSDigit = 0 ;
284
285   //we assume, that there is al least one EMC digit...
286   if(fSDigits->GetEntries() == 0) {
287     Warning("SDigits2Digits", "No SDigits !!! Do not produce Digits ") ;
288     return ;
289   }
290
291   Int_t idCurSDigit = ((AliPHOSDigit *)fSDigits->At(0))->GetId() ;
292
293   Int_t absID ;
294   for(absID = 1; absID < GetGeometry()->GetNModules()*GetGeometry()->GetNPhi()*GetGeometry()->GetNZ(); absID++){
295     Float_t noise = gRandom->Gaus(0., fPinElectronicNoise) ; 
296     if(absID < idCurSDigit ){ 
297       if(noise >fDigitThreshold ){
298         new((*fDigits)[fNdigits]) AliPHOSDigit( -1,absID,Digitize(noise) ) ;
299         fNdigits++ ;  
300       }
301     }
302     else{ //add noise and may be remove the true hit
303       Float_t signal = noise + Calibrate(((AliPHOSDigit *)fSDigits->At(iCurSDigit))->GetAmp()) ;
304       if( signal >fDigitThreshold ){
305         AliPHOSDigit * digit = (AliPHOSDigit*) fSDigits->At(iCurSDigit) ;
306         new((*fDigits)[fNdigits]) AliPHOSDigit( *digit ) ;
307         ((AliPHOSDigit *)fDigits->At(fNdigits))->SetAmp(Digitize(signal));
308         fNdigits++ ;  
309       }
310
311       if(iCurSDigit < fSDigits->GetEntries()-1){
312         iCurSDigit++ ;
313         idCurSDigit = ((AliPHOSDigit*)fSDigits->At(iCurSDigit))->GetId() ;
314       }
315       else
316         idCurSDigit = 10000000; //no real hits left    
317     }
318     
319   }  
320
321   //remove PPSD/CPV digits below thresholds
322   Int_t idigit ;
323   for(idigit = iCurSDigit; idigit < fSDigits->GetEntries() ; idigit++){  //loop over CPV/PPSD digits
324     
325     AliPHOSDigit * digit = (AliPHOSDigit *) fSDigits->At(idigit) ; 
326     Float_t ene = Calibrate(digit->GetAmp()) ;
327     
328     Int_t relid[4] ; 
329     GetGeometry()->AbsToRelNumbering(digit->GetId(), relid) ; 
330     if ( relid[0] > GetGeometry()->GetNCPVModules() ){ //ppsd
331       if ( ( (relid[1] > 0) && (ene > fPpsdEnergyThreshold)) ||    //PPSD digit
332            ( (relid[1] < 0) && (ene > fCpvEnergyThreshold ) ) )    //CPV digit 
333         new((*fDigits)[fNdigits]) AliPHOSDigit( *digit ) ;
334         fNdigits++ ;
335     }
336   }    
337     
338   fDigits->Compress() ;  
339   
340   fNdigits = fDigits->GetEntries() ;
341   fDigits->Expand(fNdigits) ;
342
343   Int_t i ;
344   for (i = 0 ; i < fNdigits ; i++) { 
345     AliPHOSDigit * digit = (AliPHOSDigit *) fDigits->At(i) ; 
346     digit->SetIndexInList(i) ;     
347   }
348
349   gAlice->TreeD()->Fill() ;
350
351   gAlice->TreeD()->Write(0,TObject::kOverwrite) ;
352  
353 }
354
355 //___________________________________________________________________________
356 void AliPHOSv1::MakeBranch(Option_t* opt, char *file)
357
358   // Called by AliRun
359
360   char *cH ; 
361   // Create new branche in the current Root Tree in the digit Tree
362   AliDetector::MakeBranch(opt) ;
363
364
365   cH = strstr(opt,"S");
366   //Create a branch for SDigits
367   if( cH ){
368     char branchname[20];
369     sprintf(branchname,"%s",GetName());  
370     if(fSDigits)
371       fSDigits->Clear();
372
373     fnSdigits = 0 ;
374     gAlice->MakeBranchInTree(gAlice->TreeS(),branchname,&fSDigits,fBufferSize,file);  
375   }
376   
377   cH = strstr(opt,"D");
378   //Create a branch for Digits
379   if( cH ){
380     char branchname[20];
381     sprintf(branchname,"%s",GetName());  
382
383     if(fDigits)
384       fDigits->Clear();
385     
386     gAlice->MakeBranchInTree(gAlice->TreeD(),branchname,&fDigits,fBufferSize,file);  
387   }
388
389   cH = strstr(opt,"R");
390   //Create a branch for Reconstruction
391   if( cH ){
392     char branchname[20];
393
394     Int_t splitlevel = 0 ; 
395
396     if(fEmcRecPoints)
397       fEmcRecPoints->Delete() ; 
398
399     if ( fEmcRecPoints && gAlice->TreeR() ) {
400       sprintf(branchname,"%sEmcRP",GetName()) ;
401       gAlice->MakeBranchInTree(gAlice->TreeR(),branchname,"TObjArray",&fEmcRecPoints, fBufferSize, splitlevel,file); 
402     }
403
404     if(fPpsdRecPoints)
405       fPpsdRecPoints->Delete() ; 
406
407     if ( fPpsdRecPoints && gAlice->TreeR() ) {
408       sprintf(branchname,"%sPpsdRP",GetName()) ;
409       gAlice->MakeBranchInTree(gAlice->TreeR(),branchname,"TObjArray",&fPpsdRecPoints, fBufferSize, splitlevel,file); 
410     }
411
412     if(fTrackSegments)
413       fTrackSegments->Clear() ; 
414     
415     if ( fTrackSegments && gAlice->TreeR() ) { 
416       sprintf(branchname,"%sTS",GetName()) ;
417       gAlice->MakeBranchInTree(gAlice->TreeR(),branchname,&fTrackSegments,fBufferSize,file);
418     }
419     
420     if(fRecParticles)
421       fRecParticles->Clear() ; 
422     
423     if ( fRecParticles && gAlice->TreeR() ) { 
424       sprintf(branchname,"%sRP",GetName()) ;
425       gAlice->MakeBranchInTree(gAlice->TreeR(),branchname,&fRecParticles,fBufferSize,file); 
426     }
427     
428   }
429
430 }
431
432 //_____________________________________________________________________________
433 void AliPHOSv1::Reconstruction(AliPHOSReconstructor * Reconstructioner)
434
435   // 1. Reinitializes the existing RecPoint, TrackSegment, and RecParticles Lists and 
436   // 2. Creates TreeR with a branch for each list
437   // 3. Steers the reconstruction processes
438   // 4. Saves the 3 lists in TreeR
439   // 5. Write the Tree to File
440   
441   fReconstructioner = Reconstructioner ;
442     
443   // 1.
444
445   //  gAlice->MakeTree("R") ; 
446   
447   MakeBranch("R") ;
448   
449   // 3.
450
451   fReconstructioner->Make(fDigits, fEmcRecPoints, fPpsdRecPoints, fTrackSegments, fRecParticles);
452
453   printf("Reconstruction: %d %d %d %d\n",
454          fEmcRecPoints->GetEntries(),fPpsdRecPoints->GetEntries(),
455          fTrackSegments->GetEntries(),fRecParticles->GetEntries());
456
457   // 4. Expand or Shrink the arrays to the proper size
458   
459   Int_t size ;
460   
461   size = fEmcRecPoints->GetEntries() ;
462   fEmcRecPoints->Expand(size) ;
463
464   size = fPpsdRecPoints->GetEntries() ;
465   fPpsdRecPoints->Expand(size) ;
466
467   size = fTrackSegments->GetEntries() ;
468   fTrackSegments->Expand(size) ;
469
470   size = fRecParticles->GetEntries() ;
471   fRecParticles->Expand(size) ;
472
473   gAlice->TreeR()->Fill() ;
474   // 5.
475
476   gAlice->TreeR()->Write(0,TObject::kOverwrite) ;
477  
478   // Deleting reconstructed objects
479   ResetReconstruction();
480   
481 }
482
483 //____________________________________________________________________________
484 void AliPHOSv1::ResetReconstruction() 
485
486   // Deleting reconstructed objects
487
488   if ( fEmcRecPoints  )  fEmcRecPoints ->Delete();
489   if ( fPpsdRecPoints )  fPpsdRecPoints->Delete();
490   if ( fTrackSegments )  fTrackSegments->Delete();
491   if ( fRecParticles  )  fRecParticles ->Delete();
492   
493 }
494
495 //____________________________________________________________________________
496
497 void AliPHOSv1::StepManager(void)
498 {
499   // Accumulates hits as long as the track stays in a single crystal or PPSD gas Cell
500
501   Int_t          relid[4] ;           // (box, layer, row, column) indices
502   Int_t          absid    ;           // absolute cell ID number
503   Float_t        xyze[4]={0,0,0,0}  ; // position wrt MRS and energy deposited
504   TLorentzVector pos      ;           // Lorentz vector of the track current position
505   TLorentzVector pmom     ;        //momentum of the particle initiated hit
506   Float_t        xyd[3]={0,0,0}   ;   //local posiiton of the entering
507   Bool_t         entered = kFALSE ;  
508   Int_t          copy     ;
509
510   Int_t tracknumber =  gAlice->GetCurrentTrackNumber() ; 
511   Int_t primary     =  gAlice->GetPrimary( gAlice->GetCurrentTrackNumber() ); 
512   TString name      =  GetGeometry()->GetName() ; 
513   Int_t trackpid    =  0  ; 
514
515   if( gMC->IsTrackEntering() ){ // create hit with position and momentum of new particle, 
516                                 // but may be without energy deposition
517
518     // Current position of the hit in the local ref. system
519       gMC -> TrackPosition(pos);
520       Float_t xyzm[3], xyzd[3] ;
521       Int_t i;
522       for (i=0; i<3; i++) xyzm[i] = pos[i];
523       gMC -> Gmtod (xyzm, xyzd, 1);    // transform coordinate from master to daughter system
524       xyd[0]  = xyzd[0];
525       xyd[1]  =-xyzd[1];
526       xyd[2]  =-xyzd[2];
527
528       
529       // Current momentum of the hit's track in the local ref. system
530       gMC -> TrackMomentum(pmom);
531       Float_t pm[3], pd[3];
532       for (i=0; i<3; i++) pm[i]   = pmom[i];
533       gMC -> Gmtod (pm, pd, 2);        // transform 3-momentum from master to daughter system
534       pmom[0] = pd[0];
535       pmom[1] =-pd[1];
536       pmom[2] =-pd[2];
537
538       trackpid = gMC->TrackPid();
539       entered = kTRUE ;      // Mark to create hit even withou energy deposition
540
541   }
542
543
544   if ( name == "IHEP" ) {       // ======> CPV is a IHEP's one
545
546     // Yuri Kharlov, 28 September 2000
547
548     static Int_t idPCPQ = gMC->VolId("PCPQ");
549     if( gMC->CurrentVolID(copy) == idPCPQ &&
550         entered &&
551         gMC->TrackCharge() != 0) {      
552       
553       // Digitize the current CPV hit:
554
555       // 1. find pad response and
556       
557       Int_t moduleNumber;
558       gMC->CurrentVolOffID(3,moduleNumber);
559       moduleNumber--;
560
561
562       TClonesArray *cpvDigits = new TClonesArray("AliPHOSCPVDigit",0);   // array of digits for current hit
563       CPVDigitize(pmom,xyd,moduleNumber,cpvDigits);
564       
565       Float_t xmean = 0;
566       Float_t zmean = 0;
567       Float_t qsum  = 0;
568       Int_t   idigit,ndigits;
569
570       // 2. go through the current digit list and sum digits in pads
571
572       ndigits = cpvDigits->GetEntriesFast();
573       for (idigit=0; idigit<ndigits-1; idigit++) {
574         AliPHOSCPVDigit  *cpvDigit1 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit);
575         Float_t x1 = cpvDigit1->GetXpad() ;
576         Float_t z1 = cpvDigit1->GetYpad() ;
577         for (Int_t jdigit=idigit+1; jdigit<ndigits; jdigit++) {
578           AliPHOSCPVDigit  *cpvDigit2 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(jdigit);
579           Float_t x2 = cpvDigit2->GetXpad() ;
580           Float_t z2 = cpvDigit2->GetYpad() ;
581           if (x1==x2 && z1==z2) {
582             Float_t qsum = cpvDigit1->GetQpad() + cpvDigit2->GetQpad() ;
583             cpvDigit2->SetQpad(qsum) ;
584             cpvDigits->RemoveAt(idigit) ;
585           }
586         }
587       }
588       cpvDigits->Compress() ;
589
590       // 3. add digits to temporary hit list fTmpHits
591
592       ndigits = cpvDigits->GetEntriesFast();
593       for (idigit=0; idigit<ndigits; idigit++) {
594         AliPHOSCPVDigit  *cpvDigit = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit);
595         relid[0] = moduleNumber + 1 ;                             // CPV (or PHOS) module number
596         relid[1] =-1 ;                                            // means CPV
597         relid[2] = cpvDigit->GetXpad() ;                          // column number of a pad
598         relid[3] = cpvDigit->GetYpad() ;                          // row    number of a pad
599         
600         // get the absolute Id number
601         GetGeometry()->RelToAbsNumbering(relid, absid) ; 
602
603         // add current digit to the temporary hit list
604         xyze[0] = 0. ;
605         xyze[1] = 0. ;
606         xyze[2] = 0. ;
607         xyze[3] = cpvDigit->GetQpad() ;                           // amplitude in a pad
608         primary = -1;                                             // No need in primary for CPV
609         AddHit(fIshunt, primary, tracknumber, absid, xyze, trackpid, pmom, xyd);
610
611         if (cpvDigit->GetQpad() > 0.02) {
612           xmean += cpvDigit->GetQpad() * (cpvDigit->GetXpad() + 0.5);
613           zmean += cpvDigit->GetQpad() * (cpvDigit->GetYpad() + 0.5);
614           qsum  += cpvDigit->GetQpad();
615         }
616       }
617       delete cpvDigits;
618     }
619   } // end of IHEP configuration
620   
621
622   if(gMC->CurrentVolID(copy) == gMC->VolId("PXTL") ) { //  We are inside a PBWO crystal
623     gMC->TrackPosition(pos) ;
624     xyze[0] = pos[0] ;
625     xyze[1] = pos[1] ;
626     xyze[2] = pos[2] ;
627     xyze[3] = gMC->Edep() ;
628
629   
630     if ( (xyze[3] != 0) || entered ) {  // Track is inside the crystal and deposits some energy or just entered 
631
632       gMC->CurrentVolOffID(10, relid[0]) ; // get the PHOS module number ;
633
634       relid[1] = 0   ;                    // means PBW04
635       gMC->CurrentVolOffID(4, relid[2]) ; // get the row number inside the module
636       gMC->CurrentVolOffID(3, relid[3]) ; // get the cell number inside the module
637       
638       // get the absolute Id number
639       GetGeometry()->RelToAbsNumbering(relid, absid) ; 
640
641       // add current hit to the hit list
642         AddHit(fIshunt, primary,tracknumber, absid, xyze, trackpid,pmom, xyd);
643
644
645     } // there is deposited energy
646   } // we are inside a PHOS Xtal
647
648
649 }
650
651 //____________________________________________________________________________
652 void AliPHOSv1::CPVDigitize (TLorentzVector p, Float_t *zxhit, Int_t moduleNumber, TClonesArray *cpvDigits)
653 {
654   // ------------------------------------------------------------------------
655   // Digitize one CPV hit:
656   // On input take exact 4-momentum p and position zxhit of the hit,
657   // find the pad response around this hit and
658   // put the amplitudes in the pads into array digits
659   //
660   // Author: Yuri Kharlov (after Serguei Sadovsky)
661   // 2 October 2000
662   // ------------------------------------------------------------------------
663
664   const Float_t kCelWr  = GetGeometry()->GetPadSizePhi()/2;  // Distance between wires (2 wires above 1 pad)
665   const Float_t kDetR   = 0.1;     // Relative energy fluctuation in track for 100 e-
666   const Float_t kdEdx   = 4.0;     // Average energy loss in CPV;
667   const Int_t   kNgamz  = 5;       // Ionization size in Z
668   const Int_t   kNgamx  = 9;       // Ionization size in Phi
669   const Float_t kNoise = 0.03;    // charge noise in one pad
670
671   Float_t rnor1,rnor2;
672
673   // Just a reminder on axes notation in the CPV module:
674   // axis Z goes along the beam
675   // axis X goes across the beam in the module plane
676   // axis Y is a normal to the module plane showing from the IP
677
678   Float_t hitX  = zxhit[0];
679   Float_t hitZ  =-zxhit[1];
680   Float_t pX    = p.Px();
681   Float_t pZ    =-p.Pz();
682   Float_t pNorm = p.Py();
683   Float_t eloss = kdEdx;
684
685   Float_t dZY   = pZ/pNorm * GetGeometry()->GetCPVGasThickness();
686   Float_t dXY   = pX/pNorm * GetGeometry()->GetCPVGasThickness();
687   gRandom->Rannor(rnor1,rnor2);
688   eloss *= (1 + kDetR*rnor1) *
689            TMath::Sqrt((1 + ( pow(dZY,2) + pow(dXY,2) ) / pow(GetGeometry()->GetCPVGasThickness(),2)));
690   Float_t zhit1 = hitZ + GetGeometry()->GetCPVActiveSize(1)/2 - dZY/2;
691   Float_t xhit1 = hitX + GetGeometry()->GetCPVActiveSize(0)/2 - dXY/2;
692   Float_t zhit2 = zhit1 + dZY;
693   Float_t xhit2 = xhit1 + dXY;
694
695   Int_t   iwht1 = (Int_t) (xhit1 / kCelWr);           // wire (x) coordinate "in"
696   Int_t   iwht2 = (Int_t) (xhit2 / kCelWr);           // wire (x) coordinate "out"
697
698   Int_t   nIter;
699   Float_t zxe[3][5];
700   if (iwht1==iwht2) {                      // incline 1-wire hit
701     nIter = 2;
702     zxe[0][0] = (zhit1 + zhit2 - dZY*0.57735) / 2;
703     zxe[1][0] = (iwht1 + 0.5) * kCelWr;
704     zxe[2][0] =  eloss/2;
705     zxe[0][1] = (zhit1 + zhit2 + dZY*0.57735) / 2;
706     zxe[1][1] = (iwht1 + 0.5) * kCelWr;
707     zxe[2][1] =  eloss/2;
708   }
709   else if (TMath::Abs(iwht1-iwht2) != 1) { // incline 3-wire hit
710     nIter = 3;
711     Int_t iwht3 = (iwht1 + iwht2) / 2;
712     Float_t xwht1 = (iwht1 + 0.5) * kCelWr; // wire 1
713     Float_t xwht2 = (iwht2 + 0.5) * kCelWr; // wire 2
714     Float_t xwht3 = (iwht3 + 0.5) * kCelWr; // wire 3
715     Float_t xwr13 = (xwht1 + xwht3) / 2;   // center 13
716     Float_t xwr23 = (xwht2 + xwht3) / 2;   // center 23
717     Float_t dxw1  = xhit1 - xwr13;
718     Float_t dxw2  = xhit2 - xwr23;
719     Float_t egm1  = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
720     Float_t egm2  = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
721     Float_t egm3  =           kCelWr / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
722     zxe[0][0] = (dXY*(xwr13-xwht1)/dXY + zhit1 + zhit1) / 2;
723     zxe[1][0] =  xwht1;
724     zxe[2][0] =  eloss * egm1;
725     zxe[0][1] = (dXY*(xwr23-xwht1)/dXY + zhit1 + zhit2) / 2;
726     zxe[1][1] =  xwht2;
727     zxe[2][1] =  eloss * egm2;
728     zxe[0][2] =  dXY*(xwht3-xwht1)/dXY + zhit1;
729     zxe[1][2] =  xwht3;
730     zxe[2][2] =  eloss * egm3;
731   }
732   else {                                   // incline 2-wire hit
733     nIter = 2;
734     Float_t xwht1 = (iwht1 + 0.5) * kCelWr;
735     Float_t xwht2 = (iwht2 + 0.5) * kCelWr;
736     Float_t xwr12 = (xwht1 + xwht2) / 2;
737     Float_t dxw1  = xhit1 - xwr12;
738     Float_t dxw2  = xhit2 - xwr12;
739     Float_t egm1  = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) );
740     Float_t egm2  = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) );
741     zxe[0][0] = (zhit1 + zhit2 - dZY*egm1) / 2;
742     zxe[1][0] =  xwht1;
743     zxe[2][0] =  eloss * egm1;
744     zxe[0][1] = (zhit1 + zhit2 + dZY*egm2) / 2;
745     zxe[1][1] =  xwht2;
746     zxe[2][1] =  eloss * egm2;
747   }
748
749   // Finite size of ionization region
750
751   Int_t nCellZ  = GetGeometry()->GetNumberOfCPVPadsZ();
752   Int_t nCellX  = GetGeometry()->GetNumberOfCPVPadsPhi();
753   Int_t nz3     = (kNgamz+1)/2;
754   Int_t nx3     = (kNgamx+1)/2;
755   cpvDigits->Expand(nIter*kNgamx*kNgamz);
756   TClonesArray &ldigits = *(TClonesArray *)cpvDigits;
757
758   for (Int_t iter=0; iter<nIter; iter++) {
759
760     Float_t zhit = zxe[0][iter];
761     Float_t xhit = zxe[1][iter];
762     Float_t qhit = zxe[2][iter];
763     Float_t zcell = zhit / GetGeometry()->GetPadSizeZ();
764     Float_t xcell = xhit / GetGeometry()->GetPadSizePhi();
765     if ( zcell<=0      || xcell<=0 ||
766          zcell>=nCellZ || xcell>=nCellX) return;
767     Int_t izcell = (Int_t) zcell;
768     Int_t ixcell = (Int_t) xcell;
769     Float_t zc = zcell - izcell - 0.5;
770     Float_t xc = xcell - ixcell - 0.5;
771     for (Int_t iz=1; iz<=kNgamz; iz++) {
772       Int_t kzg = izcell + iz - nz3;
773       if (kzg<=0 || kzg>nCellZ) continue;
774       Float_t zg = (Float_t)(iz-nz3) - zc;
775       for (Int_t ix=1; ix<=kNgamx; ix++) {
776         Int_t kxg = ixcell + ix - nx3;
777         if (kxg<=0 || kxg>nCellX) continue;
778         Float_t xg = (Float_t)(ix-nx3) - xc;
779         
780         // Now calculate pad response
781         Float_t qpad = CPVPadResponseFunction(qhit,zg,xg);
782         qpad += kNoise*rnor2;
783         if (qpad<0) continue;
784         
785         // Fill the array with pad response ID and amplitude
786         new(ldigits[cpvDigits->GetEntriesFast()]) AliPHOSCPVDigit(kxg,kzg,qpad);
787       }
788     }
789   }
790 }
791
792 //____________________________________________________________________________
793 Float_t AliPHOSv1::CPVPadResponseFunction(Float_t qhit, Float_t zhit, Float_t xhit) {
794   // ------------------------------------------------------------------------
795   // Calculate the amplitude in one CPV pad using the
796   // cumulative pad response function
797   // Author: Yuri Kharlov (after Serguei Sadovski)
798   // 3 October 2000
799   // ------------------------------------------------------------------------
800
801   Double_t dz = GetGeometry()->GetPadSizeZ()   / 2;
802   Double_t dx = GetGeometry()->GetPadSizePhi() / 2;
803   Double_t z  = zhit * GetGeometry()->GetPadSizeZ();
804   Double_t x  = xhit * GetGeometry()->GetPadSizePhi();
805   Double_t amplitude = qhit *
806     (CPVCumulPadResponse(z+dz,x+dx) - CPVCumulPadResponse(z+dz,x-dx) -
807      CPVCumulPadResponse(z-dz,x+dx) + CPVCumulPadResponse(z-dz,x-dx));
808   return (Float_t)amplitude;
809 }
810
811 //____________________________________________________________________________
812 Double_t AliPHOSv1::CPVCumulPadResponse(Double_t x, Double_t y) {
813   // ------------------------------------------------------------------------
814   // Cumulative pad response function
815   // It includes several terms from the CF decomposition in electrostatics
816   // Note: this cumulative function is wrong since omits some terms
817   //       but the cell amplitude obtained with it is correct because
818   //       these omitting terms cancel
819   // Author: Yuri Kharlov (after Serguei Sadovski)
820   // 3 October 2000
821   // ------------------------------------------------------------------------
822
823   const Double_t kA=1.0;
824   const Double_t kB=0.7;
825
826   Double_t r2       = x*x + y*y;
827   Double_t xy       = x*y;
828   Double_t cumulPRF = 0;
829   for (Int_t i=0; i<=4; i++) {
830     Double_t b1 = (2*i + 1) * kB;
831     cumulPRF += TMath::Power(-1,i) * TMath::ATan( xy / (b1*TMath::Sqrt(b1*b1 + r2)) );
832   }
833   cumulPRF *= kA/(2*TMath::Pi());
834   return cumulPRF;
835 }
836