]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSv1.cxx
Pointer to RPlists etc changed to pointer-to-pointer to listRP etc -- to allow readin...
[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 // Layout EMC + PPSD has name GPS2  
21 // Produces cumulated hits (no hits) and digits                  
22 //*-- Author: Yves Schutz (SUBATECH)
23
24
25 // --- ROOT system ---
26
27 #include "TBRIK.h"
28 #include "TNode.h"
29 #include "TRandom.h"
30
31
32 // --- Standard library ---
33
34 #include <stdio.h>
35 #include <string.h>
36 #include <stdlib.h>
37 #include <strstream.h>
38
39 // --- AliRoot header files ---
40
41 #include "AliPHOSv1.h"
42 #include "AliPHOSHit.h"
43 #include "AliPHOSDigit.h"
44 #include "AliPHOSReconstructioner.h"
45 #include "AliRun.h"
46 #include "AliConst.h"
47
48 ClassImp(AliPHOSv1)
49
50 //____________________________________________________________________________
51 AliPHOSv1::AliPHOSv1()
52 {
53   // ctor
54   fNTmpHits = 0 ; 
55   fTmpHits  = 0 ; 
56
57
58 }
59
60 //____________________________________________________________________________
61 AliPHOSv1::AliPHOSv1(const char *name, const char *title):
62 AliPHOSv0(name,title) 
63 {
64   // ctor : title is used to identify the layout
65   //        GPS2 = 5 modules (EMC + PPSD)   
66   // We use 2 arrays of hits :
67   //
68   //   - fHits (the "normal" one), which retains the hits associated with
69   //     the current primary particle being tracked
70   //     (this array is reset after each primary has been tracked).
71   //
72   //   - fTmpHits, which retains all the hits of the current event. It 
73   //     is used for the digitization part.
74  
75   fPinElectronicNoise = 0.010 ;
76   fDigitThreshold      = 0.1 ;   // 1 GeV 
77
78   // We do not want to save in TreeH the raw hits
79   // But save the cumulated hits instead (need to create the branch myself)
80   // It is put in the Digit Tree because the TreeH is filled after each primary
81  // and the TreeD at the end of the event (branch is set in FinishEvent() ).
82   
83   fTmpHits= new TClonesArray("AliPHOSHit",1000) ;
84
85   fNTmpHits = fNhits = 0 ;
86
87   fDigits = new TClonesArray("AliPHOSDigit",1000) ;
88
89
90   fIshunt     =  1 ; // All hits are associated with primary particles
91  
92 }
93
94 //____________________________________________________________________________
95 AliPHOSv1::AliPHOSv1(AliPHOSReconstructioner * Reconstructioner, const char *name, const char *title):
96   AliPHOSv0(name,title)
97 {
98   // ctor : title is used to identify the layout
99   //        GPS2 = 5 modules (EMC + PPSD)   
100   // We use 2 arrays of hits :
101   //
102   //   - fHits (the "normal" one), which retains the hits associated with
103   //     the current primary particle being tracked
104   //     (this array is reset after each primary has been tracked).
105   //
106   //   - fTmpHits, which retains all the hits of the current event. It 
107   //     is used for the digitization part.
108
109   fPinElectronicNoise = 0.010 ;
110
111   // We do not want to save in TreeH the raw hits
112   //fHits   = new TClonesArray("AliPHOSHit",100) ;
113
114   fDigits = new TClonesArray("AliPHOSDigit",1000) ;
115   fTmpHits= new TClonesArray("AliPHOSHit",1000) ;
116
117   fNTmpHits = fNhits = 0 ;
118
119   fIshunt     =  1 ; // All hits are associated with primary particles
120  
121   // gets an instance of the geometry parameters class  
122   fGeom =  AliPHOSGeometry::GetInstance(title, "") ; 
123
124   if (fGeom->IsInitialized() ) 
125     cout << "AliPHOS" << Version() << " : PHOS geometry intialized for " << fGeom->GetName() << endl ;
126   else
127     cout << "AliPHOS" << Version() << " : PHOS geometry initialization failed !" << endl ;   
128
129   // Defining the PHOS Reconstructioner
130  
131  fReconstructioner = Reconstructioner ;
132
133 }
134
135 //____________________________________________________________________________
136 AliPHOSv1::~AliPHOSv1()
137 {
138   // dtor
139
140   if ( fTmpHits) {
141     fTmpHits->Delete() ; 
142     delete fTmpHits ;
143     fTmpHits = 0 ; 
144   }
145
146 //   if ( fEmcRecPoints ) {
147 //     fEmcRecPoints->Delete() ; 
148 //     delete fEmcRecPoints ; 
149 //     fEmcRecPoints = 0 ; 
150 //   }
151
152 //   if ( fPpsdRecPoints ) { 
153 //     fPpsdRecPoints->Delete() ;
154 //     delete fPpsdRecPoints ;
155 //     fPpsdRecPoints = 0 ; 
156 //   }
157   
158 //   if ( fTrackSegments ) {
159 //     fTrackSegments->Delete() ; 
160 //     delete fTrackSegments ;
161 //     fTrackSegments = 0 ; 
162 //   }
163
164 }
165
166 //____________________________________________________________________________
167 void AliPHOSv1::AddHit(Int_t shunt, Int_t primary, Int_t tracknumber, Int_t Id, Float_t * hits)
168 {
169   // Add a hit to the hit list.
170   // A PHOS hit is the sum of all hits in a single crystal
171   //   or in a single PPSD gas cell
172
173   Int_t hitCounter ;
174   TClonesArray &ltmphits = *fTmpHits ;
175   AliPHOSHit *newHit ;
176   AliPHOSHit *curHit ;
177   Bool_t deja = kFALSE ;
178
179   // In any case, fills the fTmpHit TClonesArray (with "accumulated hits")
180
181   newHit = new AliPHOSHit(shunt, primary, tracknumber, Id, hits) ;
182
183   // We do not want to save in TreeH the raw hits 
184   //  TClonesArray &lhits = *fHits;
185
186   for ( hitCounter = 0 ; hitCounter < fNTmpHits && !deja ; hitCounter++ ) {
187     curHit = (AliPHOSHit*) ltmphits[hitCounter] ;
188   if( *curHit == *newHit ) {
189     *curHit = *curHit + *newHit ;
190     deja = kTRUE ;
191     }
192   }
193          
194   if ( !deja ) {
195     new(ltmphits[fNTmpHits]) AliPHOSHit(*newHit) ;
196     fNTmpHits++ ;
197   }
198
199   // We do not want to save in TreeH the raw hits 
200   //   new(lhits[fNhits]) AliPHOSHit(*newHit) ;    
201   //   fNhits++ ;
202
203   // Please note that the fTmpHits array must survive up to the
204   // end of the events, so it does not appear e.g. in ResetHits() (
205   // which is called at the end of each primary).  
206
207   delete newHit;
208
209 }
210
211 //___________________________________________________________________________
212 Int_t AliPHOSv1::Digitize(Float_t Energy)
213 {
214   // Applies the energy calibration
215   
216   Float_t fB = 100000000. ;
217   Float_t fA = 0. ;
218   Int_t chan = Int_t(fA + Energy*fB ) ;
219   return chan ;
220 }
221
222 //___________________________________________________________________________
223 void AliPHOSv1::FinishEvent()
224 {
225   // Makes the digits from the sum of summed hit in a single crystal or PPSD gas cell
226   // Adds to the energy the electronic noise
227   // Keeps digits with energy above fDigitThreshold
228
229   // Save the cumulated hits instead of raw hits (need to create the branch myself)
230   // It is put in the Digit Tree because the TreeH is filled after each primary
231   // and the TreeD at the end of the event.
232   
233   
234   Int_t i ;
235   Int_t relid[4];
236   Int_t j ; 
237   TClonesArray &lDigits = *fDigits ;
238   AliPHOSHit  * hit ;
239   AliPHOSDigit * newdigit ;
240   AliPHOSDigit * curdigit ;
241   Bool_t deja = kFALSE ; 
242   
243   for ( i = 0 ; i < fNTmpHits ; i++ ) {
244     hit = (AliPHOSHit*)fTmpHits->At(i) ;
245
246     // Assign primary number only if contribution is significant
247     if( hit->GetEnergy() > fDigitThreshold)
248       newdigit = new AliPHOSDigit( hit->GetPrimary(), hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
249     else
250       newdigit = new AliPHOSDigit( -1 , hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
251     deja =kFALSE ;
252     for ( j = 0 ; j < fNdigits ;  j++) { 
253       curdigit = (AliPHOSDigit*) lDigits[j] ;
254       if ( *curdigit == *newdigit) {
255         *curdigit = *curdigit + *newdigit ; 
256         deja = kTRUE ; 
257       }
258     }
259     if ( !deja ) {
260       new(lDigits[fNdigits]) AliPHOSDigit(* newdigit) ;
261       fNdigits++ ;  
262     }
263  
264     delete newdigit ;    
265   } 
266   
267   // Noise induced by the PIN diode of the PbWO crystals
268
269   Float_t energyandnoise ;
270   for ( i = 0 ; i < fNdigits ; i++ ) {
271     newdigit =  (AliPHOSDigit * ) fDigits->At(i) ;
272     fGeom->AbsToRelNumbering(newdigit->GetId(), relid) ;
273
274     if (relid[1]==0){   // Digits belong to EMC (PbW0_4 crystals)
275       energyandnoise = newdigit->GetAmp() + Digitize(gRandom->Gaus(0., fPinElectronicNoise)) ;
276
277       if (energyandnoise < 0 ) 
278         energyandnoise = 0 ;
279
280       if ( newdigit->GetAmp() < fDigitThreshold ) // if threshold not surpassed, remove digit from list
281         fDigits->RemoveAt(i) ; 
282     }
283   }
284   
285   fDigits->Compress() ;  
286
287   fNdigits =  fDigits->GetEntries() ; 
288   for (i = 0 ; i < fNdigits ; i++) { 
289     newdigit = (AliPHOSDigit *) fDigits->At(i) ; 
290     newdigit->SetIndexInList(i) ; 
291   }
292  
293 }
294
295 //___________________________________________________________________________
296 void AliPHOSv1::MakeBranch(Option_t* opt)
297 {  
298   // Create new branche in the current Root Tree in the digit Tree
299   AliDetector::MakeBranch(opt) ;
300   
301   char branchname[10];
302   sprintf(branchname,"%s",GetName());
303   char *cdD = strstr(opt,"D");
304   if (fDigits && gAlice->TreeD() && cdD) {
305     gAlice->TreeD()->Branch(branchname, &fDigits, fBufferSize);
306   }
307
308   // Create new branche PHOSCH in the current Root Tree in the digit Tree for accumulated Hits
309   if ( ! (gAlice->IsLegoRun()) ) { // only when not in lego plot mode 
310     if ( fTmpHits && gAlice->TreeD()  && cdD) {
311       char branchname[10] ;
312       sprintf(branchname, "%sCH", GetName()) ;
313       gAlice->TreeD()->Branch(branchname, &fTmpHits, fBufferSize) ;
314     }   
315   }
316
317 }
318
319 //_____________________________________________________________________________
320 void AliPHOSv1::Reconstruction(AliPHOSReconstructioner * Reconstructioner)
321
322   // 1. Reinitializes the existing RecPoint, TrackSegment, and RecParticles Lists and 
323   // 2. Creates TreeR with a branch for each list
324   // 3. Steers the reconstruction processes
325   // 4. Saves the 3 lists in TreeR
326   // 5. Write the Tree to File
327   
328   fReconstructioner = Reconstructioner ;
329   
330   char branchname[10] ;
331   
332   // 1.
333
334   //  gAlice->MakeTree("R") ; 
335   Int_t splitlevel = 0 ; 
336   
337   fEmcRecPoints->Delete() ; 
338
339   if ( fEmcRecPoints && gAlice->TreeR() ) {
340     sprintf(branchname,"%sEmcRP",GetName()) ;
341     
342     gAlice->TreeR()->Branch(branchname, "TObjArray", &fEmcRecPoints, fBufferSize, splitlevel) ; 
343   }
344
345     fPpsdRecPoints->Delete() ; 
346
347
348   if ( fPpsdRecPoints && gAlice->TreeR() ) {
349     sprintf(branchname,"%sPpsdRP",GetName()) ;
350      
351     gAlice->TreeR()->Branch(branchname, "TObjArray", &fPpsdRecPoints, fBufferSize, splitlevel) ;
352   }
353
354   fTrackSegments->Delete() ; 
355
356   if ( fTrackSegments && gAlice->TreeR() ) { 
357     sprintf(branchname,"%sTS",GetName()) ;
358     gAlice->TreeR()->Branch(branchname, &fTrackSegments, fBufferSize) ;
359   }
360
361     fRecParticles->Delete() ; 
362
363   if ( fRecParticles && gAlice->TreeR() ) { 
364      sprintf(branchname,"%sRP",GetName()) ;
365      gAlice->TreeR()->Branch(branchname, &fRecParticles, fBufferSize) ;
366   }
367
368   
369   // 3.
370   fReconstructioner->Make(fDigits, fEmcRecPoints, fPpsdRecPoints, fTrackSegments, fRecParticles);
371
372   // 4. Expand or Shrink the arrays to the proper size
373   
374   Int_t size ;
375   
376   size = fEmcRecPoints->GetEntries() ;
377   fEmcRecPoints->Expand(size) ;
378  
379   size = fPpsdRecPoints->GetEntries() ;
380   fPpsdRecPoints->Expand(size) ;
381
382   size = fTrackSegments->GetEntries() ;
383   fTrackSegments->Expand(size) ;
384
385   size = fRecParticles->GetEntries() ;
386   fRecParticles->Expand(size) ;
387
388   gAlice->TreeR()->Fill() ;
389   // 5.
390
391   gAlice->TreeR()->Write(0,TObject::kOverwrite) ;
392  
393   // Deleting reconstructed objects
394   ResetReconstruction();
395   
396 }
397
398 //____________________________________________________________________________
399 void AliPHOSv1::ResetDigits() 
400
401   // May sound strange, but cumulative hits are store in digits Tree
402   AliDetector::ResetDigits();
403   if(  fTmpHits ) {
404     fTmpHits->Delete();
405     fNTmpHits = 0 ;
406   }
407 }  
408 //____________________________________________________________________________
409 void AliPHOSv1::ResetReconstruction() 
410
411   // Deleting reconstructed objects
412
413   if ( fEmcRecPoints )   fEmcRecPoints->Delete();
414   if ( fPpsdRecPoints )  fPpsdRecPoints->Delete();
415   if ( fTrackSegments )  fTrackSegments->Delete();
416   if ( fRecParticles )   fRecParticles->Delete();
417   
418 }
419
420 //____________________________________________________________________________
421 void AliPHOSv1::SetTreeAddress()
422
423   //  TBranch *branch;
424   AliPHOS::SetTreeAddress();
425
426  //  //Branch address for TreeR: RecPpsdRecPoint
427 //   TTree *treeR = gAlice->TreeR();
428 //   if ( treeR && fPpsdRecPoints ) {
429 //     branch = treeR->GetBranch("PHOSPpsdRP");
430 //     if (branch) branch->SetAddress(&fPpsdRecPoints) ;
431   //  }
432 }
433
434 //____________________________________________________________________________
435
436 void AliPHOSv1::StepManager(void)
437 {
438   // Accumulates hits as long as the track stays in a single crystal or PPSD gas Cell
439
440   Int_t          relid[4] ;      // (box, layer, row, column) indices
441   Float_t        xyze[4] ;       // position wrt MRS and energy deposited
442   TLorentzVector pos ;
443   Int_t copy ;
444
445   Int_t tracknumber =  gAlice->CurrentTrack() ; 
446   Int_t primary =  gAlice->GetPrimary( gAlice->CurrentTrack() ); 
447   TString name = fGeom->GetName() ; 
448   if ( name == "GPS2" ) { // the CPV is a PPSD
449     if( gMC->CurrentVolID(copy) == gMC->VolId("GCEL") ) // We are inside a gas cell 
450     {
451       gMC->TrackPosition(pos) ;
452       xyze[0] = pos[0] ;
453       xyze[1] = pos[1] ;
454       xyze[2] = pos[2] ;
455       xyze[3] = gMC->Edep() ; 
456
457       if ( xyze[3] != 0 ) { // there is deposited energy 
458         gMC->CurrentVolOffID(5, relid[0]) ;  // get the PHOS Module number
459         gMC->CurrentVolOffID(3, relid[1]) ;  // get the Micromegas Module number 
460       // 1-> Geom->GetNumberOfModulesPhi() *  fGeom->GetNumberOfModulesZ() upper                         
461       //  >  fGeom->GetNumberOfModulesPhi()  *  fGeom->GetNumberOfModulesZ() lower
462         gMC->CurrentVolOffID(1, relid[2]) ;  // get the row number of the cell
463         gMC->CurrentVolID(relid[3]) ;        // get the column number 
464
465         // get the absolute Id number
466
467         Int_t absid ; 
468         fGeom->RelToAbsNumbering(relid, absid) ; 
469
470         // add current hit to the hit list      
471         AddHit(fIshunt, primary, tracknumber, absid, xyze);
472
473       } // there is deposited energy 
474      } // We are inside the gas of the CPV  
475    } // GPS2 configuration
476   
477    if(gMC->CurrentVolID(copy) == gMC->VolId("PXTL") )  //  We are inside a PBWO crystal
478      {
479        gMC->TrackPosition(pos) ;
480        xyze[0] = pos[0] ;
481        xyze[1] = pos[1] ;
482        xyze[2] = pos[2] ;
483        xyze[3] = gMC->Edep() ;
484
485        if ( xyze[3] != 0 ) {
486           gMC->CurrentVolOffID(10, relid[0]) ; // get the PHOS module number ;
487           relid[1] = 0   ;                    // means PBW04
488           gMC->CurrentVolOffID(4, relid[2]) ; // get the row number inside the module
489           gMC->CurrentVolOffID(3, relid[3]) ; // get the cell number inside the module
490
491           // get the absolute Id number
492           
493           Int_t absid ; 
494           fGeom->RelToAbsNumbering(relid, absid) ; 
495  
496           // add current hit to the hit list
497
498           AddHit(fIshunt, primary,tracknumber, absid, xyze);
499           
500        } // there is deposited energy
501      } // we are inside a PHOS Xtal
502 }
503