Geant4
v4-10.4-release
메인 페이지
관련된 페이지
모듈
네임스페이스
클래스
파일들
파일 목록
파일 멤버
모두
클래스
네임스페이스들
파일들
함수
변수
타입정의
열거형 타입
열거형 멤버
Friends
매크로
그룹들
페이지들
examples
extended
biasing
GB06
src
GB06BOptrSplitAndKillByImportance.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
//
28
//
29
#include "
GB06BOptrSplitAndKillByImportance.hh
"
30
#include "
GB06BOptnSplitAndKillByImportance.hh
"
31
32
#include "
G4BiasingProcessInterface.hh
"
33
#include "
G4ParallelGeometriesLimiterProcess.hh
"
34
#include "
G4BiasingProcessSharedData.hh
"
35
36
37
#include "
G4ParticleDefinition.hh
"
38
#include "
G4ParticleTable.hh
"
39
#include "
G4VProcess.hh
"
40
41
#include "
G4TransportationManager.hh
"
42
#include "
G4TouchableHistoryHandle.hh
"
43
44
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
45
46
GB06BOptrSplitAndKillByImportance::
47
GB06BOptrSplitAndKillByImportance
(
G4String
particleName,
48
G4String
name
)
49
:
G4VBiasingOperator
(name),
50
fParallelWorld ( nullptr ),
51
fParallelWorldIndex ( -1 ),
52
fBiasingLimiterProcess ( nullptr )
53
{
54
fParticleToBias
=
G4ParticleTable::GetParticleTable
()->
FindParticle
(particleName);
55
56
if
(
fParticleToBias
== 0 )
57
{
58
G4ExceptionDescription
ed;
59
ed <<
"Particle `"
<< particleName <<
"' not found !"
<<
G4endl
;
60
G4Exception
(
"GB06BOptrSplitAndKillByImportance(...)"
,
61
"exGB06.01"
,
62
JustWarning
,
63
ed);
64
}
65
66
fSplitAndKillByImportance
=
67
new
GB06BOptnSplitAndKillByImportance
(
"splitterFor_"
+particleName);
68
69
}
70
71
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72
73
GB06BOptrSplitAndKillByImportance::~GB06BOptrSplitAndKillByImportance
()
74
{
75
delete
fSplitAndKillByImportance
;
76
}
77
78
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
79
80
void
GB06BOptrSplitAndKillByImportance::StartRun
()
81
{
82
// ---------------
83
// -- Setup stage:
84
// ---------------
85
// -- get the particle process manager...
86
const
G4ProcessManager
* processManager =
fParticleToBias
->
GetProcessManager
();
87
// -- ... to obtain the biasing information shared among this particle processes:
88
fBiasingSharedData
=
G4BiasingProcessInterface::GetSharedData
( processManager );
89
// -- Remember the index of the parallel world:
90
fBiasingLimiterProcess
=
fBiasingSharedData
->
GetParallelGeometriesLimiterProcess
();
91
fParallelWorldIndex
=
fBiasingLimiterProcess
->
GetParallelWorldIndex
(
fParallelWorld
);
92
93
// -- Setup the biasing operation:
94
fSplitAndKillByImportance
-> SetBiasingSharedData(
fBiasingSharedData
);
95
fSplitAndKillByImportance
->
SetParallelWorldIndex
(
fParallelWorldIndex
);
96
fSplitAndKillByImportance
-> SetImportanceMap( &
fImportanceMap
);
97
98
}
99
100
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
101
102
G4VBiasingOperation
*
103
GB06BOptrSplitAndKillByImportance::
104
ProposeNonPhysicsBiasingOperation
(
const
G4Track
*
track
,
105
const
G4BiasingProcessInterface
*
/*callingProcess*/
)
106
{
107
// -- Check if current particle type is the one to bias:
108
if
( track->
GetDefinition
() !=
fParticleToBias
)
return
0;
109
110
// -- if so, request biasing:
111
return
fSplitAndKillByImportance
;
112
113
}
114
G4ParticleTable::FindParticle
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
Definition:
G4ParticleTable.cc:577
G4BiasingProcessInterface
Definition:
G4BiasingProcessInterface.hh:70
name
const XML_Char * name
Definition:
expat.h:151
GB06BOptrSplitAndKillByImportance::~GB06BOptrSplitAndKillByImportance
virtual ~GB06BOptrSplitAndKillByImportance()
Definition:
GB06BOptrSplitAndKillByImportance.cc:73
G4ExceptionDescription
std::ostringstream G4ExceptionDescription
Definition:
G4Exception.hh:45
GB06BOptrSplitAndKillByImportance::fImportanceMap
std::map< G4int, G4int > fImportanceMap
Definition:
GB06BOptrSplitAndKillByImportance.hh:115
G4ParticleTable::GetParticleTable
static G4ParticleTable * GetParticleTable()
Definition:
G4ParticleTable.cc:105
G4BiasingProcessInterface.hh
G4ParallelGeometriesLimiterProcess::GetParallelWorldIndex
G4int GetParallelWorldIndex(const G4VPhysicalVolume *parallelWorld) const
Definition:
G4ParallelGeometriesLimiterProcess.cc:407
G4BiasingProcessSharedData::GetParallelGeometriesLimiterProcess
const G4ParallelGeometriesLimiterProcess * GetParallelGeometriesLimiterProcess() const
Definition:
G4BiasingProcessSharedData.hh:76
G4endl
#define G4endl
Definition:
G4ios.hh:61
GB06BOptrSplitAndKillByImportance::fParallelWorldIndex
G4int fParallelWorldIndex
Definition:
GB06BOptrSplitAndKillByImportance.hh:112
G4String
Definition:
examples/extended/parallel/TopC/ParN02/AnnotatedFiles/G4String.hh:45
G4VProcess.hh
GB06BOptrSplitAndKillByImportance::GB06BOptrSplitAndKillByImportance
GB06BOptrSplitAndKillByImportance(G4String particleToBias, G4String name="SplitAndKillByImportance")
Definition:
GB06BOptrSplitAndKillByImportance.cc:47
G4Track::GetDefinition
G4ParticleDefinition * GetDefinition() const
GB06BOptrSplitAndKillByImportance.hh
Definition of the GB06BOptrSplitAndKillByImportance class.
GB06BOptrSplitAndKillByImportance::ProposeNonPhysicsBiasingOperation
virtual G4VBiasingOperation * ProposeNonPhysicsBiasingOperation(const G4Track *, const G4BiasingProcessInterface *) final
Definition:
GB06BOptrSplitAndKillByImportance.cc:104
GB06BOptnSplitAndKillByImportance::SetParallelWorldIndex
void SetParallelWorldIndex(G4int parallelWorldIndex)
Definition:
GB06BOptnSplitAndKillByImportance.hh:87
G4TransportationManager.hh
JustWarning
Definition:
G4ExceptionSeverity.hh:64
GB06BOptrSplitAndKillByImportance::fBiasingSharedData
const G4BiasingProcessSharedData * fBiasingSharedData
Definition:
GB06BOptrSplitAndKillByImportance.hh:113
GB06BOptrSplitAndKillByImportance::StartRun
virtual void StartRun()
Definition:
GB06BOptrSplitAndKillByImportance.cc:80
G4Track
Definition:
G4Track.hh:76
G4BiasingProcessSharedData.hh
GB06BOptrSplitAndKillByImportance::fParticleToBias
const G4ParticleDefinition * fParticleToBias
Definition:
GB06BOptrSplitAndKillByImportance.hh:110
G4ParticleTable.hh
G4VBiasingOperator
Definition:
G4VBiasingOperator.hh:181
GB06BOptnSplitAndKillByImportance.hh
Definition of the GB06BOptnSplitAndKillByImportance class.
G4VBiasingOperation
Definition:
G4VBiasingOperation.hh:77
G4Exception
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition:
G4Exception.hh:65
G4ParticleDefinition::GetProcessManager
G4ProcessManager * GetProcessManager() const
Definition:
G4ParticleDefinition.cc:258
GB06BOptrSplitAndKillByImportance::fBiasingLimiterProcess
const G4ParallelGeometriesLimiterProcess * fBiasingLimiterProcess
Definition:
GB06BOptrSplitAndKillByImportance.hh:114
GB06BOptrSplitAndKillByImportance::fSplitAndKillByImportance
GB06BOptnSplitAndKillByImportance * fSplitAndKillByImportance
Definition:
GB06BOptrSplitAndKillByImportance.hh:109
G4ParticleDefinition.hh
GB06BOptnSplitAndKillByImportance
Definition:
GB06BOptnSplitAndKillByImportance.hh:42
G4TouchableHistoryHandle.hh
G4BiasingProcessInterface::GetSharedData
const G4BiasingProcessSharedData * GetSharedData() const
Definition:
G4BiasingProcessInterface.hh:115
track
Float_t track
Definition:
extended/medical/dna/range/plot.C:35
GB06BOptrSplitAndKillByImportance::fParallelWorld
G4VPhysicalVolume * fParallelWorld
Definition:
GB06BOptrSplitAndKillByImportance.hh:111
G4ParallelGeometriesLimiterProcess.hh
G4ProcessManager
Definition:
G4ProcessManager.hh:106
다음에 의해 생성됨 :
1.8.5