source: trunk/source/processes/electromagnetic/utils/src/G4VEmModel.cc@ 968

Last change on this file since 968 was 961, checked in by garnier, 17 years ago

update processes

File size: 8.1 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// $Id: G4VEmModel.cc,v 1.25 2009/02/19 09:57:36 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name: G4VEmModel
35//
36// Author: Vladimir Ivanchenko
37//
38// Creation date: 25.07.2005
39//
40// Modifications:
41// 25.10.2005 Set default highLimit=100.TeV (V.Ivanchenko)
42// 06.02.2006 add method ComputeMeanFreePath() (mma)
43// 16.02.2009 Move implementations of virtual methods to source (VI)
44//
45//
46// Class Description:
47//
48// Abstract interface to energy loss models
49
50// -------------------------------------------------------------------
51//
52
53#include "G4VEmModel.hh"
54#include "G4LossTableManager.hh"
55#include "G4ProductionCutsTable.hh"
56
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59
60G4VEmModel::G4VEmModel(const G4String& nam):
61 fluc(0), name(nam), lowLimit(0.1*keV), highLimit(100.0*TeV),
62 polarAngleLimit(0.0),secondaryThreshold(DBL_MAX),theLPMflag(false),
63 pParticleChange(0),nuclearStopping(false),
64 currentCouple(0),currentElement(0),
65 nsec(5),flagDeexcitation(false)
66{
67 xsec.resize(nsec);
68 nSelectors = 0;
69 G4LossTableManager::Instance()->Register(this);
70}
71
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73
74G4VEmModel::~G4VEmModel()
75{
76 G4LossTableManager::Instance()->DeRegister(this);
77 G4int n = elmSelectors.size();
78 if(n > 0) {
79 for(G4int i=0; i<n; i++) {
80 delete elmSelectors[i];
81 }
82 }
83}
84
85//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
86
87void G4VEmModel::InitialiseElementSelectors(const G4ParticleDefinition* p,
88 const G4DataVector& cuts)
89{
90 // initialise before run
91 flagDeexcitation = false;
92
93 G4int nbins = G4int(std::log10(highLimit/lowLimit) + 0.5);
94 if(nbins < 3) nbins = 3;
95 G4bool spline = G4LossTableManager::Instance()->SplineFlag();
96
97 G4ProductionCutsTable* theCoupleTable=
98 G4ProductionCutsTable::GetProductionCutsTable();
99 G4int numOfCouples = theCoupleTable->GetTableSize();
100
101 // prepare vector
102 if(numOfCouples > nSelectors) {
103 elmSelectors.resize(numOfCouples);
104 nSelectors = numOfCouples;
105 }
106
107 // initialise vector
108 for(G4int i=0; i<numOfCouples; i++) {
109 const G4MaterialCutsCouple* couple =
110 theCoupleTable->GetMaterialCutsCouple(i);
111 const G4Material* material = couple->GetMaterial();
112 G4int idx = couple->GetIndex();
113
114 // selector already exist check if should be deleted
115 G4bool create = true;
116 if(elmSelectors[i]) {
117 if(material == elmSelectors[i]->GetMaterial()) create = false;
118 else delete elmSelectors[i];
119 }
120 if(create) {
121 elmSelectors[i] = new G4EmElementSelector(this,material,nbins,
122 lowLimit,highLimit,spline);
123 }
124 elmSelectors[i]->Initialise(p, cuts[idx]);
125 //elmSelectors[i]->Dump(p);
126 }
127}
128
129//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
130
131G4double G4VEmModel::ComputeDEDXPerVolume(const G4Material*,
132 const G4ParticleDefinition*,
133 G4double,G4double)
134{
135 return 0.0;
136}
137
138//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
139
140G4double G4VEmModel::CrossSectionPerVolume(const G4Material* material,
141 const G4ParticleDefinition* p,
142 G4double ekin,
143 G4double emin,
144 G4double emax)
145{
146 SetupForMaterial(p, material, ekin);
147 G4double cross = 0.0;
148 const G4ElementVector* theElementVector = material->GetElementVector();
149 const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
150 G4int nelm = material->GetNumberOfElements();
151 if(nelm > nsec) {
152 xsec.resize(nelm);
153 nsec = nelm;
154 }
155 for (G4int i=0; i<nelm; i++) {
156 cross += theAtomNumDensityVector[i]*
157 ComputeCrossSectionPerAtom(p,(*theElementVector)[i],ekin,emin,emax);
158 xsec[i] = cross;
159 }
160 return cross;
161}
162
163//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
164
165G4double G4VEmModel::ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
166 G4double, G4double, G4double,
167 G4double, G4double)
168{
169 return 0.0;
170}
171
172//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
173
174G4double G4VEmModel::MinEnergyCut(const G4ParticleDefinition*,
175 const G4MaterialCutsCouple*)
176{
177 return 0.0;
178}
179
180//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
181
182G4double G4VEmModel::GetChargeSquareRatio(const G4ParticleDefinition* p,
183 const G4Material*, G4double)
184{
185 G4double q = p->GetPDGCharge()/CLHEP::eplus;
186 return q*q;
187}
188
189//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
190
191G4double G4VEmModel::GetParticleCharge(const G4ParticleDefinition* p,
192 const G4Material*, G4double)
193{
194 return p->GetPDGCharge();
195}
196
197//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
198
199void G4VEmModel::CorrectionsAlongStep(const G4MaterialCutsCouple*,
200 const G4DynamicParticle*,
201 G4double&,G4double&,G4double)
202{}
203
204//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
205
206void G4VEmModel::SampleDeexcitationAlongStep(const G4Material*,
207 const G4Track&,
208 G4double& )
209{}
210
211//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
212
213G4double G4VEmModel::MaxSecondaryEnergy(const G4ParticleDefinition*,
214 G4double kineticEnergy)
215{
216 return kineticEnergy;
217}
218
219//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
220
221void G4VEmModel::SampleScattering(const G4DynamicParticle*, G4double)
222{}
223
224//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
225
226G4double G4VEmModel::ComputeTruePathLengthLimit(const G4Track&,
227 G4PhysicsTable*,
228 G4double)
229{
230 return DBL_MAX;
231}
232
233//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
234
235G4double G4VEmModel::ComputeGeomPathLength(G4double truePathLength)
236{
237 return truePathLength;
238}
239
240//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
241
242G4double G4VEmModel::ComputeTrueStepLength(G4double geomPathLength)
243{
244 return geomPathLength;
245}
246
247//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
248
249void G4VEmModel::DefineForRegion(const G4Region*)
250{}
251
252//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
253
254void G4VEmModel::SetupForMaterial(const G4ParticleDefinition*,
255 const G4Material*, G4double)
256{}
257
258//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.