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