source: trunk/source/processes/electromagnetic/standard/include/G4MultipleScattering71.hh@ 1192

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

update to geant4.9.2

File size: 8.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: G4MultipleScattering71.hh,v 1.6 2008/07/16 11:27:41 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02 $
28//
29//
30//------------- G4MultipleScattering71 physics process --------------------------
31// by Laszlo Urban, March 2001
32//
33// 07-08-01 new methods Store/Retrieve PhysicsTable
34// 23-08-01 new angle and z distribution,energy dependence reduced,
35// Store,Retrieve methods commented out temporarily, L.Urban
36// 11-09-01 G4MultipleScatteringx put as default: G4MultipleScattering
37// Store,Retrieve methods reactived (mma)
38// 13-09-01 Unused TrueToGeomTransformation method deleted,
39// class description (L.Urban)
40// 19-09-01 come back to previous process name msc
41// 17-04-02 NEW angle distribution + boundary algorithm modified, L.Urban
42// 22-04-02 boundary algorithm modified -> important improvement in timing !!!!
43// (L.Urban)
44// 24-05-02 changes in data members, L.Urban
45// 30-10-02 changes in data members, L.Urban
46// 20-01-03 Migrade to cut per region (V.Ivanchenko)
47// 05-02-03 changes in data members, L.Urban
48// 28-03-03 Move to model design (V.Ivanchenko)
49// 18-04-03 Change name (V.Ivanchenko)
50// 16-06-03: ShortLived are not applicable any more (V.Ivanchenko)
51// 17-08-04 name of data member facxsi changed to factail together
52// with the corresponding set function (L.Urban)
53// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
54// 15-04-05 optimize internal interfaces (V.Ivanchenko)
55// 03-10-05 Process is freezed with the name 71 (V.Ivanchenko)
56// 07-03-06 Create G4UrbanMscModel and move there step limit calculation (V.Ivanchenko)
57//
58//------------------------------------------------------------------------------
59//
60// $Id: G4MultipleScattering71.hh,v 1.6 2008/07/16 11:27:41 vnivanch Exp $
61// GEANT4 tag $Name: geant4-09-02 $
62
63// class description
64//
65// The class simulates the multiple scattering for any kind
66// of charged particle.
67//
68// class description - end
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72
73#ifndef G4MultipleScattering71_h
74#define G4MultipleScattering71_h 1
75
76#include "G4VMultipleScattering.hh"
77#include "G4MscModel71.hh"
78
79//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
80
81class G4MultipleScattering71 : public G4VMultipleScattering
82
83{
84public: // with description
85
86 G4MultipleScattering71(const G4String& processName="msc");
87
88 virtual ~G4MultipleScattering71();
89
90 // returns true for charged particles, false otherwise
91 G4bool IsApplicable (const G4ParticleDefinition& p);
92
93 virtual G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&);
94
95 virtual G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&);
96
97 G4double TruePathLengthLimit(const G4Track& track,
98 G4double& lambda,
99 G4double currentMinimalStep);
100
101 // Print few lines of informations about the process: validity range,
102 void PrintInfo();
103
104 // geom. step length distribution should be sampled or not
105 void Setsamplez(G4bool value);
106
107 // activate boundary algorithm
108 void SetBoundary(G4bool value);
109
110 // to reduce the energy/step dependence
111 void Setdtrl(G4double value);
112
113 void Setfactail(G4double value);
114
115 // Steplimit after boundary crossing = facrange*range
116 // estimated nb of steps at boundary nsmallstep = 1/facrange
117 void SetFacrange(G4double val);
118
119 // corrs to transport cross section for high energy
120 void SetNuclCorrPar(G4double val);
121 void SetFactPar(G4double val);
122
123protected:
124
125 // This function initialise models
126 void InitialiseProcess(const G4ParticleDefinition*);
127
128 // This method is used for tracking, it returns step limit
129 virtual G4double GetContinuousStepLimit(const G4Track& track,
130 G4double previousStepSize,
131 G4double currentMinimalStep,
132 G4double& currentSafety);
133
134private:
135
136 // hide assignment operator as private
137 G4MultipleScattering71 & operator = (const G4MultipleScattering71 &right);
138 G4MultipleScattering71 ( const G4MultipleScattering71 &);
139
140private: // data members
141
142 G4MscModel71* model;
143
144 G4double lowKineticEnergy;
145 G4double highKineticEnergy;
146 G4int totBins;
147
148 G4double truePathLength;
149 G4double geomPathLength;
150 G4double trueStepLength;
151 G4double range;
152
153 G4double facrange;
154 G4double tlimit;
155 G4double tlimitmin;
156 G4double dtrl;
157 G4double NuclCorrPar;
158 G4double FactPar;
159 G4double factail;
160 G4double cf;
161
162 G4int stepnolastmsc;
163 G4int nsmallstep;
164
165 G4bool samplez;
166 G4bool boundary;
167 G4bool isInitialized;
168
169};
170
171//--------------------------------------------------------------------
172// inline methods
173//--------------------------------------------------------------------
174
175inline G4bool G4MultipleScattering71::IsApplicable (const G4ParticleDefinition& p)
176{
177 return (p.GetPDGCharge() != 0.0 && !p.IsShortLived());
178}
179
180//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
181
182inline G4VParticleChange* G4MultipleScattering71::AlongStepDoIt(
183 const G4Track&,
184 const G4Step& step)
185{
186 G4double geomStepLength = step.GetStepLength();
187 if((geomStepLength == geomPathLength) && (truePathLength <= range))
188 trueStepLength = truePathLength;
189 else
190 trueStepLength = model->TrueStepLength(geomStepLength);
191 fParticleChange.ProposeTrueStepLength(trueStepLength);
192 return &fParticleChange;
193}
194
195//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
196
197inline G4VParticleChange* G4MultipleScattering71::PostStepDoIt(const G4Track& track,
198 const G4Step& step)
199{
200 fParticleChange.Initialize(track);
201 std::vector<G4DynamicParticle*>* p = 0;
202 model->SampleSecondaries(p, CurrentMaterialCutsCouple(),track.GetDynamicParticle(),
203 step.GetStepLength(),step.GetPostStepPoint()->GetSafety());
204 return &fParticleChange;
205}
206
207// geom. step length distribution should be sampled or not
208inline void G4MultipleScattering71::Setsamplez(G4bool value)
209{
210 samplez = value;
211}
212
213// activate boundary algorithm
214inline void G4MultipleScattering71::SetBoundary(G4bool value)
215{
216 boundary = value;
217}
218
219// to reduce the energy/step dependence
220inline void G4MultipleScattering71::Setdtrl(G4double value)
221{
222 dtrl = value;
223}
224
225inline void G4MultipleScattering71::Setfactail(G4double value)
226{
227 factail = value;
228}
229
230// Steplimit after boundary crossing = facrange*range
231// estimated nb of steps at boundary nsmallstep = 1/facrange
232inline void G4MultipleScattering71::SetFacrange(G4double val)
233{
234 facrange = val;
235 nsmallstep = G4int(std::log((cf+facrange-1.)/facrange)/std::log(cf))+1;
236}
237
238// corrs to transport cross section for high energy
239inline void G4MultipleScattering71::SetNuclCorrPar(G4double val)
240{
241 NuclCorrPar = val;
242}
243
244inline void G4MultipleScattering71::SetFactPar(G4double val)
245{
246 FactPar = val;
247}
248
249#endif
250
251
252
Note: See TracBrowser for help on using the repository browser.