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