Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4ParticleChangeForGamma.hh
이 파일의 문서화 페이지로 가기
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: G4ParticleChangeForGamma.hh 68795 2013-04-05 13:24:46Z gcosmo $
28 //
29 //
30 // ------------------------------------------------------------
31 // GEANT 4 class header file
32 //
33 //
34 // ------------------------------------------------------------
35 // 15 April 2005 V.Ivanchenko for gamma EM processes
36 //
37 // Modified:
38 // 30.05.05 : add UpdateStepForAtRest (V.Ivanchenko)
39 // 04.12.05 : apply UpdateStepForPostStep in any case (mma)
40 // 26.08.06 : Add->Set polarization;
41 // add const method to access track;
42 // add weight modification (V.Ivanchenko)
43 //
44 // ------------------------------------------------------------
45 //
46 // Class Description
47 // This class is a concrete class for ParticleChange for gamma processes
48 //
49 #ifndef G4ParticleChangeForGamma_h
50 #define G4ParticleChangeForGamma_h 1
51 
52 #include "globals.hh"
53 #include "G4ios.hh"
54 #include "G4VParticleChange.hh"
55 
56 class G4DynamicParticle;
57 
59 {
60 public:
61  // default constructor
63 
64  // destructor
65  virtual ~G4ParticleChangeForGamma();
66 
67  // with description
68  // ----------------------------------------------------
69  // --- the following methods are for updating G4Step -----
70 
73  // A physics process gives the final state of the particle
74  // based on information of G4Track
75 
76  inline void InitializeForPostStep(const G4Track&);
77  //Initialize all propoerties by using G4Track information
78 
79  void AddSecondary(G4DynamicParticle* aParticle);
80  // Add next secondary
81 
82  inline G4double GetProposedKineticEnergy() const;
84  // Get/Set the final kinetic energy of the current particle.
85 
86  inline const G4ThreeVector& GetProposedMomentumDirection() const;
87  inline void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz);
88  inline void ProposeMomentumDirection(const G4ThreeVector& Pfinal);
89  // Get/Propose the MomentumDirection vector: it is the final momentum direction.
90 
91  inline const G4ThreeVector& GetProposedPolarization() const;
92  inline void ProposePolarization(const G4ThreeVector& dir);
93  inline void ProposePolarization(G4double Px, G4double Py, G4double Pz);
94 
95  inline const G4Track* GetCurrentTrack() const;
96 
97  virtual void DumpInfo() const;
98 
99  // for Debug
100  virtual G4bool CheckIt(const G4Track&);
101 
102 protected:
103  // hide copy constructor and assignment operaor as protected
106 
107 private:
108 
110  // The pointer to G4Track
111 
113  // The final kinetic energy of the current particle.
114 
116  // The final momentum direction of the current particle.
117 
119  // The final polarization of the current particle.
120 };
121 
122 // ------------------------------------------------------------
123 
125 {
126  return proposedKinEnergy;
127 }
128 
130 {
132 }
133 
134 inline
136 {
138 }
139 
140 inline
142 {
144 }
145 
146 inline
148 {
152 }
153 
155 {
156  return currentTrack;
157 }
158 
159 inline
161 {
162  return proposedPolarization;
163 }
164 
165 inline
167 {
169 }
170 
171 inline
173 {
177 }
178 
180 {
182  theLocalEnergyDeposit = 0.0;
184  InitializeSecondaries(track);
185  theParentWeight = track.GetWeight();
186  isParentWeightProposed = false;
190  currentTrack = &track;
191 }
192 
193 #endif
194 
const G4ThreeVector & GetProposedPolarization() const
G4double GetKineticEnergy() const
void AddSecondary(G4DynamicParticle *aParticle)
const G4ThreeVector & GetPolarization() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposePolarization(const G4ThreeVector &dir)
G4Step * UpdateStepForPostStep(G4Step *Step)
void setX(double)
virtual G4bool CheckIt(const G4Track &)
G4double GetWeight() const
G4ParticleChangeForGamma & operator=(const G4ParticleChangeForGamma &right)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4double theNonIonizingEnergyDeposit
G4TrackStatus theStatusChange
void setZ(double)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
double energy
Definition: plottest35.C:25
void InitializeSecondaries(const G4Track &)
const G4Track * GetCurrentTrack() const
Definition: G4Step.hh:76
void InitializeForPostStep(const G4Track &)
G4double GetProposedKineticEnergy() const
TDirectory * dir
const G4ThreeVector & GetProposedMomentumDirection() const
G4Step * UpdateStepForAtRest(G4Step *pStep)
G4TrackStatus GetTrackStatus() const
const G4ThreeVector & GetMomentumDirection() const
Definition: Step.hh:41
void setY(double)