]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSv1.cxx
8de0cdbe9256b2a8503accf8b825c5da2924e66d
[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 {
63   // ctor
64
65   // Create an empty array of AliPHOSCPVModule to satisfy
66   // AliPHOSv1::Streamer when reading root file
67  
68   fReconstructioner  = 0;
69   fTrackSegmentMaker = 0;
70  
71   if ( 0==(fEMCModules=new TClonesArray("AliPHOSCPVModule",0)) ) {
72     Error("AliPHOSv1","Can not create array of EMC modules");
73     exit(1);
74   }
75
76   if ( 0==(fCPVModules=new TClonesArray("AliPHOSCPVModule",0)) ) {
77     Error("AliPHOSv1","Can not create array of CPV modules");
78     exit(1);
79   }
80
81 }
82
83 //____________________________________________________________________________
84 AliPHOSv1::AliPHOSv1(const char *name, const char *title):
85 AliPHOSv0(name,title) 
86 {
87   // ctor : title is used to identify the layout
88   //        GPS2 = 5 modules (EMC + PPSD)
89   //        IHEP = 5 modules (EMC + CPV )
90   //        MIXT = 4 modules (EMC + CPV ) and 1 module (EMC + PPSD)
91   //
92   // We store hits :
93   //   - fHits (the "normal" one), which retains the hits associated with
94   //     the current primary particle being tracked
95   //     (this array is reset after each primary has been tracked).
96   //
97
98   fPinElectronicNoise = 0.010 ;
99   fDigitThreshold      = 0.1 ;   // 1 GeV 
100
101   // We do not want to save in TreeH the raw hits
102   // But save the cumulated hits instead (need to create the branch myself)
103   // It is put in the Digit Tree because the TreeH is filled after each primary
104   // and the TreeD at the end of the event (branch is set in FinishEvent() ).
105   
106   fHits= new TClonesArray("AliPHOSHit",1000) ;
107   gAlice->AddHitList(fHits);
108
109   fNhits = 0 ;
110
111   fReconstructioner  = 0;
112   fTrackSegmentMaker = 0;
113
114   fIshunt     =  1 ; // All hits are associated with primary particles
115   
116   // Create array of EMC modules of the size of PHOS modules number
117   
118   if ( 0==(fEMCModules=new TClonesArray("AliPHOSCPVModule",fGeom->GetNModules())) ) {
119     Error("AliPHOSv1","Can not create array of EMC modules");
120     exit(1);
121   }
122   TClonesArray &lemcmodule = *fEMCModules;
123   for (Int_t i=0; i<fGeom->GetNModules(); i++) new(lemcmodule[i]) AliPHOSCPVModule();
124
125   // Create array of CPV modules for the IHEP's version of CPV
126
127   if ( strcmp(fGeom->GetName(),"IHEP") == 0 || strcmp(fGeom->GetName(),"MIXT") == 0 ) {
128     // Create array of CPV modules of the size of PHOS modules number
129
130     if ( 0==(fCPVModules=new TClonesArray("AliPHOSCPVModule",fGeom->GetNCPVModules())) ) {
131       Error("AliPHOSv1","Can not create array of CPV modules");
132       exit(1);
133     }
134     TClonesArray &lcpvmodule = *fCPVModules;
135     for (Int_t i=0; i<fGeom->GetNCPVModules(); i++) new(lcpvmodule[i]) AliPHOSCPVModule();
136   }
137   else {
138     // Create an empty array of AliPHOSCPVModule to satisfy
139     // AliPHOSv1::Streamer when writing root file
140
141     fCPVModules=new TClonesArray("AliPHOSCPVModule",0);
142
143   }
144 }
145
146 //____________________________________________________________________________
147 AliPHOSv1::AliPHOSv1(AliPHOSReconstructioner * Reconstructioner, const char *name, const char *title):
148   AliPHOSv0(name,title)
149 {
150   // ctor : title is used to identify the layout
151   //        GPS2 = 5 modules (EMC + PPSD)   
152   // We use 2 arrays of hits :
153   //
154   //   - fHits (the "normal" one), which retains the hits associated with
155   //     the current primary particle being tracked
156   //     (this array is reset after each primary has been tracked).
157   //
158   //   - fTmpHits, which retains all the hits of the current event. It 
159   //     is used for the digitization part.
160
161   fPinElectronicNoise = 0.010 ;
162
163   // We do not want to save in TreeH the raw hits
164   //fHits   = new TClonesArray("AliPHOSHit",100) ;
165
166   fDigits = 0 ;
167   fHits= new TClonesArray("AliPHOSHit",1000) ;
168   gAlice->AddHitList(fHits);
169   fNhits = 0 ;
170
171   fIshunt     =  1 ; // All hits are associated with primary particles
172  
173   // gets an instance of the geometry parameters class  
174   fGeom =  AliPHOSGeometry::GetInstance(title, "") ; 
175
176   if (fGeom->IsInitialized() ) 
177     cout << "AliPHOS" << Version() << " : PHOS geometry intialized for " << fGeom->GetName() << endl ;
178   else
179     cout << "AliPHOS" << Version() << " : PHOS geometry initialization failed !" << endl ;   
180
181   // Defining the PHOS Reconstructioner
182  
183  fReconstructioner = Reconstructioner ;
184
185 }
186
187 //____________________________________________________________________________
188 AliPHOSv1::~AliPHOSv1()
189 {
190   // dtor
191
192   if ( fHits) {
193     fHits->Delete() ; 
194     delete fHits ;
195     fHits = 0 ; 
196   }
197
198   if ( fDigits) {
199     fDigits->Delete() ; 
200     delete fDigits ;
201     fDigits = 0 ; 
202   }
203
204   if ( fEmcRecPoints ) {
205     fEmcRecPoints->Delete() ; 
206     delete fEmcRecPoints ; 
207     fEmcRecPoints = 0 ; 
208   }
209   
210   if ( fPpsdRecPoints ) { 
211     fPpsdRecPoints->Delete() ;
212     delete fPpsdRecPoints ;
213     fPpsdRecPoints = 0 ; 
214   }
215   
216   if ( fTrackSegments ) {
217     fTrackSegments->Delete() ; 
218     delete fTrackSegments ;
219     fTrackSegments = 0 ; 
220   }
221   
222 }
223
224 //____________________________________________________________________________
225 void AliPHOSv1::AddHit(Int_t shunt, Int_t primary, Int_t tracknumber, Int_t Id, Float_t * hits, Int_t trackpid)
226 {
227   // Add a hit to the hit list.
228   // A PHOS hit is the sum of all hits in a single crystal
229   //   or in a single PPSD gas cell
230
231   Int_t hitCounter ;
232   AliPHOSHit *newHit ;
233   AliPHOSHit *curHit ;
234   Bool_t deja = kFALSE ;
235
236   newHit = new AliPHOSHit(shunt, primary, tracknumber, Id, hits, trackpid) ;
237
238   for ( hitCounter = 0 ; hitCounter < fNhits && !deja ; hitCounter++ ) {
239     curHit = (AliPHOSHit*) (*fHits)[hitCounter] ;
240     if( *curHit == *newHit ) {
241       *curHit = *curHit + *newHit ;
242       deja = kTRUE ;
243     }
244   }
245          
246   if ( !deja ) {
247     new((*fHits)[fNhits]) AliPHOSHit(*newHit) ;
248     fNhits++ ;
249   }
250
251   delete newHit;
252 }
253
254 //___________________________________________________________________________
255 Int_t AliPHOSv1::Digitize(Float_t Energy)
256 {
257   // Applies the energy calibration
258   
259   Float_t fB = 100000000. ;
260   Float_t fA = 0. ;
261   Int_t chan = Int_t(fA + Energy*fB ) ;
262   return chan ;
263 }
264
265 //____________________________________________________________________________
266 void AliPHOSv1::Hit2Digit(Int_t ntracks){
267   //Collects all hits in the same active volume into digits
268   
269   if(fDigits!= 0)
270     fDigits->Clear() ;
271   else
272     fDigits = new TClonesArray("AliPHOSDigit",1000) ;
273   
274   // Branch address for digit tree
275   char branchname[20];
276   sprintf(branchname,"%s",GetName());
277   gAlice->TreeD()->Branch(branchname,&fDigits,fBufferSize);  
278   
279   gAlice->TreeD()->GetEvent(0);
280
281   
282   Int_t i ;
283   Int_t relid[4];
284   Int_t j ; 
285   AliPHOSHit   * hit ;
286   AliPHOSDigit * newdigit ;
287   AliPHOSDigit * curdigit ;
288   Bool_t deja = kFALSE ; 
289   
290   Int_t itrack ;
291   for (itrack=0; itrack<ntracks; itrack++){
292     
293     //=========== Get the Hits Tree for the Primary track itrack
294     gAlice->ResetHits();    
295     gAlice->TreeH()->GetEvent(itrack);
296       
297     for ( i = 0 ; i < fHits->GetEntries() ; i++ ) {
298       hit = (AliPHOSHit*)fHits->At(i) ;
299     
300       // Assign primary number only if contribution is significant
301       if( hit->GetEnergy() > fDigitThreshold)
302         newdigit = new AliPHOSDigit( hit->GetPrimary(), hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
303       else
304         newdigit = new AliPHOSDigit( -1               , hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
305       deja =kFALSE ;
306
307
308       for ( j = 0 ; j < fNdigits ;  j++) { 
309         curdigit = (AliPHOSDigit*) fDigits->At(j) ;
310         if ( *curdigit == *newdigit) {
311           *curdigit = *curdigit + *newdigit ; 
312           deja = kTRUE ; 
313         }
314       }
315
316       if ( !deja ) {
317         new((*fDigits)[fNdigits]) AliPHOSDigit(* newdigit) ;
318         fNdigits++ ;  
319       }
320       
321       delete newdigit ;    
322     } 
323
324   } // loop over tracks
325     
326   // Noise induced by the PIN diode of the PbWO crystals
327   
328   Float_t energyandnoise ;
329   for ( i = 0 ; i < fNdigits ; i++ ) {
330     newdigit =  (AliPHOSDigit * ) fDigits->At(i) ;
331     
332     fGeom->AbsToRelNumbering(newdigit->GetId(), relid) ;
333     
334     if (relid[1]==0){   // Digits belong to EMC (PbW0_4 crystals)
335       energyandnoise = newdigit->GetAmp() + Digitize(gRandom->Gaus(0., fPinElectronicNoise)) ;
336       
337       if (energyandnoise < 0 ) 
338         energyandnoise = 0 ;
339         
340       if ( newdigit->GetAmp() < fDigitThreshold ) // if threshold not surpassed, remove digit from list
341         fDigits->RemoveAt(i) ; 
342     }
343   }
344     
345   fDigits->Compress() ;  
346   
347   fNdigits = fDigits->GetEntries() ;
348   fDigits->Expand(fNdigits) ;
349
350   for (i = 0 ; i < fNdigits ; i++) { 
351     newdigit = (AliPHOSDigit *) fDigits->At(i) ; 
352     newdigit->SetIndexInList(i) ;     
353   }
354
355   gAlice->TreeD()->Fill() ;
356
357   gAlice->TreeD()->Write(0,TObject::kOverwrite) ;
358  
359 }
360 //___________________________________________________________________________
361 void AliPHOSv1::MakeBranch(Option_t* opt)
362 {  
363   // Create new branche in the current Root Tree in the digit Tree
364   AliDetector::MakeBranch(opt) ;
365   
366   // Create new branches EMC<i> for hits in EMC modules
367   
368   for( Int_t i=0; i<fGeom->GetNModules(); i++ ) GetEMCModule(i).MakeBranch("EMC",i+1);
369   
370   // Create new branches CPV<i> for hits in CPV modules for IHEP geometry
371   
372   if ( strcmp(fGeom->GetName(),"IHEP") == 0 || strcmp(fGeom->GetName(),"MIXT") == 0 ) {
373     for( Int_t i=0; i<fGeom->GetNCPVModules(); i++ ) GetCPVModule(i).MakeBranch("CPV",i+1);
374   }
375   
376 }
377
378 //_____________________________________________________________________________
379 void AliPHOSv1::Reconstruction(AliPHOSReconstructioner * Reconstructioner)
380
381   // 1. Reinitializes the existing RecPoint, TrackSegment, and RecParticles Lists and 
382   // 2. Creates TreeR with a branch for each list
383   // 3. Steers the reconstruction processes
384   // 4. Saves the 3 lists in TreeR
385   // 5. Write the Tree to File
386   
387   fReconstructioner = Reconstructioner ;
388   
389   char branchname[10] ;
390   
391   // 1.
392
393   //  gAlice->MakeTree("R") ; 
394   Int_t splitlevel = 0 ; 
395   
396   fEmcRecPoints->Delete() ; 
397
398   if ( fEmcRecPoints && gAlice->TreeR() ) {
399     sprintf(branchname,"%sEmcRP",GetName()) ;
400     gAlice->TreeR()->Branch(branchname, "TObjArray", &fEmcRecPoints, fBufferSize, splitlevel) ; 
401   }
402
403   fPpsdRecPoints->Delete() ; 
404
405   if ( fPpsdRecPoints && gAlice->TreeR() ) {
406     sprintf(branchname,"%sPpsdRP",GetName()) ;
407     gAlice->TreeR()->Branch(branchname, "TObjArray", &fPpsdRecPoints, fBufferSize, splitlevel) ;
408   }
409
410   fTrackSegments->Delete() ; 
411
412   if ( fTrackSegments && gAlice->TreeR() ) { 
413     sprintf(branchname,"%sTS",GetName()) ;
414     gAlice->TreeR()->Branch(branchname, &fTrackSegments, fBufferSize) ;
415   }
416
417   fRecParticles->Delete() ; 
418
419   if      (strcmp(fGeom->GetName(),"GPS2") == 0 || strcmp(fGeom->GetName(),"MIXT") == 0) {
420     if ( fRecParticles && gAlice->TreeR() ) { 
421       sprintf(branchname,"%sRP",GetName()) ;
422       gAlice->TreeR()->Branch(branchname, &fRecParticles, fBufferSize) ;
423     }
424   }
425   
426   // 3.
427   if      (strcmp(fGeom->GetName(),"GPS2") == 0 || strcmp(fGeom->GetName(),"MIXT") == 0)
428     fReconstructioner->MakePPSD(fDigits, fEmcRecPoints, fPpsdRecPoints, fTrackSegments, fRecParticles);
429   if (strcmp(fGeom->GetName(),"IHEP") == 0 || strcmp(fGeom->GetName(),"MIXT") == 0)
430     fReconstructioner->MakeCPV (fDigits, fEmcRecPoints, fPpsdRecPoints);
431
432   // 4. Expand or Shrink the arrays to the proper size
433   
434   Int_t size ;
435   
436   size = fEmcRecPoints->GetEntries() ;
437   fEmcRecPoints->Expand(size) ;
438
439   size = fPpsdRecPoints->GetEntries() ;
440   fPpsdRecPoints->Expand(size) ;
441
442   size = fTrackSegments->GetEntries() ;
443   fTrackSegments->Expand(size) ;
444
445   size = fRecParticles->GetEntries() ;
446   fRecParticles->Expand(size) ;
447
448   gAlice->TreeR()->Fill() ;
449   // 5.
450
451   gAlice->TreeR()->Write(0,TObject::kOverwrite) ;
452  
453   // Deleting reconstructed objects
454   ResetReconstruction();
455   
456 }
457
458 //____________________________________________________________________________
459 void AliPHOSv1::ResetHits() 
460 {              
461   // Reset hit tree for EMC and CPV
462   // Yuri Kharlov, 28 September 2000
463
464   AliDetector::ResetHits();
465   for (Int_t i=0; i<fGeom->GetNModules(); i++) ((AliPHOSCPVModule*)(*fEMCModules)[i]) -> Clear();
466   if ( strcmp(fGeom->GetName(),"IHEP") == 0 || strcmp(fGeom->GetName(),"MIXT") == 0 ) {
467     for (Int_t i=0; i<fGeom->GetNCPVModules(); i++) ((AliPHOSCPVModule*)(*fCPVModules)[i]) -> Clear();
468   }
469   
470 }  
471 //____________________________________________________________________________
472 void AliPHOSv1::ResetReconstruction() 
473
474   // Deleting reconstructed objects
475
476   if ( fEmcRecPoints )   fEmcRecPoints->Delete();
477   if ( fPpsdRecPoints )  fPpsdRecPoints->Delete();
478   if ( fTrackSegments )  fTrackSegments->Delete();
479   if ( fRecParticles )   fRecParticles->Delete();
480   
481 }
482
483 //____________________________________________________________________________
484 void AliPHOSv1::SetTreeAddress()
485
486   //  TBranch *branch;
487   AliPHOS::SetTreeAddress();
488
489 //  //Branch address for TreeR: RecPpsdRecPoint
490 //   TTree *treeR = gAlice->TreeR();
491 //   if ( treeR && fPpsdRecPoints ) {
492 //     branch = treeR->GetBranch("PHOSPpsdRP");
493 //     if (branch) branch->SetAddress(&fPpsdRecPoints) ;
494 //  }
495
496   // Set branch address for the Hits Tree for hits in EMC modules
497   // Yuri Kharlov, 23 November 2000.
498
499   for( Int_t i=0; i<fGeom->GetNModules(); i++ ) GetEMCModule(i).SetTreeAddress("EMC",i+1);
500
501   // Set branch address for the Hits Tree for hits in CPV modules for IHEP geometry
502   // Yuri Kharlov, 28 September 2000.
503
504   if ( strcmp(fGeom->GetName(),"IHEP") == 0 || strcmp(fGeom->GetName(),"MIXT") == 0 ) {
505     for( Int_t i=0; i<fGeom->GetNCPVModules(); i++ ) GetCPVModule(i).SetTreeAddress("CPV",i+1);
506   }
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 //    if (gMC->IsTrackEntering())
517 //      cout << "Track enters the volume " << gMC->CurrentVolName() << endl;
518 //    if (gMC->IsTrackExiting())
519 //      cout << "Track leaves the volume " << gMC->CurrentVolName() << endl;
520
521   Int_t          relid[4] ;      // (box, layer, row, column) indices
522   Int_t          absid    ;      // absolute cell ID number
523   Float_t        xyze[4]={0,0,0,0}  ; // position wrt MRS and energy deposited
524   TLorentzVector pos      ;      // Lorentz vector of the track current position
525   Int_t          copy     ;
526   Int_t          i        ;
527
528   Int_t tracknumber =  gAlice->CurrentTrack() ; 
529   Int_t primary     =  gAlice->GetPrimary( gAlice->CurrentTrack() ); 
530   TString name      =  fGeom->GetName() ; 
531   Int_t trackpid    =  gMC->TrackPid() ; 
532   if ( name == "GPS2" || name == "MIXT" ) {            // ======> CPV is a GPS' PPSD
533
534     if( gMC->CurrentVolID(copy) == gMC->VolId("GCEL") ) // We are inside a gas cell 
535     {
536       gMC->TrackPosition(pos) ;
537       xyze[0] = pos[0] ;
538       xyze[1] = pos[1] ;
539       xyze[2] = pos[2] ;
540       xyze[3] = gMC->Edep() ; 
541
542       if ( xyze[3] != 0 ) { // there is deposited energy 
543         gMC->CurrentVolOffID(5, relid[0]) ;  // get the PHOS Module number
544         if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(5),"PHO1") == 0 )
545           relid[0] += fGeom->GetNModules() - fGeom->GetNPPSDModules();
546         gMC->CurrentVolOffID(3, relid[1]) ;  // get the Micromegas Module number 
547       // 1-> fGeom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() upper
548       //   > fGeom->GetNumberOfModulesPhi() * fGeom->GetNumberOfModulesZ() lower
549         gMC->CurrentVolOffID(1, relid[2]) ;  // get the row number of the cell
550         gMC->CurrentVolID(relid[3]) ;        // get the column number 
551
552         // get the absolute Id number
553
554         fGeom->RelToAbsNumbering(relid, absid) ; 
555
556         // add current hit to the hit list      
557         AddHit(fIshunt, primary, tracknumber, absid, xyze, trackpid);
558
559       } // there is deposited energy 
560     } // We are inside the gas of the CPV  
561   } // GPS2 configuration
562
563   if ( name == "IHEP" || name == "MIXT" ) {       // ======> CPV is a IHEP's one
564
565     // Yuri Kharlov, 28 September 2000
566
567     if( gMC->CurrentVolID(copy) == gMC->VolId("CPVQ") &&
568         gMC->IsTrackEntering() &&
569         gMC->TrackCharge() != 0) {
570
571       // Charged track has just entered to the CPV sensitive plane
572       
573       AliPHOSv1 &phos = *(AliPHOSv1*)gAlice->GetModule("PHOS");
574       
575       Int_t moduleNumber;
576       gMC->CurrentVolOffID(3,moduleNumber);
577       moduleNumber--;
578       
579       // Current position of the hit in the CPV module ref. system
580
581       gMC -> TrackPosition(pos);
582       Float_t xyzm[3], xyzd[3], xyd[3]={0,0,0};
583       Int_t i;
584       for (i=0; i<3; i++) xyzm[i] = pos[i];
585       gMC -> Gmtod (xyzm, xyzd, 1);    // transform coordinate from master to daughter system
586       xyd[0]  = xyzd[0];
587       xyd[1]  =-xyzd[2];
588       
589       // Current momentum of the hit's track in the CPV module ref. system
590       
591       TLorentzVector  pmom;
592       gMC -> TrackMomentum(pmom);
593       Float_t pm[3], pd[3];
594       for (i=0; i<3; i++) pm[i]   = pmom[i];
595       gMC -> Gmtod (pm, pd, 2);        // transform 3-momentum from master to daughter system
596       pmom[0] = pd[0];
597       pmom[1] =-pd[1];
598       pmom[2] =-pd[2];
599
600       // Current particle type of the hit's track
601
602       Int_t ipart = gMC->TrackPid();
603
604       // Add the current particle in the list of the CPV hits.
605
606       phos.GetCPVModule(moduleNumber).AddHit(fIshunt,primary,pmom,xyd,ipart);
607
608       if (fDebug == 1) {
609         printf("CPV hit added to module #%2d: p = (% .4f, % .4f, % .4f, % .4f) GeV,\n",
610                moduleNumber+1,pmom.Px(),pmom.Py(),pmom.Pz(),pmom.E());
611         printf( "                            xy = (%8.4f, %8.4f) cm, ipart = %d, primary = %d\n",
612                 xyd[0],xyd[1],ipart,primary);
613       }
614
615       // Digitize the current CPV hit:
616
617       // 1. find pad response and
618       
619       TClonesArray *cpvDigits = new TClonesArray("AliPHOSCPVDigit",0);   // array of digits for current hit
620       CPVDigitize(pmom,xyd,moduleNumber,cpvDigits);
621       
622       Float_t xmean = 0;
623       Float_t zmean = 0;
624       Float_t qsum  = 0;
625       Int_t   idigit,ndigits;
626
627       // 2. go through the current digit list and sum digits in pads
628
629       ndigits = cpvDigits->GetEntriesFast();
630       for (idigit=0; idigit<ndigits-1; idigit++) {
631         AliPHOSCPVDigit  *cpvDigit1 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit);
632         Float_t x1 = cpvDigit1->GetXpad() ;
633         Float_t z1 = cpvDigit1->GetYpad() ;
634         for (Int_t jdigit=idigit+1; jdigit<ndigits; jdigit++) {
635           AliPHOSCPVDigit  *cpvDigit2 = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(jdigit);
636           Float_t x2 = cpvDigit2->GetXpad() ;
637           Float_t z2 = cpvDigit2->GetYpad() ;
638           if (x1==x2 && z1==z2) {
639             Float_t qsum = cpvDigit1->GetQpad() + cpvDigit2->GetQpad() ;
640             cpvDigit2->SetQpad(qsum) ;
641             cpvDigits->RemoveAt(idigit) ;
642           }
643         }
644       }
645       cpvDigits->Compress() ;
646
647       // 3. add digits to temporary hit list fTmpHits
648
649       ndigits = cpvDigits->GetEntriesFast();
650       for (idigit=0; idigit<ndigits; idigit++) {
651         AliPHOSCPVDigit  *cpvDigit = (AliPHOSCPVDigit*) cpvDigits->UncheckedAt(idigit);
652         relid[0] = moduleNumber + 1 ;                             // CPV (or PHOS) module number
653         relid[1] =-1 ;                                            // means CPV
654         relid[2] = cpvDigit->GetXpad() ;                          // column number of a pad
655         relid[3] = cpvDigit->GetYpad() ;                          // row    number of a pad
656         
657         // get the absolute Id number
658         fGeom->RelToAbsNumbering(relid, absid) ; 
659
660         // add current digit to the temporary hit list
661         xyze[0] = 0. ;
662         xyze[1] = 0. ;
663         xyze[2] = 0. ;
664         xyze[3] = cpvDigit->GetQpad() ;                           // amplitude in a pad
665         primary = -1;                                             // No need in primary for CPV
666         AddHit(fIshunt, primary, tracknumber, absid, xyze, trackpid);
667
668         if (cpvDigit->GetQpad() > 0.02) {
669           xmean += cpvDigit->GetQpad() * (cpvDigit->GetXpad() + 0.5);
670           zmean += cpvDigit->GetQpad() * (cpvDigit->GetYpad() + 0.5);
671           qsum  += cpvDigit->GetQpad();
672         }
673       }
674       delete cpvDigits;
675     }
676   } // end of IHEP configuration
677   
678   if(gMC->CurrentVolID(copy) == gMC->VolId("PXTL") ) { //  We are inside a PBWO crystal
679     gMC->TrackPosition(pos) ;
680     xyze[0] = pos[0] ;
681     xyze[1] = pos[1] ;
682     xyze[2] = pos[2] ;
683     xyze[3] = gMC->Edep() ;
684
685     // Track enters to the crystal from the top edge
686
687     if (gMC->IsTrackEntering()) {
688       Float_t posloc[3];
689       gMC -> Gmtod (xyze, posloc, 1);
690       if (posloc[1] > fGeom->GetCrystalSize(1)/2-0.01) {
691         Int_t row,cel;
692         Float_t xyd[3]={0,0,0};
693         AliPHOSv1 &phos = *(AliPHOSv1*)gAlice->GetModule("PHOS");
694
695         Int_t moduleNumber;
696         gMC->CurrentVolOffID(10,moduleNumber);
697         if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(10),"PHO1") == 0 )
698           moduleNumber += fGeom->GetNModules() - fGeom->GetNPPSDModules();
699         moduleNumber--;
700
701         gMC->CurrentVolOffID(4, row) ; // get the row  number inside the module
702         gMC->CurrentVolOffID(3, cel) ; // get the cell number inside the module
703         xyd[0] = -(posloc[2] + (cel-0.5-fGeom->GetNZ()  /2) *
704           (fGeom->GetCrystalSize(2) + 2 * fGeom->GetGapBetweenCrystals()));
705         xyd[1] =   posloc[0] + (row-0.5-fGeom->GetNPhi()/2) *
706           (fGeom->GetCrystalSize(0) + 2 * fGeom->GetGapBetweenCrystals());
707
708         // Current momentum of the hit's track in the CPV module ref. system
709         
710         TLorentzVector  pmom;
711         gMC -> TrackMomentum(pmom);
712         Float_t pm[3], pd[3];
713         for (i=0; i<3; i++) pm[i]   = pmom[i];
714         gMC -> Gmtod (pm, pd, 2);        // transform 3-momentum from master to daughter system
715         pmom[0] = pd[0];
716         pmom[1] =-pd[1];
717         pmom[2] =-pd[2];
718         
719         // Current particle type of the hit's track
720         
721         Int_t ipart = gMC->TrackPid();
722
723         // Add the current particle in the list of the EMC hits.
724
725         phos.GetEMCModule(moduleNumber).AddHit(fIshunt,primary,pmom,xyd,ipart);
726
727         if (fDebug == 1) {
728           printf("EMC hit added to module #%2d: p = (% .4f, % .4f, % .4f, % .4f) GeV,\n",
729                  moduleNumber+1,pmom.Px(),pmom.Py(),pmom.Pz(),pmom.E());
730           printf( "                            xy = (%8.4f, %8.4f) cm, ipart = %d, primary = %d\n",
731                   xyd[0],xyd[1],ipart,primary);
732         }
733       }
734     }
735
736     // Track is inside the crystal and deposits some energy
737
738     if ( xyze[3] != 0 ) {
739       gMC->CurrentVolOffID(10, relid[0]) ; // get the PHOS module number ;
740       if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(5),"PHO1") == 0 )
741         relid[0] += fGeom->GetNModules() - fGeom->GetNPPSDModules();
742       relid[1] = 0   ;                    // means PBW04
743       gMC->CurrentVolOffID(4, relid[2]) ; // get the row number inside the module
744       gMC->CurrentVolOffID(3, relid[3]) ; // get the cell number inside the module
745       
746       // get the absolute Id number
747       
748       fGeom->RelToAbsNumbering(relid, absid) ; 
749
750       // add current hit to the hit list
751       
752       AddHit(fIshunt, primary,tracknumber, absid, xyze, trackpid);
753
754     } // there is deposited energy
755   } // we are inside a PHOS Xtal
756 }
757
758 //____________________________________________________________________________
759 void AliPHOSv1::CPVDigitize (TLorentzVector p, Float_t *zxhit, Int_t moduleNumber, TClonesArray *cpvDigits)
760 {
761   // ------------------------------------------------------------------------
762   // Digitize one CPV hit:
763   // On input take exact 4-momentum p and position zxhit of the hit,
764   // find the pad response around this hit and
765   // put the amplitudes in the pads into array digits
766   //
767   // Author: Yuri Kharlov (after Serguei Sadovsky)
768   // 2 October 2000
769   // ------------------------------------------------------------------------
770
771   const Float_t kCelWr  = fGeom->GetPadSizePhi()/2;  // Distance between wires (2 wires above 1 pad)
772   const Float_t kDetR   = 0.1;     // Relative energy fluctuation in track for 100 e-
773   const Float_t kdEdx   = 4.0;     // Average energy loss in CPV;
774   const Int_t   kNgamz  = 5;       // Ionization size in Z
775   const Int_t   kNgamx  = 9;       // Ionization size in Phi
776   const Float_t kNoise = 0.03;    // charge noise in one pad
777
778   Float_t rnor1,rnor2;
779
780   // Just a reminder on axes notation in the CPV module:
781   // axis Z goes along the beam
782   // axis X goes across the beam in the module plane
783   // axis Y is a normal to the module plane showing from the IP
784
785   Float_t hitX  = zxhit[0];
786   Float_t hitZ  =-zxhit[1];
787   Float_t pX    = p.Px();
788   Float_t pZ    =-p.Pz();
789   Float_t pNorm = p.Py();
790   Float_t eloss = kdEdx;
791
792 //    cout << "CPVDigitize: YVK : "<<hitX<<" "<<hitZ<<" | "<<pX<<" "<<pZ<<" "<<pNorm<<endl;
793
794   Float_t dZY   = pZ/pNorm * fGeom->GetCPVGasThickness();
795   Float_t dXY   = pX/pNorm * fGeom->GetCPVGasThickness();
796   gRandom->Rannor(rnor1,rnor2);
797   eloss *= (1 + kDetR*rnor1) *
798            TMath::Sqrt((1 + ( pow(dZY,2) + pow(dXY,2) ) / pow(fGeom->GetCPVGasThickness(),2)));
799   Float_t zhit1 = hitZ + fGeom->GetCPVActiveSize(1)/2 - dZY/2;
800   Float_t xhit1 = hitX + fGeom->GetCPVActiveSize(0)/2 - dXY/2;
801   Float_t zhit2 = zhit1 + dZY;
802   Float_t xhit2 = xhit1 + dXY;
803
804   Int_t   iwht1 = (Int_t) (xhit1 / kCelWr);           // wire (x) coordinate "in"
805   Int_t   iwht2 = (Int_t) (xhit2 / kCelWr);           // wire (x) coordinate "out"
806
807   Int_t   nIter;
808   Float_t zxe[3][5];
809   if (iwht1==iwht2) {                      // incline 1-wire hit
810     nIter = 2;
811     zxe[0][0] = (zhit1 + zhit2 - dZY*0.57735) / 2;
812     zxe[1][0] = (iwht1 + 0.5) * kCelWr;
813     zxe[2][0] =  eloss/2;
814     zxe[0][1] = (zhit1 + zhit2 + dZY*0.57735) / 2;
815     zxe[1][1] = (iwht1 + 0.5) * kCelWr;
816     zxe[2][1] =  eloss/2;
817   }
818   else if (TMath::Abs(iwht1-iwht2) != 1) { // incline 3-wire hit
819     nIter = 3;
820     Int_t iwht3 = (iwht1 + iwht2) / 2;
821     Float_t xwht1 = (iwht1 + 0.5) * kCelWr; // wire 1
822     Float_t xwht2 = (iwht2 + 0.5) * kCelWr; // wire 2
823     Float_t xwht3 = (iwht3 + 0.5) * kCelWr; // wire 3
824     Float_t xwr13 = (xwht1 + xwht3) / 2;   // center 13
825     Float_t xwr23 = (xwht2 + xwht3) / 2;   // center 23
826     Float_t dxw1  = xhit1 - xwr13;
827     Float_t dxw2  = xhit2 - xwr23;
828     Float_t egm1  = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
829     Float_t egm2  = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
830     Float_t egm3  =           kCelWr / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) + kCelWr );
831     zxe[0][0] = (dXY*(xwr13-xwht1)/dXY + zhit1 + zhit1) / 2;
832     zxe[1][0] =  xwht1;
833     zxe[2][0] =  eloss * egm1;
834     zxe[0][1] = (dXY*(xwr23-xwht1)/dXY + zhit1 + zhit2) / 2;
835     zxe[1][1] =  xwht2;
836     zxe[2][1] =  eloss * egm2;
837     zxe[0][2] =  dXY*(xwht3-xwht1)/dXY + zhit1;
838     zxe[1][2] =  xwht3;
839     zxe[2][2] =  eloss * egm3;
840   }
841   else {                                   // incline 2-wire hit
842     nIter = 2;
843     Float_t xwht1 = (iwht1 + 0.5) * kCelWr;
844     Float_t xwht2 = (iwht2 + 0.5) * kCelWr;
845     Float_t xwr12 = (xwht1 + xwht2) / 2;
846     Float_t dxw1  = xhit1 - xwr12;
847     Float_t dxw2  = xhit2 - xwr12;
848     Float_t egm1  = TMath::Abs(dxw1) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) );
849     Float_t egm2  = TMath::Abs(dxw2) / ( TMath::Abs(dxw1) + TMath::Abs(dxw2) );
850     zxe[0][0] = (zhit1 + zhit2 - dZY*egm1) / 2;
851     zxe[1][0] =  xwht1;
852     zxe[2][0] =  eloss * egm1;
853     zxe[0][1] = (zhit1 + zhit2 + dZY*egm2) / 2;
854     zxe[1][1] =  xwht2;
855     zxe[2][1] =  eloss * egm2;
856   }
857
858   // Finite size of ionization region
859
860   Int_t nCellZ  = fGeom->GetNumberOfCPVPadsZ();
861   Int_t nCellX  = fGeom->GetNumberOfCPVPadsPhi();
862   Int_t nz3     = (kNgamz+1)/2;
863   Int_t nx3     = (kNgamx+1)/2;
864   cpvDigits->Expand(nIter*kNgamx*kNgamz);
865   TClonesArray &ldigits = *(TClonesArray *)cpvDigits;
866
867   for (Int_t iter=0; iter<nIter; iter++) {
868
869     Float_t zhit = zxe[0][iter];
870     Float_t xhit = zxe[1][iter];
871     Float_t qhit = zxe[2][iter];
872     Float_t zcell = zhit / fGeom->GetPadSizeZ();
873     Float_t xcell = xhit / fGeom->GetPadSizePhi();
874     if ( zcell<=0      || xcell<=0 ||
875          zcell>=nCellZ || xcell>=nCellX) return;
876     Int_t izcell = (Int_t) zcell;
877     Int_t ixcell = (Int_t) xcell;
878     Float_t zc = zcell - izcell - 0.5;
879     Float_t xc = xcell - ixcell - 0.5;
880     for (Int_t iz=1; iz<=kNgamz; iz++) {
881       Int_t kzg = izcell + iz - nz3;
882       if (kzg<=0 || kzg>nCellZ) continue;
883       Float_t zg = (Float_t)(iz-nz3) - zc;
884       for (Int_t ix=1; ix<=kNgamx; ix++) {
885         Int_t kxg = ixcell + ix - nx3;
886         if (kxg<=0 || kxg>nCellX) continue;
887         Float_t xg = (Float_t)(ix-nx3) - xc;
888         
889         // Now calculate pad response
890         Float_t qpad = CPVPadResponseFunction(qhit,zg,xg);
891         qpad += kNoise*rnor2;
892         if (qpad<0) continue;
893         
894         // Fill the array with pad response ID and amplitude
895         new(ldigits[cpvDigits->GetEntriesFast()]) AliPHOSCPVDigit(kxg,kzg,qpad);
896       }
897     }
898   }
899 }
900
901 //____________________________________________________________________________
902 Float_t AliPHOSv1::CPVPadResponseFunction(Float_t qhit, Float_t zhit, Float_t xhit) {
903   // ------------------------------------------------------------------------
904   // Calculate the amplitude in one CPV pad using the
905   // cumulative pad response function
906   // Author: Yuri Kharlov (after Serguei Sadovski)
907   // 3 October 2000
908   // ------------------------------------------------------------------------
909
910   Double_t dz = fGeom->GetPadSizeZ()   / 2;
911   Double_t dx = fGeom->GetPadSizePhi() / 2;
912   Double_t z  = zhit * fGeom->GetPadSizeZ();
913   Double_t x  = xhit * fGeom->GetPadSizePhi();
914   Double_t amplitude = qhit *
915     (CPVCumulPadResponse(z+dz,x+dx) - CPVCumulPadResponse(z+dz,x-dx) -
916      CPVCumulPadResponse(z-dz,x+dx) + CPVCumulPadResponse(z-dz,x-dx));
917   return (Float_t)amplitude;
918 }
919
920 //____________________________________________________________________________
921 Double_t AliPHOSv1::CPVCumulPadResponse(Double_t x, Double_t y) {
922   // ------------------------------------------------------------------------
923   // Cumulative pad response function
924   // It includes several terms from the CF decomposition in electrostatics
925   // Note: this cumulative function is wrong since omits some terms
926   //       but the cell amplitude obtained with it is correct because
927   //       these omitting terms cancel
928   // Author: Yuri Kharlov (after Serguei Sadovski)
929   // 3 October 2000
930   // ------------------------------------------------------------------------
931
932   const Double_t kA=1.0;
933   const Double_t kB=0.7;
934
935   Double_t r2       = x*x + y*y;
936   Double_t xy       = x*y;
937   Double_t cumulPRF = 0;
938   for (Int_t i=0; i<=4; i++) {
939     Double_t b1 = (2*i + 1) * kB;
940     cumulPRF += TMath::Power(-1,i) * TMath::ATan( xy / (b1*TMath::Sqrt(b1*b1 + r2)) );
941   }
942   cumulPRF *= kA/(2*TMath::Pi());
943   return cumulPRF;
944 }
945