]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSDebug.cxx
Adding event generator for e+e- pair production
[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   AliPHOSGeometry::GetInstance(title, "") ; 
128
129   if (GetGeometry()->IsInitialized() ) 
130     Info("AliPHOSv1", "AliPHOS %d : PHOS geometry intialized for %s", Version(), GetGeometry()->GetName() );
131   else
132     Info("AliPHOSv1", "AliPHOS %d : PHOS geometry initialization failed !", Version() ) ;   
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 {
216   // Collects all hits in the same active volume into digit
217   // OBSOLETE replace by SDigitizer
218
219   Int_t i ;
220   Int_t j ; 
221   AliPHOSHit   * hit ;
222   AliPHOSDigit * newdigit ;
223   AliPHOSDigit * curdigit ;
224   Bool_t deja = kFALSE ; 
225   
226
227   Int_t itrack ;
228   for (itrack=0; itrack<gAlice->GetNtrack(); itrack++){
229         
230     //=========== Get the Hits Tree for the Primary track itrack
231     gAlice->ResetHits();    
232     gAlice->TreeH()->GetEvent(itrack);
233       
234
235     for ( i = 0 ; i < fHits->GetEntries() ; i++ ) {
236       hit = (AliPHOSHit*)fHits->At(i) ;
237     
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() ) ) ;
241       else
242         newdigit = new AliPHOSDigit( -1               , hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
243       deja =kFALSE ;
244
245       for ( j = 0 ; j < fnSdigits ;  j++) { 
246         curdigit = (AliPHOSDigit*) fSDigits->At(j) ;
247         if ( *curdigit == *newdigit) {
248           *curdigit = *curdigit + *newdigit ; 
249           deja = kTRUE ; 
250         }
251       }
252
253       if ( !deja ) {
254         new((*fSDigits)[fnSdigits]) AliPHOSDigit(* newdigit) ;
255         fnSdigits++ ;  
256       }
257       
258       delete newdigit ;    
259     } 
260
261   } // loop over tracks
262
263   fSDigits->Sort() ;
264
265   fnSdigits = fSDigits->GetEntries() ;
266   fSDigits->Expand(fnSdigits) ;
267
268   for (i = 0 ; i < fnSdigits ; i++) { 
269     AliPHOSDigit * digit = (AliPHOSDigit *) fSDigits->At(i) ; 
270     digit->SetIndexInList(i) ;     
271   }
272
273   gAlice->TreeS()->Fill() ;
274   gAlice->TreeS()->Write(0,TObject::kOverwrite) ;
275  
276
277 }
278 //____________________________________________________________________________
279 void AliPHOSv1::SDigits2Digits()
280 {
281   // Adds noise to the summable digits and removes everething below thresholds
282   // Note, that sDigits should be SORTED in accordance with abs ID.
283   // OBSOLETE Replaced by Digitzer
284
285   gAlice->TreeS()->GetEvent(0) ;
286
287   // First calculate noise induced by the PIN diode of the PbWO crystals
288   Int_t iCurSDigit = 0 ;
289
290   //we assume, that there is al least one EMC digit...
291   if(fSDigits->GetEntries() == 0) {
292     Warning("SDigits2Digits", "No SDigits !!! Do not produce Digits ") ;
293     return ;
294   }
295
296   Int_t idCurSDigit = ((AliPHOSDigit *)fSDigits->At(0))->GetId() ;
297
298   Int_t absID ;
299   for(absID = 1; absID < GetGeometry()->GetNModules()*GetGeometry()->GetNPhi()*GetGeometry()->GetNZ(); absID++){
300     Float_t noise = gRandom->Gaus(0., fPinElectronicNoise) ; 
301     if(absID < idCurSDigit ){ 
302       if(noise >fDigitThreshold ){
303         new((*fDigits)[fNdigits]) AliPHOSDigit( -1,absID,Digitize(noise) ) ;
304         fNdigits++ ;  
305       }
306     }
307     else{ //add noise and may be remove the true hit
308       Float_t signal = noise + Calibrate(((AliPHOSDigit *)fSDigits->At(iCurSDigit))->GetAmp()) ;
309       if( signal >fDigitThreshold ){
310         AliPHOSDigit * digit = (AliPHOSDigit*) fSDigits->At(iCurSDigit) ;
311         new((*fDigits)[fNdigits]) AliPHOSDigit( *digit ) ;
312         ((AliPHOSDigit *)fDigits->At(fNdigits))->SetAmp(Digitize(signal));
313         fNdigits++ ;  
314       }
315
316       if(iCurSDigit < fSDigits->GetEntries()-1){
317         iCurSDigit++ ;
318         idCurSDigit = ((AliPHOSDigit*)fSDigits->At(iCurSDigit))->GetId() ;
319       }
320       else
321         idCurSDigit = 10000000; //no real hits left    
322     }
323     
324   }  
325
326   //remove PPSD/CPV digits below thresholds
327   Int_t idigit ;
328   for(idigit = iCurSDigit; idigit < fSDigits->GetEntries() ; idigit++){  //loop over CPV/PPSD digits
329     
330     AliPHOSDigit * digit = (AliPHOSDigit *) fSDigits->At(idigit) ; 
331     Float_t ene = Calibrate(digit->GetAmp()) ;
332     
333     Int_t relid[4] ; 
334     GetGeometry()->AbsToRelNumbering(digit->GetId(), relid) ; 
335     if ( relid[0] > GetGeometry()->GetNCPVModules() ){ //ppsd
336       if ( ( (relid[1] > 0) && (ene > fPpsdEnergyThreshold)) ||    //PPSD digit
337            ( (relid[1] < 0) && (ene > fCpvEnergyThreshold ) ) )    //CPV digit 
338         new((*fDigits)[fNdigits]) AliPHOSDigit( *digit ) ;
339         fNdigits++ ;
340     }
341   }    
342     
343   fDigits->Compress() ;  
344   
345   fNdigits = fDigits->GetEntries() ;
346   fDigits->Expand(fNdigits) ;
347
348   Int_t i ;
349   for (i = 0 ; i < fNdigits ; i++) { 
350     AliPHOSDigit * digit = (AliPHOSDigit *) fDigits->At(i) ; 
351     digit->SetIndexInList(i) ;     
352   }
353
354   gAlice->TreeD()->Fill() ;
355
356   gAlice->TreeD()->Write(0,TObject::kOverwrite) ;
357  
358 }
359
360 //___________________________________________________________________________
361 void AliPHOSv1::MakeBranch(Option_t* opt, char *file)
362
363   // Called by AliRun
364
365   char *cH ; 
366   // Create new branche in the current Root Tree in the digit Tree
367   AliDetector::MakeBranch(opt) ;
368
369
370   cH = strstr(opt,"S");
371   //Create a branch for SDigits
372   if( cH ){
373     char branchname[20];
374     sprintf(branchname,"%s",GetName());  
375     if(fSDigits)
376       fSDigits->Clear();
377
378     fnSdigits = 0 ;
379     gAlice->MakeBranchInTree(gAlice->TreeS(),branchname,&fSDigits,fBufferSize,file);  
380   }
381   
382   cH = strstr(opt,"D");
383   //Create a branch for Digits
384   if( cH ){
385     char branchname[20];
386     sprintf(branchname,"%s",GetName());  
387
388     if(fDigits)
389       fDigits->Clear();
390     
391     gAlice->MakeBranchInTree(gAlice->TreeD(),branchname,&fDigits,fBufferSize,file);  
392   }
393
394   cH = strstr(opt,"R");
395   //Create a branch for Reconstruction
396   if( cH ){
397     char branchname[20];
398
399     Int_t splitlevel = 0 ; 
400
401     if(fEmcRecPoints)
402       fEmcRecPoints->Delete() ; 
403
404     if ( fEmcRecPoints && gAlice->TreeR() ) {
405       sprintf(branchname,"%sEmcRP",GetName()) ;
406       gAlice->MakeBranchInTree(gAlice->TreeR(),branchname,"TObjArray",&fEmcRecPoints, fBufferSize, splitlevel,file); 
407     }
408
409     if(fPpsdRecPoints)
410       fPpsdRecPoints->Delete() ; 
411
412     if ( fPpsdRecPoints && gAlice->TreeR() ) {
413       sprintf(branchname,"%sPpsdRP",GetName()) ;
414       gAlice->MakeBranchInTree(gAlice->TreeR(),branchname,"TObjArray",&fPpsdRecPoints, fBufferSize, splitlevel,file); 
415     }
416
417     if(fTrackSegments)
418       fTrackSegments->Clear() ; 
419     
420     if ( fTrackSegments && gAlice->TreeR() ) { 
421       sprintf(branchname,"%sTS",GetName()) ;
422       gAlice->MakeBranchInTree(gAlice->TreeR(),branchname,&fTrackSegments,fBufferSize,file);
423     }
424     
425     if(fRecParticles)
426       fRecParticles->Clear() ; 
427     
428     if ( fRecParticles && gAlice->TreeR() ) { 
429       sprintf(branchname,"%sRP",GetName()) ;
430       gAlice->MakeBranchInTree(gAlice->TreeR(),branchname,&fRecParticles,fBufferSize,file); 
431     }
432     
433   }
434
435 }
436
437 //_____________________________________________________________________________
438 void AliPHOSv1::Reconstruction(AliPHOSReconstructioner * Reconstructioner)
439
440   // 1. Reinitializes the existing RecPoint, TrackSegment, and RecParticles Lists and 
441   // 2. Creates TreeR with a branch for each list
442   // 3. Steers the reconstruction processes
443   // 4. Saves the 3 lists in TreeR
444   // 5. Write the Tree to File
445   
446   fReconstructioner = Reconstructioner ;
447     
448   // 1.
449
450   //  gAlice->MakeTree("R") ; 
451   
452   MakeBranch("R") ;
453   
454   // 3.
455
456   fReconstructioner->Make(fDigits, fEmcRecPoints, fPpsdRecPoints, fTrackSegments, fRecParticles);
457
458   printf("Reconstruction: %d %d %d %d\n",
459          fEmcRecPoints->GetEntries(),fPpsdRecPoints->GetEntries(),
460          fTrackSegments->GetEntries(),fRecParticles->GetEntries());
461
462   // 4. Expand or Shrink the arrays to the proper size
463   
464   Int_t size ;
465   
466   size = fEmcRecPoints->GetEntries() ;
467   fEmcRecPoints->Expand(size) ;
468
469   size = fPpsdRecPoints->GetEntries() ;
470   fPpsdRecPoints->Expand(size) ;
471
472   size = fTrackSegments->GetEntries() ;
473   fTrackSegments->Expand(size) ;
474
475   size = fRecParticles->GetEntries() ;
476   fRecParticles->Expand(size) ;
477
478   gAlice->TreeR()->Fill() ;
479   // 5.
480
481   gAlice->TreeR()->Write(0,TObject::kOverwrite) ;
482  
483   // Deleting reconstructed objects
484   ResetReconstruction();
485   
486 }
487
488 //____________________________________________________________________________
489 void AliPHOSv1::ResetReconstruction() 
490
491   // Deleting reconstructed objects
492
493   if ( fEmcRecPoints  )  fEmcRecPoints ->Delete();
494   if ( fPpsdRecPoints )  fPpsdRecPoints->Delete();
495   if ( fTrackSegments )  fTrackSegments->Delete();
496   if ( fRecParticles  )  fRecParticles ->Delete();
497   
498 }
499
500 //____________________________________________________________________________
501
502 void AliPHOSv1::StepManager(void)
503 {
504   // Accumulates hits as long as the track stays in a single crystal or PPSD gas Cell
505
506   Int_t          relid[4] ;           // (box, layer, row, column) indices
507   Int_t          absid    ;           // absolute cell ID number
508   Float_t        xyze[4]={0,0,0,0}  ; // position wrt MRS and energy deposited
509   TLorentzVector pos      ;           // Lorentz vector of the track current position
510   TLorentzVector pmom     ;        //momentum of the particle initiated hit
511   Float_t        xyd[3]={0,0,0}   ;   //local posiiton of the entering
512   Bool_t         entered = kFALSE ;  
513   Int_t          copy     ;
514
515   Int_t tracknumber =  gAlice->CurrentTrack() ; 
516   Int_t primary     =  gAlice->GetPrimary( gAlice->CurrentTrack() ); 
517   TString name      =  GetGeometry()->GetName() ; 
518   Int_t trackpid    =  0  ; 
519
520   if( gMC->IsTrackEntering() ){ // create hit with position and momentum of new particle, 
521                                 // but may be without energy deposition
522
523     // Current position of the hit in the local ref. system
524       gMC -> TrackPosition(pos);
525       Float_t xyzm[3], xyzd[3] ;
526       Int_t i;
527       for (i=0; i<3; i++) xyzm[i] = pos[i];
528       gMC -> Gmtod (xyzm, xyzd, 1);    // transform coordinate from master to daughter system
529       xyd[0]  = xyzd[0];
530       xyd[1]  =-xyzd[1];
531       xyd[2]  =-xyzd[2];
532
533       
534       // Current momentum of the hit's track in the local ref. system
535       gMC -> TrackMomentum(pmom);
536       Float_t pm[3], pd[3];
537       for (i=0; i<3; i++) pm[i]   = pmom[i];
538       gMC -> Gmtod (pm, pd, 2);        // transform 3-momentum from master to daughter system
539       pmom[0] = pd[0];
540       pmom[1] =-pd[1];
541       pmom[2] =-pd[2];
542
543       trackpid = gMC->TrackPid();
544       entered = kTRUE ;      // Mark to create hit even withou energy deposition
545
546   }
547
548
549   if ( name == "GPS2" || name == "MIXT" ) {            // ======> CPV is a GPS' PPSD
550
551     if( gMC->CurrentVolID(copy) == gMC->VolId("PPCE") ) // We are inside a gas cell 
552     {
553       gMC->TrackPosition(pos) ;
554       xyze[0] = pos[0] ;
555       xyze[1] = pos[1] ;
556       xyze[2] = pos[2] ;
557       xyze[3] = gMC->Edep() ; 
558
559       if ( (xyze[3] != 0) || entered ) { // there is deposited energy or new particle entering  PPSD
560         gMC->CurrentVolOffID(5, relid[0]) ;  // get the PHOS Module number
561         if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(5),"PHO1") == 0 ){
562           relid[0] += GetGeometry()->GetNModules() - GetGeometry()->GetNPPSDModules();
563         }
564         gMC->CurrentVolOffID(3, relid[1]) ;  // get the Micromegas Module number 
565       // 1-> GetGeometry()->GetNumberOfModulesPhi() * GetGeometry()->GetNumberOfModulesZ() upper
566       //   > GetGeometry()->GetNumberOfModulesPhi() * GetGeometry()->GetNumberOfModulesZ() lower
567         gMC->CurrentVolOffID(1, relid[2]) ;  // get the row number of the cell
568         gMC->CurrentVolID(relid[3]) ;        // get the column number 
569
570         // get the absolute Id number
571
572         GetGeometry()->RelToAbsNumbering(relid, absid) ; 
573
574         // add current hit to the hit list      
575           AddHit(fIshunt, primary, tracknumber, absid, xyze, trackpid, pmom, xyd);
576
577
578       } // there is deposited energy 
579     } // We are inside the gas of the CPV  
580   } // GPS2 configuration
581
582   if ( name == "IHEP" || name == "MIXT" ) {       // ======> CPV is a IHEP's one
583
584     // Yuri Kharlov, 28 September 2000
585
586     if( gMC->CurrentVolID(copy) == gMC->VolId("PCPQ") &&
587         entered &&
588         gMC->TrackCharge() != 0) {      
589       
590       // Digitize the current CPV hit:
591
592       // 1. find pad response and
593       
594       Int_t moduleNumber;
595       gMC->CurrentVolOffID(3,moduleNumber);
596       moduleNumber--;
597
598
599       TClonesArray *cpvDigits = new TClonesArray("AliPHOSCPVDigit",0);   // array of digits for current hit
600       CPVDigitize(pmom,xyd,moduleNumber,cpvDigits);
601       
602       Float_t xmean = 0;
603       Float_t zmean = 0;
604       Float_t qsum  = 0;
605       Int_t   idigit,ndigits;
606
607       // 2. go through the current digit list and sum digits in pads
608
609       ndigits = cpvDigits->GetEntriesFast();
610       for (idigit=0; idigit<ndigits-1; idigit++) {
611         AliPHOSCPVDigit  *cpvDigit1 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit);
612         Float_t x1 = cpvDigit1->GetXpad() ;
613         Float_t z1 = cpvDigit1->GetYpad() ;
614         for (Int_t jdigit=idigit+1; jdigit<ndigits; jdigit++) {
615           AliPHOSCPVDigit  *cpvDigit2 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(jdigit);
616           Float_t x2 = cpvDigit2->GetXpad() ;
617           Float_t z2 = cpvDigit2->GetYpad() ;
618           if (x1==x2 && z1==z2) {
619             Float_t qsum = cpvDigit1->GetQpad() + cpvDigit2->GetQpad() ;
620             cpvDigit2->SetQpad(qsum) ;
621             cpvDigits->RemoveAt(idigit) ;
622           }
623         }
624       }
625       cpvDigits->Compress() ;
626
627       // 3. add digits to temporary hit list fTmpHits
628
629       ndigits = cpvDigits->GetEntriesFast();
630       for (idigit=0; idigit<ndigits; idigit++) {
631         AliPHOSCPVDigit  *cpvDigit = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit);
632         relid[0] = moduleNumber + 1 ;                             // CPV (or PHOS) module number
633         relid[1] =-1 ;                                            // means CPV
634         relid[2] = cpvDigit->GetXpad() ;                          // column number of a pad
635         relid[3] = cpvDigit->GetYpad() ;                          // row    number of a pad
636         
637         // get the absolute Id number
638         GetGeometry()->RelToAbsNumbering(relid, absid) ; 
639
640         // add current digit to the temporary hit list
641         xyze[0] = 0. ;
642         xyze[1] = 0. ;
643         xyze[2] = 0. ;
644         xyze[3] = cpvDigit->GetQpad() ;                           // amplitude in a pad
645         primary = -1;                                             // No need in primary for CPV
646         AddHit(fIshunt, primary, tracknumber, absid, xyze, trackpid, pmom, xyd);
647
648         if (cpvDigit->GetQpad() > 0.02) {
649           xmean += cpvDigit->GetQpad() * (cpvDigit->GetXpad() + 0.5);
650           zmean += cpvDigit->GetQpad() * (cpvDigit->GetYpad() + 0.5);
651           qsum  += cpvDigit->GetQpad();
652         }
653       }
654       delete cpvDigits;
655     }
656   } // end of IHEP configuration
657   
658
659   if(gMC->CurrentVolID(copy) == gMC->VolId("PXTL") ) { //  We are inside a PBWO crystal
660     gMC->TrackPosition(pos) ;
661     xyze[0] = pos[0] ;
662     xyze[1] = pos[1] ;
663     xyze[2] = pos[2] ;
664     xyze[3] = gMC->Edep() ;
665
666   
667     if ( (xyze[3] != 0) || entered ) {  // Track is inside the crystal and deposits some energy or just entered 
668
669       gMC->CurrentVolOffID(10, relid[0]) ; // get the PHOS module number ;
670
671       if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(10),"PHO1") == 0 )
672         relid[0] += GetGeometry()->GetNModules() - GetGeometry()->GetNPPSDModules();      
673
674       relid[1] = 0   ;                    // means PBW04
675       gMC->CurrentVolOffID(4, relid[2]) ; // get the row number inside the module
676       gMC->CurrentVolOffID(3, relid[3]) ; // get the cell number inside the module
677       
678       // get the absolute Id number
679       GetGeometry()->RelToAbsNumbering(relid, absid) ; 
680
681       // add current hit to the hit list
682         AddHit(fIshunt, primary,tracknumber, absid, xyze, trackpid,pmom, xyd);
683
684
685     } // there is deposited energy
686   } // we are inside a PHOS Xtal
687
688
689 }
690
691 //____________________________________________________________________________
692 void AliPHOSv1::CPVDigitize (TLorentzVector p, Float_t *zxhit, Int_t moduleNumber, TClonesArray *cpvDigits)
693 {
694   // ------------------------------------------------------------------------
695   // Digitize one CPV hit:
696   // On input take exact 4-momentum p and position zxhit of the hit,
697   // find the pad response around this hit and
698   // put the amplitudes in the pads into array digits
699   //
700   // Author: Yuri Kharlov (after Serguei Sadovsky)
701   // 2 October 2000
702   // ------------------------------------------------------------------------
703
704   const Float_t kCelWr  = GetGeometry()->GetPadSizePhi()/2;  // Distance between wires (2 wires above 1 pad)
705   const Float_t kDetR   = 0.1;     // Relative energy fluctuation in track for 100 e-
706   const Float_t kdEdx   = 4.0;     // Average energy loss in CPV;
707   const Int_t   kNgamz  = 5;       // Ionization size in Z
708   const Int_t   kNgamx  = 9;       // Ionization size in Phi
709   const Float_t kNoise = 0.03;    // charge noise in one pad
710
711   Float_t rnor1,rnor2;
712
713   // Just a reminder on axes notation in the CPV module:
714   // axis Z goes along the beam
715   // axis X goes across the beam in the module plane
716   // axis Y is a normal to the module plane showing from the IP
717
718   Float_t hitX  = zxhit[0];
719   Float_t hitZ  =-zxhit[1];
720   Float_t pX    = p.Px();
721   Float_t pZ    =-p.Pz();
722   Float_t pNorm = p.Py();
723   Float_t eloss = kdEdx;
724
725   Float_t dZY   = pZ/pNorm * GetGeometry()->GetCPVGasThickness();
726   Float_t dXY   = pX/pNorm * GetGeometry()->GetCPVGasThickness();
727   gRandom->Rannor(rnor1,rnor2);
728   eloss *= (1 + kDetR*rnor1) *
729            TMath::Sqrt((1 + ( pow(dZY,2) + pow(dXY,2) ) / pow(GetGeometry()->GetCPVGasThickness(),2)));
730   Float_t zhit1 = hitZ + GetGeometry()->GetCPVActiveSize(1)/2 - dZY/2;
731   Float_t xhit1 = hitX + GetGeometry()->GetCPVActiveSize(0)/2 - dXY/2;
732   Float_t zhit2 = zhit1 + dZY;
733   Float_t xhit2 = xhit1 + dXY;
734
735   Int_t   iwht1 = (Int_t) (xhit1 / kCelWr);           // wire (x) coordinate "in"
736   Int_t   iwht2 = (Int_t) (xhit2 / kCelWr);           // wire (x) coordinate "out"
737
738   Int_t   nIter;
739   Float_t zxe[3][5];
740   if (iwht1==iwht2) {                      // incline 1-wire hit
741     nIter = 2;
742     zxe[0][0] = (zhit1 + zhit2 - dZY*0.57735) / 2;
743     zxe[1][0] = (iwht1 + 0.5) * kCelWr;
744     zxe[2][0] =  eloss/2;
745     zxe[0][1] = (zhit1 + zhit2 + dZY*0.57735) / 2;
746     zxe[1][1] = (iwht1 + 0.5) * kCelWr;
747     zxe[2][1] =  eloss/2;
748   }
749   else if (TMath::Abs(iwht1-iwht2) != 1) { // incline 3-wire hit
750     nIter = 3;
751     Int_t iwht3 = (iwht1 + iwht2) / 2;
752     Float_t xwht1 = (iwht1 + 0.5) * kCelWr; // wire 1
753     Float_t xwht2 = (iwht2 + 0.5) * kCelWr; // wire 2
754     Float_t xwht3 = (iwht3 + 0.5) * kCelWr; // wire 3
755     Float_t xwr13 = (xwht1 + xwht3) / 2;   // center 13
756     Float_t xwr23 = (xwht2 + xwht3) / 2;   // center 23
757     Float_t dxw1  = xhit1 - xwr13;
758     Float_t dxw2  = xhit2 - xwr23;
759     Float_t egm1  = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
760     Float_t egm2  = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
761     Float_t egm3  =           kCelWr / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
762     zxe[0][0] = (dXY*(xwr13-xwht1)/dXY + zhit1 + zhit1) / 2;
763     zxe[1][0] =  xwht1;
764     zxe[2][0] =  eloss * egm1;
765     zxe[0][1] = (dXY*(xwr23-xwht1)/dXY + zhit1 + zhit2) / 2;
766     zxe[1][1] =  xwht2;
767     zxe[2][1] =  eloss * egm2;
768     zxe[0][2] =  dXY*(xwht3-xwht1)/dXY + zhit1;
769     zxe[1][2] =  xwht3;
770     zxe[2][2] =  eloss * egm3;
771   }
772   else {                                   // incline 2-wire hit
773     nIter = 2;
774     Float_t xwht1 = (iwht1 + 0.5) * kCelWr;
775     Float_t xwht2 = (iwht2 + 0.5) * kCelWr;
776     Float_t xwr12 = (xwht1 + xwht2) / 2;
777     Float_t dxw1  = xhit1 - xwr12;
778     Float_t dxw2  = xhit2 - xwr12;
779     Float_t egm1  = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) );
780     Float_t egm2  = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) );
781     zxe[0][0] = (zhit1 + zhit2 - dZY*egm1) / 2;
782     zxe[1][0] =  xwht1;
783     zxe[2][0] =  eloss * egm1;
784     zxe[0][1] = (zhit1 + zhit2 + dZY*egm2) / 2;
785     zxe[1][1] =  xwht2;
786     zxe[2][1] =  eloss * egm2;
787   }
788
789   // Finite size of ionization region
790
791   Int_t nCellZ  = GetGeometry()->GetNumberOfCPVPadsZ();
792   Int_t nCellX  = GetGeometry()->GetNumberOfCPVPadsPhi();
793   Int_t nz3     = (kNgamz+1)/2;
794   Int_t nx3     = (kNgamx+1)/2;
795   cpvDigits->Expand(nIter*kNgamx*kNgamz);
796   TClonesArray &ldigits = *(TClonesArray *)cpvDigits;
797
798   for (Int_t iter=0; iter<nIter; iter++) {
799
800     Float_t zhit = zxe[0][iter];
801     Float_t xhit = zxe[1][iter];
802     Float_t qhit = zxe[2][iter];
803     Float_t zcell = zhit / GetGeometry()->GetPadSizeZ();
804     Float_t xcell = xhit / GetGeometry()->GetPadSizePhi();
805     if ( zcell<=0      || xcell<=0 ||
806          zcell>=nCellZ || xcell>=nCellX) return;
807     Int_t izcell = (Int_t) zcell;
808     Int_t ixcell = (Int_t) xcell;
809     Float_t zc = zcell - izcell - 0.5;
810     Float_t xc = xcell - ixcell - 0.5;
811     for (Int_t iz=1; iz<=kNgamz; iz++) {
812       Int_t kzg = izcell + iz - nz3;
813       if (kzg<=0 || kzg>nCellZ) continue;
814       Float_t zg = (Float_t)(iz-nz3) - zc;
815       for (Int_t ix=1; ix<=kNgamx; ix++) {
816         Int_t kxg = ixcell + ix - nx3;
817         if (kxg<=0 || kxg>nCellX) continue;
818         Float_t xg = (Float_t)(ix-nx3) - xc;
819         
820         // Now calculate pad response
821         Float_t qpad = CPVPadResponseFunction(qhit,zg,xg);
822         qpad += kNoise*rnor2;
823         if (qpad<0) continue;
824         
825         // Fill the array with pad response ID and amplitude
826         new(ldigits[cpvDigits->GetEntriesFast()]) AliPHOSCPVDigit(kxg,kzg,qpad);
827       }
828     }
829   }
830 }
831
832 //____________________________________________________________________________
833 Float_t AliPHOSv1::CPVPadResponseFunction(Float_t qhit, Float_t zhit, Float_t xhit) {
834   // ------------------------------------------------------------------------
835   // Calculate the amplitude in one CPV pad using the
836   // cumulative pad response function
837   // Author: Yuri Kharlov (after Serguei Sadovski)
838   // 3 October 2000
839   // ------------------------------------------------------------------------
840
841   Double_t dz = GetGeometry()->GetPadSizeZ()   / 2;
842   Double_t dx = GetGeometry()->GetPadSizePhi() / 2;
843   Double_t z  = zhit * GetGeometry()->GetPadSizeZ();
844   Double_t x  = xhit * GetGeometry()->GetPadSizePhi();
845   Double_t amplitude = qhit *
846     (CPVCumulPadResponse(z+dz,x+dx) - CPVCumulPadResponse(z+dz,x-dx) -
847      CPVCumulPadResponse(z-dz,x+dx) + CPVCumulPadResponse(z-dz,x-dx));
848   return (Float_t)amplitude;
849 }
850
851 //____________________________________________________________________________
852 Double_t AliPHOSv1::CPVCumulPadResponse(Double_t x, Double_t y) {
853   // ------------------------------------------------------------------------
854   // Cumulative pad response function
855   // It includes several terms from the CF decomposition in electrostatics
856   // Note: this cumulative function is wrong since omits some terms
857   //       but the cell amplitude obtained with it is correct because
858   //       these omitting terms cancel
859   // Author: Yuri Kharlov (after Serguei Sadovski)
860   // 3 October 2000
861   // ------------------------------------------------------------------------
862
863   const Double_t kA=1.0;
864   const Double_t kB=0.7;
865
866   Double_t r2       = x*x + y*y;
867   Double_t xy       = x*y;
868   Double_t cumulPRF = 0;
869   for (Int_t i=0; i<=4; i++) {
870     Double_t b1 = (2*i + 1) * kB;
871     cumulPRF += TMath::Power(-1,i) * TMath::ATan( xy / (b1*TMath::Sqrt(b1*b1 + r2)) );
872   }
873   cumulPRF *= kA/(2*TMath::Pi());
874   return cumulPRF;
875 }
876