source: trunk/source/processes/electromagnetic/utils/src/G4EmMultiModel.cc

Last change on this file was 1340, checked in by garnier, 14 years ago

update ti head

File size: 5.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: G4EmMultiModel.cc,v 1.8 2010/08/17 17:36:59 vnivanch Exp $
27// GEANT4 tag $Name: emutils-V09-03-23 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name:   G4EmMultiModel
35//
36// Author:        Vladimir Ivanchenko
37//
38// Creation date: 03.05.2004
39//
40// Modifications:
41// 15-04-05 optimize internal interface (V.Ivanchenko)
42// 04-07-10 updated interfaces according to g4 9.4 (V.Ivanchenko)
43//
44
45
46//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48
49#include "G4EmMultiModel.hh"
50#include "Randomize.hh"
51
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
53
54G4EmMultiModel::G4EmMultiModel(const G4String& nam)
55  : G4VEmModel(nam), nModels(0)
56{
57  model.clear();
58  cross_section.clear();
59}
60
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62
63G4EmMultiModel::~G4EmMultiModel()
64{}
65
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
67
68void G4EmMultiModel::AddModel(G4VEmModel* p)
69{
70  cross_section.push_back(0.0);
71  model.push_back(p);
72  ++nModels;
73}
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
76
77void G4EmMultiModel::Initialise(const G4ParticleDefinition* p, 
78                                const G4DataVector& cuts)
79{
80  if(nModels > 0) {
81    G4cout << "### Initialisation of EM MultiModel " << GetName()
82           << " including following list of models:" << G4endl;
83    for(G4int i=0; i<nModels; ++i) {
84      G4cout << "    " << (model[i])->GetName();
85      (model[i])->SetParticleChange(pParticleChange, GetModelOfFluctuations());
86      (model[i])->Initialise(p, cuts);
87    }
88    G4cout << G4endl;
89  }
90}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
93
94G4double G4EmMultiModel::ComputeDEDX(const G4MaterialCutsCouple* couple,
95                                     const G4ParticleDefinition* p,
96                                     G4double kineticEnergy,
97                                     G4double cutEnergy)
98{
99  SetCurrentCouple(couple);
100  G4double dedx = 0.0;
101
102  if(nModels > 0) {
103    for(G4int i=0; i<nModels; i++) {
104      dedx += (model[i])->ComputeDEDX(couple, p, cutEnergy, kineticEnergy);
105    }
106  } 
107
108  return dedx;
109}
110
111//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
112
113G4double G4EmMultiModel::ComputeCrossSectionPerAtom(const G4ParticleDefinition* p,
114                                                    G4double kinEnergy,
115                                                    G4double Z,
116                                                    G4double A,
117                                                    G4double cutEnergy,
118                                                    G4double maxEnergy)
119{
120  G4double cross = 0.0;
121  if(nModels>0) {
122    for(G4int i=0; i<nModels; ++i) {
123      (model[i])->SetCurrentCouple(CurrentCouple());
124      cross += (model[i])->ComputeCrossSectionPerAtom(p, kinEnergy, Z, A, 
125                                                      cutEnergy, maxEnergy);
126    }
127  } 
128  return cross;
129}
130
131//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
132
133void G4EmMultiModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
134                                       const G4MaterialCutsCouple* couple,
135                                       const G4DynamicParticle* dp,
136                                       G4double minEnergy,
137                                       G4double maxEnergy)
138{
139  SetCurrentCouple(couple);
140  if(nModels > 0) {
141    G4int i;
142    G4double cross = 0.0;
143    for(i=0; i<nModels; ++i) {
144      cross += (model[i])->CrossSection(couple, dp->GetParticleDefinition(), 
145                                        dp->GetKineticEnergy(), minEnergy, maxEnergy);
146      cross_section[i] = cross;
147    }
148
149    cross *= G4UniformRand();
150
151    for(i=0; i<nModels; ++i) {
152      if(cross <= cross_section[i]) {
153        (model[i])->SampleSecondaries(vdp, couple, dp, minEnergy, maxEnergy);
154        return;
155      }
156    }
157  } 
158}
159
160//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
161
Note: See TracBrowser for help on using the repository browser.