source: trunk/source/processes/electromagnetic/utils/src/G4VMscModel.cc@ 1181

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

maj sur la beta de geant 4.9.3

File size: 5.5 KB
RevLine 
[968]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//
[1055]26// $Id: G4VMscModel.cc,v 1.12 2009/05/10 19:36:19 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03-beta-cand-01 $
[968]28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name: G4VMscModel
35//
36// Author: Vladimir Ivanchenko
37//
38// Creation date: 07.03.2008
39//
40// Modifications:
41//
42//
43// Class Description:
44//
45// General interface to msc models
46
47// -------------------------------------------------------------------
48//
49
50#include "G4VMscModel.hh"
[1055]51#include "G4ParticleChangeForMSC.hh"
52#include "G4TransportationManager.hh"
[968]53
54//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56
57G4VMscModel::G4VMscModel(const G4String& nam):
58 G4VEmModel(nam),
[1055]59 safetyHelper(0),
60 facrange(0.04),
[968]61 facgeom(2.5),
62 facsafety(0.25),
63 skin(3.0),
64 dtrl(0.05),
65 lambdalimit(mm),
[1055]66 geommax(1.e50*mm),
[968]67 steppingAlgorithm(fUseSafety),
68 samplez(false),
69 latDisplasment(true)
70{}
71
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73
74G4VMscModel::~G4VMscModel()
75{}
76
77//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
[1055]78
79G4ParticleChangeForMSC* G4VMscModel::GetParticleChangeForMSC()
80{
81 G4ParticleChangeForMSC* p = 0;
82 if (pParticleChange) {
83 p = static_cast<G4ParticleChangeForMSC*>(pParticleChange);
84 } else {
85 p = new G4ParticleChangeForMSC();
86 }
87 return p;
88}
89
90//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
91
92void G4VMscModel::InitialiseSafetyHelper()
93{
94 if(!safetyHelper) {
95 safetyHelper = G4TransportationManager::GetTransportationManager()
96 ->GetSafetyHelper();
97 safetyHelper->InitialiseHelper();
98 }
99}
100
101//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
102
103void G4VMscModel::ComputeDisplacement(G4ParticleChangeForMSC* fParticleChange,
104 const G4ThreeVector& dir,
105 G4double displacement,
106 G4double postsafety)
107{
108 const G4ThreeVector* pos = fParticleChange->GetProposedPosition();
109 G4double r = displacement;
110 if(r > postsafety) {
111 G4double newsafety = safetyHelper->ComputeSafety(*pos);
112 if(r > newsafety) r = newsafety;
113 }
114 if(r > 0.) {
115
116 // compute new endpoint of the Step
117 G4ThreeVector newPosition = *pos + r*dir;
118
119 // definitely not on boundary
120 if(displacement == r) {
121 safetyHelper->ReLocateWithinVolume(newPosition);
122
123 } else {
124 // check safety after displacement
125 G4double postsafety = safetyHelper->ComputeSafety(newPosition);
126
127 // displacement to boundary
128 if(postsafety <= 0.0) {
129 safetyHelper->Locate(newPosition,
130 *fParticleChange->GetProposedMomentumDirection());
131
132 // not on the boundary
133 } else {
134 safetyHelper->ReLocateWithinVolume(newPosition);
135 }
136 }
137 fParticleChange->ProposePosition(newPosition);
138 }
139}
140
141//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
142
143void G4VMscModel::SampleScattering(const G4DynamicParticle*, G4double)
144{}
145
146//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
147
148G4double G4VMscModel::ComputeTruePathLengthLimit(const G4Track&,
149 G4PhysicsTable*,
150 G4double)
151{
152 return DBL_MAX;
153}
154
155//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
156
157G4double G4VMscModel::ComputeGeomPathLength(G4double truePathLength)
158{
159 return truePathLength;
160}
161
162//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
163
164G4double G4VMscModel::ComputeTrueStepLength(G4double geomPathLength)
165{
166 return geomPathLength;
167}
168
169//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
170
171void G4VMscModel::SampleSecondaries(std::vector<G4DynamicParticle*>*,
172 const G4MaterialCutsCouple*,
173 const G4DynamicParticle*,
174 G4double, G4double)
175{}
176
177//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.