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