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