Geant4  10.02.p01
G4StepLimiter.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 //
27 // $Id: G4StepLimiter.cc 68048 2013-03-13 14:34:07Z gcosmo $
28 //
29 // --------------------------------------------------------------
30 // History
31 //
32 // 23-01-04 first implementation (H.Kurashige)
33 // --------------------------------------------------------------
34 
35 #include "G4StepLimiter.hh"
37 
38 #include "G4Step.hh"
39 #include "G4UserLimits.hh"
40 #include "G4VParticleChange.hh"
41 
42 
45  : G4VProcess(aName, fGeneral )
46 {
47  // set Process Sub Type
48  SetProcessSubType(static_cast<int>(STEP_LIMITER));
49 
50  if (verboseLevel>1) {
51  G4cout << GetProcessName() << " is created "<< G4endl;
52  }
53 }
54 
55 
58 {
59 }
60 
61 
64  : G4VProcess(right)
65 {
66 }
67 
68 
70 G4double
72  const G4Track& aTrack,
73  G4double, // previousStepSize
75 {
76  // condition is set to "Not Forced"
77  *condition = NotForced;
78 
79  G4double proposedStep = DBL_MAX;
80  G4UserLimits* pUserLimits =
82  if (pUserLimits) {
83  // max allowed step length
84  //
85  proposedStep = pUserLimits->GetMaxAllowedStep(aTrack);
86  if (proposedStep < 0.) proposedStep =0.;
87  }
88  return proposedStep;
89 }
90 
94  const G4Step& )
95 // Do Nothing
96 //
97 {
99  return &aParticleChange;
100 }
101 
102 
103 
104 
105 
106 
107 
108 
109 
110 
111 
112 
113 
G4double condition(const G4ErrorSymMatrix &m)
G4int verboseLevel
Definition: G4VProcess.hh:368
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4GLOB_DLL std::ostream G4cout
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4StepLimiter(const G4String &processName="StepLimiter")
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
Definition: G4Step.hh:76
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4UserLimits * GetUserLimits() const
virtual void Initialize(const G4Track &)
G4LogicalVolume * GetLogicalVolume() const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
G4VPhysicalVolume * GetVolume() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
virtual G4double GetMaxAllowedStep(const G4Track &)
virtual ~G4StepLimiter()
G4ForceCondition
#define DBL_MAX
Definition: templates.hh:83