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