Geant4-11
G4EmMultiModel.cc
Go to the documentation of this file.
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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4EmMultiModel
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 03.05.2004
37//
38// Modifications:
39// 15-04-05 optimize internal interface (V.Ivanchenko)
40// 04-07-10 updated interfaces according to g4 9.4 (V.Ivanchenko)
41//
42
43
44//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
45//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46
47#include "G4EmMultiModel.hh"
48#include "Randomize.hh"
49#include "G4EmParameters.hh"
50
51//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
52
54 : G4VEmModel(nam)
55{
56 model.clear();
57 cross_section.clear();
58}
59
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
61
63{}
64
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66
68{
69 cross_section.push_back(0.0);
70 model.push_back(p);
71 ++nModels;
72}
73
74//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
75
77 const G4DataVector& cuts)
78{
80 G4int verb = IsMaster() ? param->Verbose() : param->WorkerVerbose();
81 if(verb > 0) {
82 G4cout << "### Initialisation of EM MultiModel " << GetName()
83 << " including following list of " << nModels << " models:" << G4endl;
84 }
85 for(G4int i=0; i<nModels; ++i) {
86 G4cout << " " << (model[i])->GetName();
88 (model[i])->Initialise(p, cuts);
89 }
90 if(verb > 0) { G4cout << G4endl; }
91}
92
93//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
94
96 const G4ParticleDefinition* p,
97 G4double kineticEnergy,
98 G4double cutEnergy)
99{
100 G4double dedx = 0.0;
101 for(G4int i=0; i<nModels; ++i) {
102 dedx += (model[i])->ComputeDEDXPerVolume(mat, p, cutEnergy, kineticEnergy);
103 }
104
105 return dedx;
106}
107
108//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
109
111 G4double kinEnergy,
112 G4double Z,
113 G4double A,
114 G4double cutEnergy,
115 G4double maxEnergy)
116{
117 G4double cross = 0.0;
118 for(G4int i=0; i<nModels; ++i) {
120 cross += (model[i])->ComputeCrossSectionPerAtom(p, kinEnergy, Z, A,
121 cutEnergy, maxEnergy);
122 }
123 return cross;
124}
125
126//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
127
128void G4EmMultiModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
129 const G4MaterialCutsCouple* couple,
130 const G4DynamicParticle* dp,
131 G4double minEnergy,
132 G4double maxEnergy)
133{
134 SetCurrentCouple(couple);
135 if(nModels > 0) {
136 G4int i;
137 G4double cross = 0.0;
138 for(i=0; i<nModels; ++i) {
139 cross += (model[i])->CrossSection(couple, dp->GetParticleDefinition(),
140 dp->GetKineticEnergy(), minEnergy, maxEnergy);
141 cross_section[i] = cross;
142 }
143
144 cross *= G4UniformRand();
145
146 for(i=0; i<nModels; ++i) {
147 if(cross <= cross_section[i]) {
148 (model[i])->SampleSecondaries(vdp, couple, dp, minEnergy, maxEnergy);
149 return;
150 }
151 }
152 }
153}
154
155//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
156
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
G4EmMultiModel(const G4String &nam="MultiModel")
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double tmax) final
std::vector< G4VEmModel * > model
std::vector< G4double > cross_section
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) final
~G4EmMultiModel() override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) final
void Initialise(const G4ParticleDefinition *, const G4DataVector &) final
void AddModel(G4VEmModel *)
static G4EmParameters * Instance()
G4int Verbose() const
G4int WorkerVerbose() const
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:614
G4bool IsMaster() const
Definition: G4VEmModel.hh:746
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:472
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:447
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:425
G4double CrossSection(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:539
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:490
const G4String & GetName() const
Definition: G4VEmModel.hh:837