Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4UnknownDecay.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: G4UnknownDecay.hh 105727 2017-08-16 12:47:05Z gcosmo $
28 //
29 //
30 // ------------------------------------------------------------
31 // GEANT 4 class header file
32 //
33 //
34 #ifndef G4UnknownDecay_h
35 #define G4UnknownDecay_h 1
36 
38 
39 #include "G4ios.hh"
40 #include "globals.hh"
41 #include "G4VDiscreteProcess.hh"
43 
45 {
46  // Class Description
47  // This class is a decay process for "unknown" particle.
48 
49  public:
50  // Constructors
51  G4UnknownDecay(const G4String& processName ="UnknownDecay");
52 
53  // Destructor
54  virtual ~G4UnknownDecay();
55 
56  private:
57  // copy constructor
59 
60  // Assignment Operation (generated)
61  G4UnknownDecay & operator=(const G4UnknownDecay &right);
62 
63  public: //With Description
64 
66  const G4Track& aTrack,
67  const G4Step& aStep
68  ) override;
69 
70  virtual void BuildPhysicsTable(const G4ParticleDefinition&) override;
71  // In G4UnknownDecay, thePhysicsTable stores values of
72  // beta * std::sqrt( 1 - beta*beta)
73  // as a function of normalized kinetic enregy (=Ekin/mass),
74  // becasuse this table is universal for all particle types,
75 
76 
77  virtual G4bool IsApplicable(const G4ParticleDefinition&) override;
78  // returns "true" if the decay process can be applied to
79  // the particle type.
80 
81  virtual void ProcessDescription(std::ostream& outFile) const override;
82  //
83 
84  protected: // With Description
85  virtual G4VParticleChange* DecayIt(
86  const G4Track& aTrack,
87  const G4Step& aStep
88  ) ;
89  // The DecayIt() method returns by pointer a particle-change object,
90  // which has information of daughter particles.
91 
92  public:
93 
95  const G4Track& track,
96  G4double previousStepSize,
98  ) override;
99 
100 
101  protected: // With Description
102  // GetMeanFreePath returns ctau*beta*gamma for decay in flight
103  virtual G4double GetMeanFreePath(const G4Track& aTrack,
104  G4double previousStepSize,
105  G4ForceCondition* condition
106  ) override;
107 
108  private:
110  // controle flag for output message
111  // 0: Silent
112  // 1: Warning message
113  // 2: More
114 
115  private:
116  // HighestValue.
118 
119  // ParticleChange for decay process
121 
122 };
123 
124 inline
126  const G4Track& track,
127  G4double /*previousStepSize*/,
129  )
130 {
131  // pre-assigned UnknownDecay time
133 
134  if (pTime < 0.) pTime = DBL_MIN;
135 
136  // condition is set to "Not Forced"
137  *condition = NotForced;
138 
139  // reminder proper time
140  G4double remainder = pTime - track.GetProperTime();
141  if (remainder <= 0.0) remainder = DBL_MIN;
142 
143  // use pre-assigned Decay time to determine PIL
144  //return GetMeanFreePath(track, previousStepSize, condition);
145  return remainder*CLHEP::c_light;
146 
147 }
148 
149 inline
151  const G4Track& aTrack,
152  const G4Step& aStep
153  )
154 {
155  return DecayIt(aTrack, aStep);
156 }
157 
158 #endif
159 
160 
161 
162 
163 
164 
165 
166 
167 
168 
virtual void BuildPhysicsTable(const G4ParticleDefinition &) override
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override
virtual G4VParticleChange * DecayIt(const G4Track &aTrack, const G4Step &aStep)
G4UnknownDecay(const G4String &processName="UnknownDecay")
G4UnknownDecay & operator=(const G4UnknownDecay &right)
G4double condition(const G4ErrorSymMatrix &m)
G4double GetProperTime() const
virtual G4bool IsApplicable(const G4ParticleDefinition &) override
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4double GetPreAssignedDecayProperTime() const
Definition: G4Step.hh:76
virtual G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
int G4int
Definition: G4Types.hh:78
static constexpr double c_light
#define DBL_MIN
Definition: templates.hh:75
G4ParticleChangeForDecay fParticleChangeForDecay
G4ForceCondition
virtual ~G4UnknownDecay()
virtual void ProcessDescription(std::ostream &outFile) const override
const G4double HighestValue
const G4DynamicParticle * GetDynamicParticle() const