Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4LowEWentzelVIModel.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: G4LowEWentzelVIModel.cc 74528 2013-10-12 17:24:24Z vnivanch $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4LowEWentzelVIModel
34 //
35 // Author: V.Ivanchenko
36 //
37 // Creation date: 11.02.2014 from G4WentzelVIModel
38 //
39 // Modifications:
40 //
41 // Class Description:
42 //
43 
44 // -------------------------------------------------------------------
45 //
46 
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
49 
50 #include "G4LowEWentzelVIModel.hh"
51 #include "G4PhysicalConstants.hh"
52 #include "G4SystemOfUnits.hh"
53 
54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55 
56 using namespace std;
57 
59  G4WentzelVIModel(false,"LowEnWentzelVI")
60 {
62 }
63 
64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
65 
67 {}
68 
69 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
70 
72  const G4Track& track,
73  G4double& currentMinimalStep)
74 {
75  G4double tlimit = currentMinimalStep;
76  const G4DynamicParticle* dp = track.GetDynamicParticle();
77  G4StepPoint* sp = track.GetStep()->GetPreStepPoint();
78  G4StepStatus stepStatus = sp->GetStepStatus();
79  singleScatteringMode = false;
80  //G4cout << "G4LowEWentzelVIModel::ComputeTruePathLengthLimit stepStatus= "
81  // << stepStatus << " " << track.GetDefinition()->GetParticleName()
82  // << G4endl;
83 
84  // initialisation for each step, lambda may be computed from scratch
90  /*
91  G4cout << "lambdaeff= " << lambdaeff << " Range= " << currentRange
92  << " tlimit= " << tlimit << " 1-cost= " << 1 - cosTetMaxNuc << G4endl;
93  */
94  // extra check for abnormal situation
95  // this check needed to run MSC with eIoni and eBrem inactivated
96  tlimit = std::min(tlimit, currentRange);
97 
98  // stop here if small range particle
99  if(tlimit < tlimitminfix) {
100  return ConvertTrueToGeom(tlimit, currentMinimalStep);
101  }
102 
103  // pre step
104  G4double presafety = sp->GetSafety();
105  // far from geometry boundary
106  if(currentRange < presafety) {
107  return ConvertTrueToGeom(tlimit, currentMinimalStep);
108  }
109 
110  // compute presafety again if presafety <= 0 and no boundary
111  // i.e. when it is needed for optimization purposes
112  if(stepStatus != fGeomBoundary && presafety < tlimitminfix) {
113  presafety = ComputeSafety(sp->GetPosition(), tlimit);
114  if(currentRange < presafety) {
115  return ConvertTrueToGeom(tlimit, currentMinimalStep);
116  }
117  }
118  /*
119  G4cout << "e(MeV)= " << preKinEnergy/MeV
120  << " " << particle->GetParticleName()
121  << " CurLimit(mm)= " << tlimit/mm <<" safety(mm)= " << presafety/mm
122  << " R(mm)= " <<currentRange/mm
123  << " L0(mm^-1)= " << lambdaeff*mm
124  <<G4endl;
125  */
126  // natural limit for high energy
127  G4double rlimit = std::max(facrange*currentRange, lambdaeff);
128  //G4double rlimit = std::max(facrange*currentRange,
129  // 0.7*(1.0 - cosTetMaxNuc)*lambdaeff);
130 
131  // low-energy e-
132  rlimit = std::max(rlimit, facsafety*presafety);
133 
134  // cut correction
135  //G4double rcut = currentCouple->GetProductionCuts()->GetProductionCut(1);
136  //G4cout << "rcut= " << rcut << " rlimit= " << rlimit << " presafety= "
137  // << presafety << " 1-cosThetaMax= " <<1-cosThetaMax
138  //<< " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc << G4endl;
139  //if(rcut > rlimit) { rlimit = std::min(rlimit, rcut*sqrt(rlimit/rcut)); }
140 
141  tlimit = std::min(tlimit, rlimit);
142  tlimit = std::max(tlimit, tlimitminfix);
143 
144  // step limit in infinite media
145  tlimit = std::min(tlimit, 50*currentMaterial->GetRadlen()/facgeom);
146 
147  //compute geomlimit and force few steps within a volume
149  && stepStatus == fGeomBoundary) {
150 
151  G4double geomlimit = ComputeGeomLimit(track, presafety, currentRange);
152  tlimit = std::min(tlimit, geomlimit/facgeom);
153  }
154  /*
155  G4cout << particle->GetParticleName() << " E(MeV)= " << preKinEnergy
156  << " L0= " << lambdaeff << " R= " << currentRange
157  << " tlimit= " << tlimit
158  << " currentMinimalStep= " << currentMinimalStep << G4endl;
159  */
160  return ConvertTrueToGeom(tlimit, currentMinimalStep);
161 }
162 
163 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
void SetSingleScatteringFactor(G4double)
G4double facgeom
Definition: G4VMscModel.hh:176
G4WentzelOKandVIxSection * wokvi
G4double GetKineticEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
G4double facrange
Definition: G4VMscModel.hh:175
G4StepStatus GetStepStatus() const
void DefineMaterial(const G4MaterialCutsCouple *)
G4double ConvertTrueToGeom(G4double &tLength, G4double &gLength)
Definition: G4VMscModel.hh:245
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
const G4MaterialCutsCouple * currentCouple
const G4Material * currentMaterial
const G4Step * GetStep() const
G4double ComputeSafety(const G4ThreeVector &position, G4double limit=DBL_MAX)
Definition: G4VMscModel.hh:237
G4StepStatus
Definition: G4StepStatus.hh:49
G4StepPoint * GetPreStepPoint() const
const G4ThreeVector & GetPosition() const
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:283
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
G4double ComputeGeomLimit(const G4Track &, G4double &presafety, G4double limit)
Definition: G4VMscModel.hh:255
G4double GetTransportMeanFreePath(const G4ParticleDefinition *part, G4double kinEnergy)
Definition: G4VMscModel.hh:340
G4double GetRadlen() const
Definition: G4Material.hh:220
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double facsafety
Definition: G4VMscModel.hh:177
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
const G4ParticleDefinition * particle
G4double GetSafety() const
G4MscStepLimitType steppingAlgorithm
Definition: G4VMscModel.hh:185
double G4double
Definition: G4Types.hh:76
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep)