Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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//
28//
29// ------------------------------------------------------------
30// GEANT 4 class header file
31//
32//
33#ifndef G4UnknownDecay_h
34#define G4UnknownDecay_h 1
35
37
38#include "G4ios.hh"
39#include "globals.hh"
40#include "G4VDiscreteProcess.hh"
42
44{
45 // Class Description
46 // This class is a decay process for "unknown" particle.
47
48 public:
49 // Constructors
50 G4UnknownDecay(const G4String& processName ="UnknownDecay");
51
52 // Destructor
53 virtual ~G4UnknownDecay();
54
55 private:
56 // copy constructor
57 G4UnknownDecay(const G4UnknownDecay &right);
58
59 // Assignment Operation (generated)
60 G4UnknownDecay & operator=(const G4UnknownDecay &right);
61
62 public: //With Description
63
65 const G4Track& aTrack,
66 const G4Step& aStep
67 ) override;
68
69 virtual void BuildPhysicsTable(const G4ParticleDefinition&) override;
70 // In G4UnknownDecay, thePhysicsTable stores values of
71 // beta * std::sqrt( 1 - beta*beta)
72 // as a function of normalized kinetic enregy (=Ekin/mass),
73 // becasuse this table is universal for all particle types,
74
75
76 virtual G4bool IsApplicable(const G4ParticleDefinition&) override;
77 // returns "true" if the decay process can be applied to
78 // the particle type.
79
80 virtual void ProcessDescription(std::ostream& outFile) const override;
81 //
82
83 protected: // With Description
85 const G4Track& aTrack,
86 const G4Step& aStep
87 ) ;
88 // The DecayIt() method returns by pointer a particle-change object,
89 // which has information of daughter particles.
90
91 public:
92
94 const G4Track& track,
95 G4double previousStepSize,
97 ) override;
98
99
100 protected: // With Description
101 // GetMeanFreePath returns ctau*beta*gamma for decay in flight
102 virtual G4double GetMeanFreePath(const G4Track& aTrack,
103 G4double previousStepSize,
105 ) override;
106
107 private:
108 G4int verboseLevel;
109 // controle flag for output message
110 // 0: Silent
111 // 1: Warning message
112 // 2: More
113
114 private:
115 // HighestValue.
116 const G4double HighestValue;
117
118 // ParticleChange for decay process
119 G4ParticleChangeForDecay fParticleChangeForDecay;
120
121};
122
123inline
125 const G4Track& track,
126 G4double /*previousStepSize*/,
128 )
129{
130 // pre-assigned UnknownDecay time
132
133 if (pTime < 0.) pTime = DBL_MIN;
134
135 // condition is set to "Not Forced"
137
138 // reminder proper time
139 G4double remainder = pTime - track.GetProperTime();
140 if (remainder <= 0.0) remainder = DBL_MIN;
141
142 // use pre-assigned Decay time to determine PIL
143 //return GetMeanFreePath(track, previousStepSize, condition);
144 return remainder*CLHEP::c_light;
145
146}
147
148inline
150 const G4Track& aTrack,
151 const G4Step& aStep
152 )
153{
154 return DecayIt(aTrack, aStep);
155}
156
157#endif
158
159
160
161
162
163
164
165
166
167
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ NotForced
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4double GetPreAssignedDecayProperTime() const
Definition: G4Step.hh:62
G4double GetProperTime() const
const G4DynamicParticle * GetDynamicParticle() const
virtual void BuildPhysicsTable(const G4ParticleDefinition &) override
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override
virtual G4bool IsApplicable(const G4ParticleDefinition &) override
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
virtual G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
virtual ~G4UnknownDecay()
virtual G4VParticleChange * DecayIt(const G4Track &aTrack, const G4Step &aStep)
virtual void ProcessDescription(std::ostream &outFile) const override
#define DBL_MIN
Definition: templates.hh:54