]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALv1.cxx
Reduced the number of daughter volumes (Sahal)
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALv1.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 EMCAL Manager class 
20 // An object of this class does not produce digits
21 // It is the one to use if you do want to produce outputs in TREEH 
22 //                  
23 //*-- Author: Sahal Yacoob (LBL /UCT)
24 //*--       : Jennifer Klay (LBL)
25
26
27 // --- ROOT system ---
28 #include "TPGON.h"
29 #include "TTUBS.h"
30 #include "TNode.h"
31 #include "TRandom.h"
32 #include "TTree.h"
33 #include "TGeometry.h"
34
35
36 // --- Standard library ---
37
38 #include <stdio.h>
39 #include <string.h>
40 #include <stdlib.h>
41 #include <strstream.h>
42 #include <iostream.h>
43 #include <math.h>
44 // --- AliRoot header files ---
45
46 #include "AliEMCALv1.h"
47 #include "AliEMCALHit.h"
48 #include "AliEMCALGeometry.h"
49 #include "AliConst.h"
50 #include "AliRun.h"
51 #include "AliMC.h"
52
53 ClassImp(AliEMCALv1)
54
55
56 //______________________________________________________________________
57 AliEMCALv1::AliEMCALv1():AliEMCALv0(){
58   // ctor
59 }
60 //______________________________________________________________________
61 AliEMCALv1::AliEMCALv1(const char *name, const char *title):
62     AliEMCALv0(name,title){
63     // Standard Creator.
64
65     fHits= new TClonesArray("AliEMCALHit",1000);
66     gAlice->AddHitList(fHits);
67
68     fNhits = 0;
69     fSamplingFraction = 12.9 ; 
70     fLayerToPreshowerRatio = 5.0/6.0; 
71     fIshunt     =  1; // All hits are associated with primary particles
72 }
73 //______________________________________________________________________
74 AliEMCALv1::~AliEMCALv1(){
75     // dtor
76
77     if ( fHits) {
78         fHits->Delete();
79         delete fHits;
80         fHits = 0;
81     }
82 }
83 //______________________________________________________________________
84 void AliEMCALv1::AddHit(Int_t shunt, Int_t primary, Int_t tracknumber, Int_t iparent, Float_t ienergy, 
85                         Int_t id, Float_t * hits,Float_t * p){
86     // Add a hit to the hit list.
87     // An EMCAL hit is the sum of all hits in a single segment 
88     //   originating from the same enterring particle 
89     Int_t hitCounter;
90     
91     AliEMCALHit *newHit;
92     AliEMCALHit *curHit;
93     Bool_t deja = kFALSE;
94
95     newHit = new AliEMCALHit(shunt, primary, tracknumber, iparent, ienergy, id, hits, p);
96     for ( hitCounter = fNhits-1; hitCounter >= 0 && !deja; hitCounter-- ) {
97         curHit = (AliEMCALHit*) (*fHits)[hitCounter];
98         // We add hits with the same tracknumber, while GEANT treats
99         // primaries succesively
100         if(curHit->GetPrimary() != primary) break;
101         if( *curHit == *newHit ) {
102             *curHit = *curHit + *newHit;
103             deja = kTRUE;
104         } // end if
105     } // end for hitCounter
106
107     if ( !deja ) {
108         new((*fHits)[fNhits]) AliEMCALHit(*newHit);
109         fNhits++;
110     } // end if
111
112     delete newHit;
113 }
114 //______________________________________________________________________
115 void AliEMCALv1::StepManager(void){
116     // Accumulates hits as long as the track stays in a single
117     // crystal or PPSD gas Cell
118     Int_t          id[2];           // (layer, phi, Eta) indices
119     Int_t          absid;
120     // position wrt MRS and energy deposited
121     Float_t        xyze[4]={0.,0.,0.,0.};
122     Float_t        pmom[4]={0.,0.,0.,0.};
123     TLorentzVector pos; // Lorentz vector of the track current position.
124     TLorentzVector mom; // Lorentz vector of the track current momentum.
125     Int_t tracknumber =  gAlice->CurrentTrack();
126     Int_t primary = 0;
127     static Int_t iparent = 0;
128     static Float_t ienergy = 0;
129     Int_t copy = 0;
130
131
132  if(gMC->IsTrackEntering() && (strcmp(gMC->CurrentVolName(),"XALU") == 0)) // This Particle in enterring the Calorimeter 
133  {
134     iparent = tracknumber;
135     gMC->TrackMomentum(mom);
136     ienergy = mom[3]; 
137 }
138  if(gMC->CurrentVolID(copy) == gMC->VolId("XPHI") ) { // We are in a Scintillator Layer 
139
140     gMC->CurrentVolOffID(1, id[0]); // get the POLY copy number;
141     gMC->CurrentVolID(id[1]); // get the phi number inside the layer
142     primary = gAlice->GetPrimary(tracknumber);
143     gMC->TrackPosition(pos);
144     gMC->TrackMomentum(mom);
145     xyze[0] = pos[0];
146     xyze[1] = pos[1];
147     xyze[2] = pos[2];
148     xyze[3] = fSamplingFraction*(gMC->Edep()); // Correct for sampling calorimeter
149     pmom[0] = mom[0];
150     pmom[1] = mom[1];
151     pmom[2] = mom[2];
152     pmom[3] = mom[3];
153     
154     if(xyze[3] > 0.){// Track is inside the crystal and deposits some energy
155         absid = (id[0]-1)*(fGeom->GetNPhi()) + id[1];
156         if((absid/fGeom->GetNPhi()) < (2*fGeom->GetNZ()))
157         {xyze[3] = fLayerToPreshowerRatio*xyze[3] ;}  //                                             Preshower readout must be scaled
158         AddHit(fIshunt, primary,tracknumber, iparent, ienergy, absid, xyze, pmom);
159     } // there is deposited energy
160 }
161 }