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