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