]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSv1.cxx
56e5a4e9b6b17db216689fc784b1b2e93f72c86e
[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 ( fRecParticles && gAlice->TreeR() ) { 
419     sprintf(branchname,"%sRP",GetName()) ;
420     gAlice->TreeR()->Branch(branchname, &fRecParticles, fBufferSize) ;
421   }
422   
423   // 3.
424
425   fReconstructioner->Make(fDigits, fEmcRecPoints, fPpsdRecPoints, fTrackSegments, fRecParticles);
426
427   printf("Reconstruction: %d %d %d %d\n",
428          fEmcRecPoints->GetEntries(),fPpsdRecPoints->GetEntries(),
429          fTrackSegments->GetEntries(),fRecParticles->GetEntries());
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         }
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
698         if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(10),"PHO1") == 0 )
699           moduleNumber += fGeom->GetNModules() - fGeom->GetNPPSDModules();
700         moduleNumber--;
701
702         gMC->CurrentVolOffID(4, row) ; // get the row  number inside the module
703         gMC->CurrentVolOffID(3, cel) ; // get the cell number inside the module
704         xyd[0] = -(posloc[2] + (cel-0.5-fGeom->GetNZ()  /2) *
705           (fGeom->GetCrystalSize(2) + 2 * fGeom->GetGapBetweenCrystals()));
706         xyd[1] =   posloc[0] + (row-0.5-fGeom->GetNPhi()/2) *
707           (fGeom->GetCrystalSize(0) + 2 * fGeom->GetGapBetweenCrystals());
708
709         // Current momentum of the hit's track in the CPV module ref. system
710         
711         TLorentzVector  pmom;
712         gMC -> TrackMomentum(pmom);
713         Float_t pm[3], pd[3];
714         for (i=0; i<3; i++) pm[i]   = pmom[i];
715         gMC -> Gmtod (pm, pd, 2);        // transform 3-momentum from master to daughter system
716         pmom[0] = pd[0];
717         pmom[1] =-pd[1];
718         pmom[2] =-pd[2];
719         
720         // Current particle type of the hit's track
721         
722         Int_t ipart = gMC->TrackPid();
723
724         // Add the current particle in the list of the EMC hits.
725
726         phos.GetEMCModule(moduleNumber).AddHit(fIshunt,primary,pmom,xyd,ipart);
727
728         if (fDebug == 1) {
729           printf("EMC hit added to module #%2d: p = (% .4f, % .4f, % .4f, % .4f) GeV,\n",
730                  moduleNumber+1,pmom.Px(),pmom.Py(),pmom.Pz(),pmom.E());
731           printf( "                            xy = (%8.4f, %8.4f) cm, ipart = %d, primary = %d\n",
732                   xyd[0],xyd[1],ipart,primary);
733         }
734       }
735     }
736
737     // Track is inside the crystal and deposits some energy
738
739     if ( xyze[3] != 0 ) {
740       gMC->CurrentVolOffID(10, relid[0]) ; // get the PHOS module number ;
741       if ( name == "MIXT" && strcmp(gMC->CurrentVolOffName(10),"PHO1") == 0 )
742         relid[0] += fGeom->GetNModules() - fGeom->GetNPPSDModules();      
743
744       relid[1] = 0   ;                    // means PBW04
745       gMC->CurrentVolOffID(4, relid[2]) ; // get the row number inside the module
746       gMC->CurrentVolOffID(3, relid[3]) ; // get the cell number inside the module
747       
748       // get the absolute Id number
749       fGeom->RelToAbsNumbering(relid, absid) ; 
750
751       // add current hit to the hit list
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