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

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

update

File size: 5.6 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.20 2008/11/13 23:13:18 schaelic Exp $
27// GEANT4 tag $Name: geant4-09-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//
44//
45// Class Description:
46//
47// Abstract interface to energy loss models
48
49// -------------------------------------------------------------------
50//
51
52#include "G4VEmModel.hh"
53#include "G4LossTableManager.hh"
54#include "G4ProductionCutsTable.hh"
55
56//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58
59G4VEmModel::G4VEmModel(const G4String& nam):
60 fluc(0), name(nam), lowLimit(0.1*keV), highLimit(100.0*TeV),
61 polarAngleLimit(0.0),secondaryThreshold(DBL_MAX),theLPMflag(false),
62 pParticleChange(0),nuclearStopping(false),nsec(5)
63{
64 xsec.resize(nsec);
65 nSelectors = 0;
66 G4LossTableManager::Instance()->Register(this);
67}
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
70
71G4VEmModel::~G4VEmModel()
72{
73 G4LossTableManager::Instance()->DeRegister(this);
74 G4int n = elmSelectors.size();
75 if(n > 0) {
76 for(G4int i=0; i<n; i++) {
77 delete elmSelectors[i];
78 }
79 }
80}
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
83
84G4double G4VEmModel::CrossSectionPerVolume(const G4Material* material,
85 const G4ParticleDefinition* p,
86 G4double ekin,
87 G4double emin,
88 G4double emax)
89{
90 SetupForMaterial(p, material, ekin);
91 G4double cross = 0.0;
92 const G4ElementVector* theElementVector = material->GetElementVector();
93 const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
94 G4int nelm = material->GetNumberOfElements();
95 if(nelm > nsec) {
96 xsec.resize(nelm);
97 nsec = nelm;
98 }
99 for (G4int i=0; i<nelm; i++) {
100 cross += theAtomNumDensityVector[i]*
101 ComputeCrossSectionPerAtom(p,(*theElementVector)[i],ekin,emin,emax);
102 xsec[i] = cross;
103 }
104 return cross;
105}
106
107//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
108
109G4double G4VEmModel::ComputeMeanFreePath(const G4ParticleDefinition* p,
110 G4double ekin,
111 const G4Material* material,
112 G4double emin,
113 G4double emax)
114{
115 G4double mfp = DBL_MAX;
116 G4double cross = CrossSectionPerVolume(material,p,ekin,emin,emax);
117 if (cross > DBL_MIN) mfp = 1./cross;
118 return mfp;
119}
120
121//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
122
123void G4VEmModel::InitialiseElementSelectors(const G4ParticleDefinition* p,
124 const G4DataVector& cuts)
125{
126 G4int nbins = G4int(std::log10(highLimit/lowLimit) + 0.5);
127 if(nbins < 3) nbins = 3;
128 G4bool spline = G4LossTableManager::Instance()->SplineFlag();
129
130 G4ProductionCutsTable* theCoupleTable=
131 G4ProductionCutsTable::GetProductionCutsTable();
132 G4int numOfCouples = theCoupleTable->GetTableSize();
133
134 // prepare vector
135 if(numOfCouples > nSelectors) {
136 elmSelectors.resize(numOfCouples);
137 nSelectors = numOfCouples;
138 }
139
140 // initialise vector
141 for(G4int i=0; i<numOfCouples; i++) {
142 const G4MaterialCutsCouple* couple =
143 theCoupleTable->GetMaterialCutsCouple(i);
144 const G4Material* material = couple->GetMaterial();
145 G4int idx = couple->GetIndex();
146
147 // selector already exist check if should be deleted
148 G4bool create = true;
149 if(elmSelectors[i]) {
150 if(material == elmSelectors[i]->GetMaterial()) create = false;
151 else delete elmSelectors[i];
152 }
153 if(create) {
154 elmSelectors[i] = new G4EmElementSelector(this,material,nbins,
155 lowLimit,highLimit,spline);
156 }
157 elmSelectors[i]->Initialise(p, cuts[idx]);
158 //elmSelectors[i]->Dump(p);
159 }
160}
161
162//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
163
164
Note: See TracBrowser for help on using the repository browser.