Geant4
v4-10.4-release
메인 페이지
관련된 페이지
모듈
네임스페이스
클래스
파일들
파일 목록
파일 멤버
모두
클래스
네임스페이스들
파일들
함수
변수
타입정의
열거형 타입
열거형 멤버
Friends
매크로
그룹들
페이지들
source
processes
transportation
src
G4NeutronKiller.cc
이 파일의 문서화 페이지로 가기
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.cc 110694 2018-06-08 06:39:34Z gcosmo $
27
//
28
//---------------------------------------------------------------------------
29
//
30
// ClassName: G4NeutronKiller
31
//
32
// Description: The process to kill particles to save CPU
33
//
34
// Author: V.Ivanchenko 26/09/00 for HARP software
35
//
36
//----------------------------------------------------------------------------
37
//
38
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
39
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
40
41
#include "
G4NeutronKiller.hh
"
42
43
#include "
G4SystemOfUnits.hh
"
44
#include "
G4NeutronKillerMessenger.hh
"
45
#include "
G4TransportationProcessType.hh
"
46
47
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48
49
G4NeutronKiller::G4NeutronKiller
(
const
G4String
& processName,
G4ProcessType
aType)
50
:
G4VDiscreteProcess
(processName, aType)
51
{
52
// set Process Sub Type
53
SetProcessSubType
(
NEUTRON_KILLER
);
54
55
kinEnergyThreshold
= 0.0;
56
timeThreshold
=
DBL_MAX
;
57
pMess
=
new
G4NeutronKillerMessenger
(
this
);
58
}
59
60
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
61
62
G4NeutronKiller::~G4NeutronKiller
()
63
{
64
delete
pMess
;
65
}
66
67
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
68
69
G4bool
G4NeutronKiller::IsApplicable
(
const
G4ParticleDefinition
& particle)
70
{
71
return
(particle.
GetParticleName
() ==
"neutron"
);
72
}
73
74
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
75
76
void
G4NeutronKiller::SetTimeLimit
(
G4double
val)
77
{
78
timeThreshold
= val;
79
if
(
verboseLevel
> 0)
80
G4cout
<<
"### G4NeutronKiller: timeLimit(ns) = "
81
<<
timeThreshold
/
ns
<<
G4endl
;
82
}
83
84
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
85
86
void
G4NeutronKiller::SetKinEnergyLimit
(
G4double
val)
87
{
88
kinEnergyThreshold
= val;
89
if
(
verboseLevel
> 0)
90
G4cout
<<
"### G4NeutronKiller: Tracking cut E(MeV) = "
91
<<
kinEnergyThreshold
/
MeV
<<
G4endl
;
92
}
93
94
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
95
96
G4double
G4NeutronKiller::PostStepGetPhysicalInteractionLength
(
97
const
G4Track
& aTrack,
98
G4double
,
G4ForceCondition
*
condition
)
99
{
100
// condition is set to "Not Forced"
101
*condition =
NotForced
;
102
103
return
(aTrack.
GetGlobalTime
() >
timeThreshold
||
104
aTrack.
GetKineticEnergy
() <
kinEnergyThreshold
) ? 0.0 :
DBL_MAX
;
105
}
106
107
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
108
109
G4double
G4NeutronKiller::GetMeanFreePath
(
const
G4Track
&,
G4double
,
110
G4ForceCondition
*)
111
{
112
return
DBL_MAX
;
113
}
114
115
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
116
117
G4VParticleChange
*
G4NeutronKiller::PostStepDoIt
(
const
G4Track
& aTrack,
118
const
G4Step
&)
119
{
120
pParticleChange
->
Initialize
(aTrack);
121
pParticleChange
->
ProposeTrackStatus
(
fStopAndKill
);
122
return
pParticleChange
;
123
}
124
125
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4Track::GetKineticEnergy
G4double GetKineticEnergy() const
G4NeutronKiller::~G4NeutronKiller
virtual ~G4NeutronKiller()
Definition:
G4NeutronKiller.cc:62
MeV
static constexpr double MeV
Definition:
G4SIunits.hh:214
G4VDiscreteProcess
Definition:
G4VDiscreteProcess.hh:58
G4endl
#define G4endl
Definition:
G4ios.hh:61
G4ParticleDefinition::GetParticleName
const G4String & GetParticleName() const
Definition:
G4ParticleDefinition.hh:121
fStopAndKill
Definition:
G4TrackStatus.hh:56
G4String
Definition:
examples/extended/parallel/TopC/ParN02/AnnotatedFiles/G4String.hh:45
G4NeutronKillerMessenger.hh
G4VParticleChange
Definition:
G4VParticleChange.hh:94
condition
G4double condition(const G4ErrorSymMatrix &m)
G4NeutronKiller.hh
G4NeutronKiller::GetMeanFreePath
virtual G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *)
Definition:
G4NeutronKiller.cc:109
G4VParticleChange::Initialize
virtual void Initialize(const G4Track &)
G4TransportationProcessType.hh
G4double
double G4double
Definition:
G4Types.hh:76
G4bool
bool G4bool
Definition:
G4Types.hh:79
G4VProcess::pParticleChange
G4VParticleChange * pParticleChange
Definition:
G4VProcess.hh:283
G4NeutronKiller::G4NeutronKiller
G4NeutronKiller(const G4String &processName="nKiller", G4ProcessType aType=fGeneral)
Definition:
G4NeutronKiller.cc:49
G4Track
Definition:
G4Track.hh:76
G4NeutronKiller::IsApplicable
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition:
G4NeutronKiller.cc:69
G4ProcessType
G4ProcessType
Definition:
G4ProcessType.hh:43
G4NeutronKiller::kinEnergyThreshold
G4double kinEnergyThreshold
Definition:
G4NeutronKiller.hh:98
NotForced
Definition:
G4ForceCondition.hh:56
G4Track::GetGlobalTime
G4double GetGlobalTime() const
G4NeutronKiller::pMess
G4NeutronKillerMessenger * pMess
Definition:
G4NeutronKiller.hh:101
G4SystemOfUnits.hh
G4ParticleDefinition
Definition:
G4ParticleDefinition.hh:73
G4NeutronKiller::PostStepGetPhysicalInteractionLength
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
Definition:
G4NeutronKiller.cc:96
G4Step
Definition:
G4Step.hh:76
G4NeutronKiller::SetTimeLimit
void SetTimeLimit(G4double)
Definition:
G4NeutronKiller.cc:76
G4NeutronKiller::SetKinEnergyLimit
void SetKinEnergyLimit(G4double)
Definition:
G4NeutronKiller.cc:86
G4VProcess::verboseLevel
G4int verboseLevel
Definition:
G4VProcess.hh:371
G4NeutronKiller::timeThreshold
G4double timeThreshold
Definition:
G4NeutronKiller.hh:99
G4VProcess::SetProcessSubType
void SetProcessSubType(G4int)
Definition:
G4VProcess.hh:435
G4ForceCondition
G4ForceCondition
Definition:
G4ForceCondition.hh:49
G4cout
G4GLOB_DLL std::ostream G4cout
NEUTRON_KILLER
Definition:
G4TransportationProcessType.hh:50
G4NeutronKillerMessenger
Definition:
G4NeutronKillerMessenger.hh:53
DBL_MAX
#define DBL_MAX
Definition:
templates.hh:83
G4VParticleChange::ProposeTrackStatus
void ProposeTrackStatus(G4TrackStatus status)
G4NeutronKiller::PostStepDoIt
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
Definition:
G4NeutronKiller.cc:117
ns
#define ns
Definition:
xmlparse.cc:614
다음에 의해 생성됨 :
1.8.5