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