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