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