]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TGeant4/TG4CutVector.cxx
Put back the two bellows in front of the absorber.
[u/mrichter/AliRoot.git] / TGeant4 / TG4CutVector.cxx
1 // $Id$
2 // Category: physics
3 //
4 // See the class description in the header file.
5
6 #include "TG4CutVector.h"
7 #include "TG4G3Defaults.h"
8
9 #include <G4Track.hh>
10 #include <G4ParticleDefinition.hh>
11 #include <G4VProcess.hh>
12
13 TG4CutVector::TG4CutVector()
14 {
15   // initialize fCutVector 
16   fCutVector = new TG4doubleVector;
17   for (G4int i=0; i<kNoG3Cuts; i++) fCutVector->insert(0.); 
18 }
19
20 TG4CutVector::TG4CutVector(const TG4CutVector& right)
21 {
22   // copy fCutVector 
23   fCutVector = new TG4doubleVector;
24   for (G4int i=0; i<kNoG3Cuts; i++) {
25     fCutVector->insert((*right.fCutVector)[i]);
26   }   
27 }
28
29 TG4CutVector::~TG4CutVector() {
30 //
31   delete fCutVector;
32 }
33
34 // operators
35
36 TG4CutVector& TG4CutVector::operator=(const TG4CutVector& right)
37 {
38   // check assignement to self
39   if (this == &right) return *this;
40
41   // initialize fCutVector 
42   fCutVector->clear();
43   for (G4int i=0; i<kNoG3Cuts; i++) {
44     fCutVector->insert((*right.fCutVector)[i]);
45   }
46   
47   return *this;   
48 }  
49
50 G4double TG4CutVector::operator[](G4int index) const 
51 {
52 //
53   if (index < kNoG3Cuts)
54     return (*fCutVector)[index];
55   else {
56     TG4Globals::Exception(
57       "TG4CutVector::operator[]: index out of the vector scope");
58     return 0.;  
59   }    
60 }  
61
62 // public methods
63
64 void TG4CutVector::SetG3Cut(TG3Cut g3Cut, G4double cutValue)
65 {
66 // Sets the cutValue for the specified cut.
67 // ---
68
69   if (g3Cut<kNoG3Cuts) {
70     (*fCutVector)[g3Cut] = cutValue;
71   }  
72   else {
73     TG4Globals::Exception(
74       "TG4CutVector::SetG3Cut: Inconsistent cut.");
75   }
76 }
77
78 void TG4CutVector::SetG3Defaults()
79 {
80 // Sets G3 default values for all cuts.
81 // ---
82
83   for (G4int i=0; i<kNoG3Cuts; i++) {
84    (*fCutVector)[i] = TG4G3Defaults::CutValue(i);
85   } 
86 }
87
88 G4double TG4CutVector::GetMinEkine(const G4Track& track) const
89 {
90 // Returns the cut value for the particle associated with
91 // specified track.
92 // ---
93
94   G4ParticleDefinition* particle = track.GetDefinition(); 
95   G4String particleName = particle->GetParticleName();
96
97   if  (particleName == "gamma") {
98     return GetMinEkineForGamma(track);
99   }    
100   else if (particleName == "e-") {    
101     return GetMinEkineForElectron(track);
102   }  
103   else if ((particle->GetParticleType() == "baryon") ||
104            (particle->GetParticleType() == "meson") ||
105            (particle->GetParticleType() == "nucleus")) { 
106     if (particle->GetPDGCharge() == 0) 
107       return GetMinEkineForHadron(track);
108     else 
109       return GetMinEkineForNeutralHadron(track);
110   }    
111   else if ((particleName == "mu-") || (particleName == "mu+")) {
112     return GetMinEkineForMuon(track);
113   }  
114   else {
115     G4String text = "TG4CutVector::GetMinEkine: \n";
116     text = text + "    The kinetic energy cut for " + particleName;
117     text = text + " is not defined.";   
118     TG4Globals::Warning(text);
119     return 0.;
120   }    
121 }
122
123 G4double TG4CutVector::GetMinEkineForGamma(const G4Track& track) const
124 {
125 // Returns the cut value for gamma.
126 // (Cut is not applied for "opticalphoton" 
127 //  as it is treated in G4 as a particle different 
128 //  from "gamma" in G4.)
129 // ---
130
131   G4cout << "TG4CutVector::GetMinEkineForGamma start" << endl;
132   const G4VProcess* kpCreatorProcess = track.GetCreatorProcess();
133   G4String processName = "";
134   if (kpCreatorProcess) processName = kpCreatorProcess->GetProcessName();
135
136   if ((processName == "eBrem") || (processName == "IeBrem")) {
137     return (*fCutVector)[kBCUTE];
138   }     
139   else if ((processName == "MuBrems") || (processName == "IMuBrems") || 
140            (processName == "//hBrems")|| (processName == "//IhBrems")) {
141            // hadron Brehmstrahlung is not defined in G4
142     return (*fCutVector)[kBCUTM];
143   }
144   else {
145     return (*fCutVector)[kCUTGAM];
146   }
147 }
148
149 G4double TG4CutVector::GetMinEkineForElectron(const G4Track& track) const
150 {
151 // Returns the cut value for e-.
152 // Should these cuts be applied to e+ too ??
153 // ---
154
155   const G4VProcess* kpCreatorProcess = track.GetCreatorProcess();
156   G4String processName = "";
157   if (kpCreatorProcess) processName = kpCreatorProcess->GetProcessName();
158
159   if ((processName == "eIoni") || (processName == "IeIoni")) {
160      // !! Geant4 treats delta rays + continuous energy loss
161      // within one process
162     return (*fCutVector)[kDCUTE];
163   }
164   else if ((processName == "MuIoni") || (processName == "IMuIoni")) {
165      // !! Geant4 treats delta rays + continuous energy loss
166      // within one process
167     return (*fCutVector)[kDCUTM];
168   }
169   else {   
170     return (*fCutVector)[kCUTELE];
171   }
172 }
173
174 G4double TG4CutVector::GetMinEkineForHadron(const G4Track& track) const
175 {
176 // Returns the cut value for charged hadron.
177 // ---
178
179   return (*fCutVector)[kCUTHAD];
180 }
181
182 G4double TG4CutVector::GetMinEkineForNeutralHadron(const G4Track& track) const
183 {
184 // Returns the cut value for neutral hadron.
185 // ---
186
187   return (*fCutVector)[kCUTNEU];
188 }
189
190 G4double TG4CutVector::GetMinEkineForMuon(const G4Track& track) const
191 {
192 // Returns the cut value for neutral muon.
193 // ---
194
195   return (*fCutVector)[kCUTMUO];
196 }
197
198 G4double TG4CutVector::GetMinEkineForOther(const G4Track& track) const
199 {
200 // Returns 0.
201 // ---
202
203   return 0.;
204 }
205