Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4UnknownDecay.hh
Go to the documentation of this file.
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$
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
58 G4UnknownDecay(const G4UnknownDecay &right);
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 );
69
70 virtual void BuildPhysicsTable(const G4ParticleDefinition&);
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
78 // returns "true" if the decay process can be applied to
79 // the particle type.
80
81 protected: // With Description
83 const G4Track& aTrack,
84 const G4Step& aStep
85 );
86 // The DecayIt() method returns by pointer a particle-change object,
87 // which has information of daughter particles.
88
89 public:
90
92 const G4Track& track,
93 G4double previousStepSize,
95 );
96
97
98 protected: // With Description
99 // GetMeanFreePath returns ctau*beta*gamma for decay in flight
100 virtual G4double GetMeanFreePath(const G4Track& aTrack,
101 G4double previousStepSize,
103 );
104
105 public:
106 void SetVerboseLevel(G4int value);
107 G4int GetVerboseLevel() const;
108
109 private:
110 G4int verboseLevel;
111 // controle flag for output message
112 // 0: Silent
113 // 1: Warning message
114 // 2: More
115
116 private:
117 // HighestValue.
118 const G4double HighestValue;
119
120 // ParticleChange for decay process
121 G4ParticleChangeForDecay fParticleChangeForDecay;
122
123};
124
125inline
127 const G4Track& track,
128 G4double /*previousStepSize*/,
130 )
131{
132 // pre-assigned UnknownDecay time
134
135 if (pTime < 0.) pTime = DBL_MIN;
136
137 // condition is set to "Not Forced"
139
140 // reminder proper time
141 G4double remainder = pTime - track.GetProperTime();
142 if (remainder <= 0.0) remainder = DBL_MIN;
143
144 // use pre-assigned Decay time to determine PIL
145 //return GetMeanFreePath(track, previousStepSize, condition);
146 return remainder*CLHEP::c_light;
147
148}
149
150inline
151 void G4UnknownDecay::SetVerboseLevel(G4int value){ verboseLevel = value; }
152
153inline
154 G4int G4UnknownDecay::GetVerboseLevel() const { return verboseLevel; }
155
156inline
158 const G4Track& aTrack,
159 const G4Step& aStep
160 )
161{
162 return DecayIt(aTrack, aStep);
163}
164
165#endif
166
167
168
169
170
171
172
173
174
175
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ NotForced
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4double GetPreAssignedDecayProperTime() const
Definition: G4Step.hh:78
G4double GetProperTime() const
const G4DynamicParticle * GetDynamicParticle() const
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
void SetVerboseLevel(G4int value)
virtual G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
G4int GetVerboseLevel() const
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
virtual G4bool IsApplicable(const G4ParticleDefinition &)
virtual ~G4UnknownDecay()
virtual G4VParticleChange * DecayIt(const G4Track &aTrack, const G4Step &aStep)
#define DBL_MIN
Definition: templates.hh:75