]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSDebug.cxx
5e0df8e245452636821b4a1c38a765473c749d30
[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 + PPSD has name GPS2:
22 // Produces cumulated hits
23 //---
24 // Layout EMC + CPV  has name IHEP:
25 // Produces hits for CPV, cumulated hits
26 //---
27 // Layout EMC + CPV + PPSD has name GPS:
28 // Produces hits for CPV, cumulated hits
29 //---
30 //*-- Author: Yves Schutz (SUBATECH)
31
32
33 // --- ROOT system ---
34
35 #include "TBRIK.h"
36 #include "TNode.h"
37 #include "TRandom.h"
38 #include "TTree.h"
39
40
41 // --- Standard library ---
42
43 #include <stdio.h>
44 #include <string.h>
45 #include <stdlib.h>
46 #include <strstream.h>
47
48 // --- AliRoot header files ---
49
50 #include "AliPHOSv1.h"
51 #include "AliPHOSHit.h"
52 #include "AliPHOSDigit.h"
53 #include "AliPHOSReconstructioner.h"
54 #include "AliRun.h"
55 #include "AliConst.h"
56 #include "AliMC.h"
57
58 ClassImp(AliPHOSv1)
59
60 //____________________________________________________________________________
61 AliPHOSv1::AliPHOSv1():
62 AliPHOSv0()
63 {
64   // ctor
65  
66   fReconstructioner  = 0;
67   fTrackSegmentMaker = 0;
68
69 }
70
71 //____________________________________________________________________________
72 AliPHOSv1::AliPHOSv1(const char *name, const char *title):
73 AliPHOSv0(name,title) 
74 {
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)
79   //
80   // We store hits :
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).
84   //
85
86   fPinElectronicNoise = 0.010 ;
87   fDigitThreshold      = 0.01 ;   // 1 GeV 
88   fDigitizeA= 0. ;             
89   fDigitizeB = 10000000. ;    
90
91
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() ).
96   
97   fHits= new TClonesArray("AliPHOSHit",1000) ;
98
99   fNhits = 0 ;
100
101   fReconstructioner  = 0;
102   fTrackSegmentMaker = 0;
103
104   fIshunt     =  1 ; // All hits are associated with primary particles
105   
106 }
107
108 //____________________________________________________________________________
109 AliPHOSv1::AliPHOSv1(AliPHOSReconstructioner * Reconstructioner, const char *name, const char *title):
110   AliPHOSv0(name,title)
111 {
112   // ctor : title is used to identify the layout
113   //        GPS2 = 5 modules (EMC + PPSD)   
114
115   fPinElectronicNoise = 0.010 ;
116
117   // We do not want to save in TreeH the raw hits
118
119   fDigits = 0 ;
120   fHits= new TClonesArray("AliPHOSHit",1000) ;
121
122   fNhits = 0 ;
123
124   fIshunt     =  1 ; // All hits are associated with primary particles
125  
126   // gets an instance of the geometry parameters class  
127   fGeom =  AliPHOSGeometry::GetInstance(title, "") ; 
128
129   if (fGeom->IsInitialized() ) 
130     cout << "AliPHOS" << Version() << " : PHOS geometry intialized for " << fGeom->GetName() << endl ;
131   else
132     cout << "AliPHOS" << Version() << " : PHOS geometry initialization failed !" << endl ;   
133
134   // Defining the PHOS Reconstructioner
135  
136  fReconstructioner = Reconstructioner ;
137
138 }
139
140 //____________________________________________________________________________
141 AliPHOSv1::~AliPHOSv1()
142 {
143   // dtor
144
145   if ( fHits) {
146     fHits->Delete() ; 
147     delete fHits ;
148     fHits = 0 ; 
149   }
150
151   if ( fSDigits) {
152     fSDigits->Delete() ; 
153     delete fSDigits ;
154     fSDigits = 0 ; 
155   }
156
157   if ( fDigits) {
158     fDigits->Delete() ; 
159     delete fDigits ;
160     fDigits = 0 ; 
161   }
162
163   if ( fEmcRecPoints ) {
164     fEmcRecPoints->Delete() ; 
165     delete fEmcRecPoints ; 
166     fEmcRecPoints = 0 ; 
167   }
168   
169   if ( fPpsdRecPoints ) { 
170     fPpsdRecPoints->Delete() ;
171     delete fPpsdRecPoints ;
172     fPpsdRecPoints = 0 ; 
173   }
174     
175   if ( fTrackSegments ) {
176     fTrackSegments->Delete() ; 
177     delete fTrackSegments ;
178     fTrackSegments = 0 ; 
179   }
180   
181 }
182
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)
185 {
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
189
190   Int_t hitCounter ;
191   AliPHOSHit *newHit ;
192   AliPHOSHit *curHit ;
193   Bool_t deja = kFALSE ;
194
195   newHit = new AliPHOSHit(shunt, primary, tracknumber, Id, hits, trackpid, p, lpos) ;
196
197   for ( hitCounter = fNhits-1 ; hitCounter >= 0 && !deja ; hitCounter-- ) {
198     curHit = (AliPHOSHit*) (*fHits)[hitCounter] ;
199     if( *curHit == *newHit ) {
200       *curHit = *curHit + *newHit ;
201       deja = kTRUE ;
202     }
203   }
204          
205   if ( !deja ) {
206     new((*fHits)[fNhits]) AliPHOSHit(*newHit) ;
207     fNhits++ ;
208   }
209
210   delete newHit;
211 }
212
213 //____________________________________________________________________________
214 void AliPHOSv1::Hits2SDigits(){
215   //Collects all hits in the same active volume into digit
216
217   Int_t i ;
218   Int_t j ; 
219   AliPHOSHit   * hit ;
220   AliPHOSDigit * newdigit ;
221   AliPHOSDigit * curdigit ;
222   Bool_t deja = kFALSE ; 
223   
224
225   Int_t itrack ;
226   for (itrack=0; itrack<gAlice->GetNtrack(); itrack++){
227         
228     //=========== Get the Hits Tree for the Primary track itrack
229     gAlice->ResetHits();    
230     gAlice->TreeH()->GetEvent(itrack);
231       
232
233     for ( i = 0 ; i < fHits->GetEntries() ; i++ ) {
234       hit = (AliPHOSHit*)fHits->At(i) ;
235     
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() ) ) ;
239       else
240         newdigit = new AliPHOSDigit( -1               , hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
241       deja =kFALSE ;
242
243       for ( j = 0 ; j < fnSdigits ;  j++) { 
244         curdigit = (AliPHOSDigit*) fSDigits->At(j) ;
245         if ( *curdigit == *newdigit) {
246           *curdigit = *curdigit + *newdigit ; 
247           deja = kTRUE ; 
248         }
249       }
250
251       if ( !deja ) {
252         new((*fSDigits)[fnSdigits]) AliPHOSDigit(* newdigit) ;
253         fnSdigits++ ;  
254       }
255       
256       delete newdigit ;    
257     } 
258
259   } // loop over tracks
260
261   fSDigits->Sort() ;
262
263   fnSdigits = fSDigits->GetEntries() ;
264   fSDigits->Expand(fnSdigits) ;
265
266   for (i = 0 ; i < fnSdigits ; i++) { 
267     AliPHOSDigit * digit = (AliPHOSDigit *) fSDigits->At(i) ; 
268     digit->SetIndexInList(i) ;     
269   }
270
271   gAlice->TreeS()->Fill() ;
272   gAlice->TreeS()->Write(0,TObject::kOverwrite) ;
273  
274
275 }
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.
280
281
282   gAlice->TreeS()->GetEvent(0) ;
283
284   // First calculate noise induced by the PIN diode of the PbWO crystals
285   Int_t iCurSDigit = 0 ;
286
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 ;
290     return ;
291   }
292
293   Int_t idCurSDigit = ((AliPHOSDigit *)fSDigits->At(0))->GetId() ;
294
295   Int_t absID ;
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) ) ;
301         fNdigits++ ;  
302       }
303     }
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));
310         fNdigits++ ;  
311       }
312
313       if(iCurSDigit < fSDigits->GetEntries()-1){
314         iCurSDigit++ ;
315         idCurSDigit = ((AliPHOSDigit*)fSDigits->At(iCurSDigit))->GetId() ;
316       }
317       else
318         idCurSDigit = 10000000; //no real hits left    
319     }
320     
321   }  
322
323   //remove PPSD/CPV digits below thresholds
324   Int_t idigit ;
325   for(idigit = iCurSDigit; idigit < fSDigits->GetEntries() ; idigit++){  //loop over CPV/PPSD digits
326     
327     AliPHOSDigit * digit = (AliPHOSDigit *) fSDigits->At(idigit) ; 
328     Float_t ene = Calibrate(digit->GetAmp()) ;
329     
330     Int_t relid[4] ; 
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 ) ;
336         fNdigits++ ;
337     }
338   }    
339     
340   fDigits->Compress() ;  
341   
342   fNdigits = fDigits->GetEntries() ;
343   fDigits->Expand(fNdigits) ;
344
345   Int_t i ;
346   for (i = 0 ; i < fNdigits ; i++) { 
347     AliPHOSDigit * digit = (AliPHOSDigit *) fDigits->At(i) ; 
348     digit->SetIndexInList(i) ;     
349   }
350
351   gAlice->TreeD()->Fill() ;
352
353   gAlice->TreeD()->Write(0,TObject::kOverwrite) ;
354  
355 }
356
357 //___________________________________________________________________________
358 void AliPHOSv1::MakeBranch(Option_t* opt, char *file)
359
360
361   char *cH ; 
362   // Create new branche in the current Root Tree in the digit Tree
363   AliDetector::MakeBranch(opt) ;
364
365
366   cH = strstr(opt,"S");
367   //Create a branch for SDigits
368   if( cH ){
369     char branchname[20];
370     sprintf(branchname,"%s",GetName());  
371     if(fSDigits)
372       fSDigits->Clear();
373
374     fnSdigits = 0 ;
375     gAlice->MakeBranchInTree(gAlice->TreeS(),branchname,&fSDigits,fBufferSize,file);  
376   }
377   
378   cH = strstr(opt,"D");
379   //Create a branch for Digits
380   if( cH ){
381     char branchname[20];
382     sprintf(branchname,"%s",GetName());  
383
384     if(fDigits)
385       fDigits->Clear();
386     
387     gAlice->MakeBranchInTree(gAlice->TreeD(),branchname,&fDigits,fBufferSize,file);  
388   }
389
390   cH = strstr(opt,"R");
391   //Create a branch for Reconstruction
392   if( cH ){
393     char branchname[20];
394
395     Int_t splitlevel = 0 ; 
396
397     if(fEmcRecPoints)
398       fEmcRecPoints->Delete() ; 
399
400     if ( fEmcRecPoints && gAlice->TreeR() ) {
401       sprintf(branchname,"%sEmcRP",GetName()) ;
402       gAlice->MakeBranchInTree(gAlice->TreeR(),branchname,"TObjArray",&fEmcRecPoints, fBufferSize, splitlevel,file); 
403     }
404
405     if(fPpsdRecPoints)
406       fPpsdRecPoints->Delete() ; 
407
408     if ( fPpsdRecPoints && gAlice->TreeR() ) {
409       sprintf(branchname,"%sPpsdRP",GetName()) ;
410       gAlice->MakeBranchInTree(gAlice->TreeR(),branchname,"TObjArray",&fPpsdRecPoints, fBufferSize, splitlevel,file); 
411     }
412
413     if(fTrackSegments)
414       fTrackSegments->Clear() ; 
415     
416     if ( fTrackSegments && gAlice->TreeR() ) { 
417       sprintf(branchname,"%sTS",GetName()) ;
418       gAlice->MakeBranchInTree(gAlice->TreeR(),branchname,&fTrackSegments,fBufferSize,file);
419     }
420     
421     if(fRecParticles)
422       fRecParticles->Clear() ; 
423     
424     if ( fRecParticles && gAlice->TreeR() ) { 
425       sprintf(branchname,"%sRP",GetName()) ;
426       gAlice->MakeBranchInTree(gAlice->TreeR(),branchname,&fRecParticles,fBufferSize,file); 
427     }
428     
429   }
430
431 }
432
433 //_____________________________________________________________________________
434 void AliPHOSv1::Reconstruction(AliPHOSReconstructioner * Reconstructioner)
435
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
441   
442   fReconstructioner = Reconstructioner ;
443     
444   // 1.
445
446   //  gAlice->MakeTree("R") ; 
447   
448   MakeBranch("R") ;
449   
450   // 3.
451
452   fReconstructioner->Make(fDigits, fEmcRecPoints, fPpsdRecPoints, fTrackSegments, fRecParticles);
453
454   printf("Reconstruction: %d %d %d %d\n",
455          fEmcRecPoints->GetEntries(),fPpsdRecPoints->GetEntries(),
456          fTrackSegments->GetEntries(),fRecParticles->GetEntries());
457
458   // 4. Expand or Shrink the arrays to the proper size
459   
460   Int_t size ;
461   
462   size = fEmcRecPoints->GetEntries() ;
463   fEmcRecPoints->Expand(size) ;
464
465   size = fPpsdRecPoints->GetEntries() ;
466   fPpsdRecPoints->Expand(size) ;
467
468   size = fTrackSegments->GetEntries() ;
469   fTrackSegments->Expand(size) ;
470
471   size = fRecParticles->GetEntries() ;
472   fRecParticles->Expand(size) ;
473
474   gAlice->TreeR()->Fill() ;
475   // 5.
476
477   gAlice->TreeR()->Write(0,TObject::kOverwrite) ;
478  
479   // Deleting reconstructed objects
480   ResetReconstruction();
481   
482 }
483
484 //____________________________________________________________________________
485 void AliPHOSv1::ResetReconstruction() 
486
487   // Deleting reconstructed objects
488
489   if ( fEmcRecPoints  )  fEmcRecPoints ->Delete();
490   if ( fPpsdRecPoints )  fPpsdRecPoints->Delete();
491   if ( fTrackSegments )  fTrackSegments->Delete();
492   if ( fRecParticles  )  fRecParticles ->Delete();
493   
494 }
495
496 //____________________________________________________________________________
497
498 void AliPHOSv1::StepManager(void)
499 {
500   // Accumulates hits as long as the track stays in a single crystal or PPSD gas Cell
501
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 ;  
509   Int_t          copy     ;
510
511   Int_t tracknumber =  gAlice->CurrentTrack() ; 
512   Int_t primary     =  gAlice->GetPrimary( gAlice->CurrentTrack() ); 
513   TString name      =  fGeom->GetName() ; 
514   Int_t trackpid    =  0  ; 
515
516   if( gMC->IsTrackEntering() ){ // create hit with position and momentum of new particle, 
517                                 // but may be without energy deposition
518
519     // Current position of the hit in the local ref. system
520       gMC -> TrackPosition(pos);
521       Float_t xyzm[3], xyzd[3] ;
522       Int_t i;
523       for (i=0; i<3; i++) xyzm[i] = pos[i];
524       gMC -> Gmtod (xyzm, xyzd, 1);    // transform coordinate from master to daughter system
525       xyd[0]  = xyzd[0];
526       xyd[1]  =-xyzd[1];
527       xyd[2]  =-xyzd[2];
528
529       
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
535       pmom[0] = pd[0];
536       pmom[1] =-pd[1];
537       pmom[2] =-pd[2];
538
539       trackpid = gMC->TrackPid();
540       entered = kTRUE ;      // Mark to create hit even withou energy deposition
541
542   }
543
544
545   if ( name == "GPS2" || name == "MIXT" ) {            // ======> CPV is a GPS' PPSD
546
547     if( gMC->CurrentVolID(copy) == gMC->VolId("GCEL") ) // We are inside a gas cell 
548     {
549       gMC->TrackPosition(pos) ;
550       xyze[0] = pos[0] ;
551       xyze[1] = pos[1] ;
552       xyze[2] = pos[2] ;
553       xyze[3] = gMC->Edep() ; 
554
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();
559         }
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 
565
566         // get the absolute Id number
567
568         fGeom->RelToAbsNumbering(relid, absid) ; 
569
570         // add current hit to the hit list      
571           AddHit(fIshunt, primary, tracknumber, absid, xyze, trackpid, pmom, xyd);
572
573
574       } // there is deposited energy 
575     } // We are inside the gas of the CPV  
576   } // GPS2 configuration
577
578   if ( name == "IHEP" || name == "MIXT" ) {       // ======> CPV is a IHEP's one
579
580     // Yuri Kharlov, 28 September 2000
581
582     if( gMC->CurrentVolID(copy) == gMC->VolId("CPVQ") &&
583         entered &&
584         gMC->TrackCharge() != 0) {      
585       
586       // Digitize the current CPV hit:
587
588       // 1. find pad response and
589       
590       Int_t moduleNumber;
591       gMC->CurrentVolOffID(3,moduleNumber);
592       moduleNumber--;
593
594
595       TClonesArray *cpvDigits = new TClonesArray("AliPHOSCPVDigit",0);   // array of digits for current hit
596       CPVDigitize(pmom,xyd,moduleNumber,cpvDigits);
597       
598       Float_t xmean = 0;
599       Float_t zmean = 0;
600       Float_t qsum  = 0;
601       Int_t   idigit,ndigits;
602
603       // 2. go through the current digit list and sum digits in pads
604
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) ;
618           }
619         }
620       }
621       cpvDigits->Compress() ;
622
623       // 3. add digits to temporary hit list fTmpHits
624
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
632         
633         // get the absolute Id number
634         fGeom->RelToAbsNumbering(relid, absid) ; 
635
636         // add current digit to the temporary hit list
637         xyze[0] = 0. ;
638         xyze[1] = 0. ;
639         xyze[2] = 0. ;
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);
643
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();
648         }
649       }
650       delete cpvDigits;
651     }
652   } // end of IHEP configuration
653   
654
655   if(gMC->CurrentVolID(copy) == gMC->VolId("PXTL") ) { //  We are inside a PBWO crystal
656     gMC->TrackPosition(pos) ;
657     xyze[0] = pos[0] ;
658     xyze[1] = pos[1] ;
659     xyze[2] = pos[2] ;
660     xyze[3] = gMC->Edep() ;
661
662   
663     if ( (xyze[3] != 0) || entered ) {  // Track is inside the crystal and deposits some energy or just entered 
664
665       gMC->CurrentVolOffID(10, relid[0]) ; // get the PHOS module number ;
666
667       if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(10),"PHO1") == 0 )
668         relid[0] += fGeom->GetNModules() - fGeom->GetNPPSDModules();      
669
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
673       
674       // get the absolute Id number
675       fGeom->RelToAbsNumbering(relid, absid) ; 
676
677       // add current hit to the hit list
678         AddHit(fIshunt, primary,tracknumber, absid, xyze, trackpid,pmom, xyd);
679
680
681     } // there is deposited energy
682   } // we are inside a PHOS Xtal
683
684
685 }
686
687 //____________________________________________________________________________
688 void AliPHOSv1::CPVDigitize (TLorentzVector p, Float_t *zxhit, Int_t moduleNumber, TClonesArray *cpvDigits)
689 {
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
695   //
696   // Author: Yuri Kharlov (after Serguei Sadovsky)
697   // 2 October 2000
698   // ------------------------------------------------------------------------
699
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
706
707   Float_t rnor1,rnor2;
708
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
713
714   Float_t hitX  = zxhit[0];
715   Float_t hitZ  =-zxhit[1];
716   Float_t pX    = p.Px();
717   Float_t pZ    =-p.Pz();
718   Float_t pNorm = p.Py();
719   Float_t eloss = kdEdx;
720
721 //    cout << "CPVDigitize: YVK : "<<hitX<<" "<<hitZ<<" | "<<pX<<" "<<pZ<<" "<<pNorm<<endl;
722
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;
732
733   Int_t   iwht1 = (Int_t) (xhit1 / kCelWr);           // wire (x) coordinate "in"
734   Int_t   iwht2 = (Int_t) (xhit2 / kCelWr);           // wire (x) coordinate "out"
735
736   Int_t   nIter;
737   Float_t zxe[3][5];
738   if (iwht1==iwht2) {                      // incline 1-wire hit
739     nIter = 2;
740     zxe[0][0] = (zhit1 + zhit2 - dZY*0.57735) / 2;
741     zxe[1][0] = (iwht1 + 0.5) * kCelWr;
742     zxe[2][0] =  eloss/2;
743     zxe[0][1] = (zhit1 + zhit2 + dZY*0.57735) / 2;
744     zxe[1][1] = (iwht1 + 0.5) * kCelWr;
745     zxe[2][1] =  eloss/2;
746   }
747   else if (TMath::Abs(iwht1-iwht2) != 1) { // incline 3-wire hit
748     nIter = 3;
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;
761     zxe[1][0] =  xwht1;
762     zxe[2][0] =  eloss * egm1;
763     zxe[0][1] = (dXY*(xwr23-xwht1)/dXY + zhit1 + zhit2) / 2;
764     zxe[1][1] =  xwht2;
765     zxe[2][1] =  eloss * egm2;
766     zxe[0][2] =  dXY*(xwht3-xwht1)/dXY + zhit1;
767     zxe[1][2] =  xwht3;
768     zxe[2][2] =  eloss * egm3;
769   }
770   else {                                   // incline 2-wire hit
771     nIter = 2;
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;
780     zxe[1][0] =  xwht1;
781     zxe[2][0] =  eloss * egm1;
782     zxe[0][1] = (zhit1 + zhit2 + dZY*egm2) / 2;
783     zxe[1][1] =  xwht2;
784     zxe[2][1] =  eloss * egm2;
785   }
786
787   // Finite size of ionization region
788
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;
795
796   for (Int_t iter=0; iter<nIter; iter++) {
797
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;
817         
818         // Now calculate pad response
819         Float_t qpad = CPVPadResponseFunction(qhit,zg,xg);
820         qpad += kNoise*rnor2;
821         if (qpad<0) continue;
822         
823         // Fill the array with pad response ID and amplitude
824         new(ldigits[cpvDigits->GetEntriesFast()]) AliPHOSCPVDigit(kxg,kzg,qpad);
825       }
826     }
827   }
828 }
829
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)
836   // 3 October 2000
837   // ------------------------------------------------------------------------
838
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;
847 }
848
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)
858   // 3 October 2000
859   // ------------------------------------------------------------------------
860
861   const Double_t kA=1.0;
862   const Double_t kB=0.7;
863
864   Double_t r2       = x*x + y*y;
865   Double_t xy       = x*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)) );
870   }
871   cumulPRF *= kA/(2*TMath::Pi());
872   return cumulPRF;
873 }
874