Geant4  10.02.p01
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 // $Id: G4EmMultiModel.cc 92921 2015-09-21 15:06:51Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4EmMultiModel
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 03.05.2004
38 //
39 // Modifications:
40 // 15-04-05 optimize internal interface (V.Ivanchenko)
41 // 04-07-10 updated interfaces according to g4 9.4 (V.Ivanchenko)
42 //
43 
44 
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47 
48 #include "G4EmMultiModel.hh"
49 #include "Randomize.hh"
50 
51 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
52 
54  : G4VEmModel(nam), nModels(0)
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 {
79  if(nModels > 0) {
80  G4cout << "### Initialisation of EM MultiModel " << GetName()
81  << " including following list of models:" << G4endl;
82  for(G4int i=0; i<nModels; ++i) {
83  G4cout << " " << (model[i])->GetName();
85  (model[i])->Initialise(p, cuts);
86  }
87  G4cout << G4endl;
88  }
89 }
90 
91 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
92 
94  const G4ParticleDefinition* p,
95  G4double kineticEnergy,
96  G4double cutEnergy)
97 {
98  SetCurrentCouple(couple);
99  G4double dedx = 0.0;
100 
101  if(nModels > 0) {
102  for(G4int i=0; i<nModels; i++) {
103  dedx += (model[i])->ComputeDEDX(couple, p, cutEnergy, kineticEnergy);
104  }
105  }
106 
107  return dedx;
108 }
109 
110 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
111 
113  G4double kinEnergy,
114  G4double Z,
115  G4double A,
116  G4double cutEnergy,
117  G4double maxEnergy)
118 {
119  G4double cross = 0.0;
120  if(nModels>0) {
121  for(G4int i=0; i<nModels; ++i) {
123  cross += (model[i])->ComputeCrossSectionPerAtom(p, kinEnergy, Z, A,
124  cutEnergy, maxEnergy);
125  }
126  }
127  return cross;
128 }
129 
130 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
131 
132 void G4EmMultiModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
133  const G4MaterialCutsCouple* couple,
134  const G4DynamicParticle* dp,
135  G4double minEnergy,
136  G4double maxEnergy)
137 {
138  SetCurrentCouple(couple);
139  if(nModels > 0) {
140  G4int i;
141  G4double cross = 0.0;
142  for(i=0; i<nModels; ++i) {
143  cross += (model[i])->CrossSection(couple, dp->GetParticleDefinition(),
144  dp->GetKineticEnergy(), minEnergy, maxEnergy);
145  cross_section[i] = cross;
146  }
147 
148  cross *= G4UniformRand();
149 
150  for(i=0; i<nModels; ++i) {
151  if(cross <= cross_section[i]) {
152  (model[i])->SampleSecondaries(vdp, couple, dp, minEnergy, maxEnergy);
153  return;
154  }
155  }
156  }
157 }
158 
159 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
160 
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double tmax)
void AddModel(G4VEmModel *)
virtual ~G4EmMultiModel()
G4double GetKineticEnergy() const
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:610
int G4int
Definition: G4Types.hh:78
G4EmMultiModel(const G4String &nam="MultiModel")
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:452
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=0)
Definition: G4VEmModel.cc:410
G4double CrossSection(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:501
virtual G4double ComputeDEDX(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
const G4ParticleDefinition * GetParticleDefinition() const
std::vector< G4double > cross_section
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:445
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:422
#define G4endl
Definition: G4ios.hh:61
const G4String & GetName() const
Definition: G4VEmModel.hh:795
double G4double
Definition: G4Types.hh:76
std::vector< G4VEmModel * > model