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