Geant4  10.02.p01
G4NeutronKiller.hh
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: G4NeutronKiller.hh 68048 2013-03-13 14:34:07Z gcosmo $
27 //
28 //---------------------------------------------------------------------------
29 //
30 // ClassName: G4NeutronKiller
31 //
32 // Description: The process to kill neutrons to save CPU
33 //
34 // Author: V.Ivanchenko 26/09/00 for HARP software
35 //
36 //----------------------------------------------------------------------------
37 //
38 // Class description:
39 //
40 // G4NeutronKiller allows to remove unwanted neutrons from simulation in
41 // order to improve CPU performance. There are two parameters:
42 //
43 // 1) low energy threshold for neutron kinetic energy (default 0)
44 // 2) time limit for neutron track (default DBL_MAX)
45 //
46 // These parameters can be changed via Set methods or by UI commands:
47 // /physics_engine/neutron/energyCut
48 // /physics_engine/neutron/timeLimit
49 //
50 // If a neutron track is killed no energy deposition is added to the step
51 //
52 
53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55 
56 #ifndef NeutronKiller_h
57 #define NeutronKiller_h 1
58 
59 #include "globals.hh"
60 #include "G4VDiscreteProcess.hh"
61 #include "G4ParticleDefinition.hh"
62 #include "G4Step.hh"
63 #include "G4Track.hh"
64 
66 
67 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
68 
70 {
71 public:
72 
73  G4NeutronKiller(const G4String& processName = "nKiller",
74  G4ProcessType aType = fGeneral );
75 
76  virtual ~G4NeutronKiller();
77 
79 
80  void SetTimeLimit(G4double);
81 
83 
85  G4double previousStepSize,
87 
89 
91 
92 private:
93 
94  // hide assignment operator as private
97 
100 
102 };
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
105 
107  const G4Track& aTrack,
109 {
110  // condition is set to "Not Forced"
111  *condition = NotForced;
112 
113  G4double limit = DBL_MAX;
114  if(aTrack.GetGlobalTime() > timeThreshold ||
115  aTrack.GetKineticEnergy() < kinEnergyThreshold) limit = 0.0;
116  return limit;
117 }
118 
119 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
120 
123 {
124  return DBL_MAX;
125 }
126 
127 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
128 
130  const G4Step&)
131 {
132  pParticleChange->Initialize(aTrack);
134  return pParticleChange;
135 }
136 
137 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
138 
139 #endif
140 
G4double condition(const G4ErrorSymMatrix &m)
virtual void Initialize(const G4Track &)
G4double timeThreshold
void SetKinEnergyLimit(G4double)
virtual ~G4NeutronKiller()
G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *)
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4bool IsApplicable(const G4ParticleDefinition &)
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4double GetKineticEnergy() const
bool G4bool
Definition: G4Types.hh:79
G4NeutronKiller(const G4String &processName="nKiller", G4ProcessType aType=fGeneral)
Definition: G4Step.hh:76
G4double GetGlobalTime() const
void SetTimeLimit(G4double)
G4NeutronKiller & operator=(const G4NeutronKiller &right)
G4double kinEnergyThreshold
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
double G4double
Definition: G4Types.hh:76
G4NeutronKillerMessenger * pMess
void ProposeTrackStatus(G4TrackStatus status)
G4ForceCondition
#define DBL_MAX
Definition: templates.hh:83
G4ProcessType