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

Last change on this file since 1199 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 8.3 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.30 2009/09/23 14:42:47 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
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#include "G4ParticleChangeForLoss.hh"
57#include "G4ParticleChangeForGamma.hh"
58
59//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
61
62G4VEmModel::G4VEmModel(const G4String& nam):
63 fluc(0), name(nam), lowLimit(0.1*keV), highLimit(100.0*TeV),
64 eMinActive(0.0),eMaxActive(DBL_MAX),
65 polarAngleLimit(0.0),secondaryThreshold(DBL_MAX),theLPMflag(false),
66 pParticleChange(0),nuclearStopping(false),
67 currentCouple(0),currentElement(0),
68 nsec(5),flagDeexcitation(false)
69{
70 xsec.resize(nsec);
71 nSelectors = 0;
72 G4LossTableManager::Instance()->Register(this);
73}
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
76
77G4VEmModel::~G4VEmModel()
78{
79 G4LossTableManager::Instance()->DeRegister(this);
80 G4int n = elmSelectors.size();
81 if(n > 0) {
82 for(G4int i=0; i<n; i++) {
83 delete elmSelectors[i];
84 }
85 }
86}
87
88//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
89
90G4ParticleChangeForLoss* G4VEmModel::GetParticleChangeForLoss()
91{
92 G4ParticleChangeForLoss* p = 0;
93 if (pParticleChange) {
94 p = static_cast<G4ParticleChangeForLoss*>(pParticleChange);
95 } else {
96 p = new G4ParticleChangeForLoss();
97 pParticleChange = p;
98 }
99 return p;
100}
101
102//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
103
104G4ParticleChangeForGamma* G4VEmModel::GetParticleChangeForGamma()
105{
106 G4ParticleChangeForGamma* p = 0;
107 if (pParticleChange) {
108 p = static_cast<G4ParticleChangeForGamma*>(pParticleChange);
109 } else {
110 p = new G4ParticleChangeForGamma();
111 pParticleChange = p;
112 }
113 return p;
114}
115
116//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
117
118void G4VEmModel::InitialiseElementSelectors(const G4ParticleDefinition* p,
119 const G4DataVector& cuts)
120{
121 // initialise before run
122 flagDeexcitation = false;
123
124 G4int nbins = G4int(std::log10(highLimit/lowLimit) + 0.5);
125 if(nbins < 3) nbins = 3;
126 G4bool spline = G4LossTableManager::Instance()->SplineFlag();
127
128 G4ProductionCutsTable* theCoupleTable=
129 G4ProductionCutsTable::GetProductionCutsTable();
130 G4int numOfCouples = theCoupleTable->GetTableSize();
131
132 // prepare vector
133 if(numOfCouples > nSelectors) elmSelectors.reserve(numOfCouples);
134
135 // initialise vector
136 for(G4int i=0; i<numOfCouples; i++) {
137 currentCouple = theCoupleTable->GetMaterialCutsCouple(i);
138 const G4Material* material = currentCouple->GetMaterial();
139 G4int idx = currentCouple->GetIndex();
140
141 // selector already exist check if should be deleted
142 G4bool create = true;
143 if(i < nSelectors) {
144 if(elmSelectors[i]) {
145 if(material == elmSelectors[i]->GetMaterial()) create = false;
146 else delete elmSelectors[i];
147 }
148 } else {
149 nSelectors++;
150 elmSelectors.push_back(0);
151 }
152 if(create) {
153 elmSelectors[i] = new G4EmElementSelector(this,material,nbins,
154 lowLimit,highLimit,spline);
155 }
156 elmSelectors[i]->Initialise(p, cuts[idx]);
157 //elmSelectors[i]->Dump(p);
158 }
159}
160
161//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
162
163G4double G4VEmModel::ComputeDEDXPerVolume(const G4Material*,
164 const G4ParticleDefinition*,
165 G4double,G4double)
166{
167 return 0.0;
168}
169
170//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
171
172G4double G4VEmModel::CrossSectionPerVolume(const G4Material* material,
173 const G4ParticleDefinition* p,
174 G4double ekin,
175 G4double emin,
176 G4double emax)
177{
178 SetupForMaterial(p, material, ekin);
179 G4double cross = 0.0;
180 const G4ElementVector* theElementVector = material->GetElementVector();
181 const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
182 G4int nelm = material->GetNumberOfElements();
183 if(nelm > nsec) {
184 xsec.resize(nelm);
185 nsec = nelm;
186 }
187 for (G4int i=0; i<nelm; i++) {
188 cross += theAtomNumDensityVector[i]*
189 ComputeCrossSectionPerAtom(p,(*theElementVector)[i],ekin,emin,emax);
190 xsec[i] = cross;
191 }
192 return cross;
193}
194
195//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
196
197G4double G4VEmModel::ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
198 G4double, G4double, G4double,
199 G4double, G4double)
200{
201 return 0.0;
202}
203
204//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
205
206void G4VEmModel::DefineForRegion(const G4Region*)
207{}
208
209//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
210
211G4double G4VEmModel::MinEnergyCut(const G4ParticleDefinition*,
212 const G4MaterialCutsCouple*)
213{
214 return 0.0;
215}
216
217//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
218
219G4double G4VEmModel::GetChargeSquareRatio(const G4ParticleDefinition* p,
220 const G4Material*, G4double)
221{
222 G4double q = p->GetPDGCharge()/CLHEP::eplus;
223 return q*q;
224}
225
226//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
227
228G4double G4VEmModel::GetParticleCharge(const G4ParticleDefinition* p,
229 const G4Material*, G4double)
230{
231 return p->GetPDGCharge();
232}
233
234//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
235
236void G4VEmModel::CorrectionsAlongStep(const G4MaterialCutsCouple*,
237 const G4DynamicParticle*,
238 G4double&,G4double&,G4double)
239{}
240
241//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
242
243void G4VEmModel::SampleDeexcitationAlongStep(const G4Material*,
244 const G4Track&,
245 G4double& )
246{}
247
248//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
249
250G4double G4VEmModel::MaxSecondaryEnergy(const G4ParticleDefinition*,
251 G4double kineticEnergy)
252{
253 return kineticEnergy;
254}
255
256//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
257
258void G4VEmModel::SetupForMaterial(const G4ParticleDefinition*,
259 const G4Material*, G4double)
260{}
261
262//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.