Geant4
v4-10.4-release
메인 페이지
관련된 페이지
모듈
네임스페이스
클래스
파일들
파일 목록
파일 멤버
모두
클래스
네임스페이스들
파일들
함수
변수
타입정의
열거형 타입
열거형 멤버
Friends
매크로
그룹들
페이지들
source
processes
electromagnetic
adjoint
src
G4eAdjointMultipleScattering.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: G4eAdjointMultipleScattering.cc 67990 2013-03-13 10:56:28Z gcosmo $
27
//
28
// -----------------------------------------------------------------------------
29
//
30
// GEANT4 Class file
31
//
32
// File name: G4eAdjointMultipleScattering
33
//
34
// Author: Vladimir Ivanchenko
35
//
36
// Creation date: 10 March 2008
37
//
38
// Modifications:
39
//
40
// -----------------------------------------------------------------------------
41
//
42
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
43
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
44
45
#include "
G4eAdjointMultipleScattering.hh
"
46
#include "
G4UrbanAdjointMscModel.hh
"
47
#include "
G4MscStepLimitType.hh
"
48
#include "
G4Electron.hh
"
49
#include "
G4Positron.hh
"
50
#include "
G4DynamicParticle.hh
"
51
52
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
53
54
using namespace
std;
55
56
G4eAdjointMultipleScattering::G4eAdjointMultipleScattering
(
const
G4String
& processName)
57
:
G4VMultipleScattering
(processName)
58
{
59
isInitialized
=
false
;
60
}
61
62
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
63
64
G4eAdjointMultipleScattering::~G4eAdjointMultipleScattering
()
65
{}
66
67
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
68
69
G4bool
G4eAdjointMultipleScattering::IsApplicable
(
const
G4ParticleDefinition
&
p
)
70
{
71
return
(p.
GetPDGCharge
() != 0.0 && !p.
IsShortLived
());
72
}
73
74
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
75
76
void
G4eAdjointMultipleScattering::InitialiseProcess
(
const
G4ParticleDefinition
*)
77
{
78
if
(
isInitialized
) {
return
; }
79
if
(!
EmModel
(0)) {
SetEmModel
(
new
G4UrbanAdjointMscModel
(), 0); }
80
AddEmModel
(1,
EmModel
(0));
81
isInitialized
=
true
;
82
}
83
84
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
85
86
void
G4eAdjointMultipleScattering::PrintInfo
()
87
{
88
G4cout
<<
" RangeFactor= "
<<
RangeFactor
()
89
<<
", stepLimitType: "
<<
StepLimitType
()
90
<<
", latDisplacement: "
<<
LateralDisplasmentFlag
();
91
if
(
StepLimitType
() ==
fUseDistanceToBoundary
) {
92
G4cout
<<
", skin= "
<<
Skin
() <<
", geomFactor= "
<<
GeomFactor
();
93
}
94
G4cout
<<
G4endl
;
95
}
96
97
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
98
99
void
G4eAdjointMultipleScattering::StartTracking
(
G4Track
* )
100
{
G4DynamicParticle
* aDynPart =
new
G4DynamicParticle
(
G4Electron::Electron
(),
G4ThreeVector
(0.,0.,1.),1.);
101
G4Track
* tempTrack =
new
G4Track
(aDynPart,0.,
G4ThreeVector
(0.,0.,0.));
102
G4VMultipleScattering::StartTracking
( tempTrack);
103
delete
tempTrack;
104
}
G4UrbanAdjointMscModel
Definition:
G4UrbanAdjointMscModel.hh:74
G4ThreeVector
CLHEP::Hep3Vector G4ThreeVector
Definition:
G4ThreeVector.hh:42
G4Positron.hh
G4Electron.hh
G4VMultipleScattering
Definition:
G4VMultipleScattering.hh:92
G4UrbanAdjointMscModel.hh
G4VMultipleScattering::StepLimitType
G4MscStepLimitType StepLimitType() const
Definition:
G4VMultipleScattering.hh:378
G4endl
#define G4endl
Definition:
G4ios.hh:61
p
const char * p
Definition:
xmltok.h:285
G4String
Definition:
examples/extended/parallel/TopC/ParN02/AnnotatedFiles/G4String.hh:45
G4ParticleDefinition::GetPDGCharge
G4double GetPDGCharge() const
Definition:
G4ParticleDefinition.hh:125
G4DynamicParticle
Definition:
G4DynamicParticle.hh:73
G4VMultipleScattering::RangeFactor
G4double RangeFactor() const
Definition:
G4VMultipleScattering.hh:347
G4VMultipleScattering::GeomFactor
G4double GeomFactor() const
Definition:
G4VMultipleScattering.hh:364
G4eAdjointMultipleScattering::PrintInfo
void PrintInfo() override
Definition:
G4eAdjointMultipleScattering.cc:86
G4VMultipleScattering::LateralDisplasmentFlag
G4bool LateralDisplasmentFlag() const
Definition:
G4VMultipleScattering.hh:319
G4bool
bool G4bool
Definition:
G4Types.hh:79
G4DynamicParticle.hh
G4Track
Definition:
G4Track.hh:76
G4eAdjointMultipleScattering::InitialiseProcess
void InitialiseProcess(const G4ParticleDefinition *) override
Definition:
G4eAdjointMultipleScattering.cc:76
G4eAdjointMultipleScattering::StartTracking
void StartTracking(G4Track *) override
Definition:
G4eAdjointMultipleScattering.cc:99
G4VMultipleScattering::EmModel
G4VMscModel * EmModel(size_t index=0) const
Definition:
G4VMultipleScattering.cc:162
G4ParticleDefinition
Definition:
G4ParticleDefinition.hh:73
G4VMultipleScattering::StartTracking
void StartTracking(G4Track *) override
Definition:
G4VMultipleScattering.cc:366
G4eAdjointMultipleScattering.hh
G4MscStepLimitType.hh
G4Electron::Electron
static G4Electron * Electron()
Definition:
G4Electron.cc:94
G4VMultipleScattering::SetEmModel
void SetEmModel(G4VMscModel *, size_t index=0)
Definition:
G4VMultipleScattering.cc:154
G4eAdjointMultipleScattering::IsApplicable
G4bool IsApplicable(const G4ParticleDefinition &p) final
Definition:
G4eAdjointMultipleScattering.cc:69
fUseDistanceToBoundary
Definition:
G4MscStepLimitType.hh:52
G4eAdjointMultipleScattering::~G4eAdjointMultipleScattering
virtual ~G4eAdjointMultipleScattering()
Definition:
G4eAdjointMultipleScattering.cc:64
G4cout
G4GLOB_DLL std::ostream G4cout
G4VMultipleScattering::AddEmModel
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=nullptr)
Definition:
G4VMultipleScattering.cc:144
G4eAdjointMultipleScattering::G4eAdjointMultipleScattering
G4eAdjointMultipleScattering(const G4String &processName="msc")
Definition:
G4eAdjointMultipleScattering.cc:56
G4ParticleDefinition::IsShortLived
G4bool IsShortLived() const
Definition:
G4ParticleDefinition.hh:158
G4eAdjointMultipleScattering::isInitialized
G4bool isInitialized
Definition:
G4eAdjointMultipleScattering.hh:85
G4VMultipleScattering::Skin
G4double Skin() const
Definition:
G4VMultipleScattering.hh:333
다음에 의해 생성됨 :
1.8.5