Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4UnknownDecay.cc
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 implementation file
31//
32// ------------------------------------------------------------
33//
34
35#include "G4UnknownDecay.hh"
36
38#include "G4SystemOfUnits.hh"
39#include "G4DynamicParticle.hh"
40#include "G4DecayProducts.hh"
41#include "G4PhysicsLogVector.hh"
43#include "G4DecayProcessType.hh"
44
45// constructor
47 :G4VDiscreteProcess(processName, fDecay),
48 verboseLevel(1),
49 HighestValue(20.0)
50{
51 // set Process Sub Type
52 SetProcessSubType(static_cast<int>(DECAY_Unknown));
53
54#ifdef G4VERBOSE
55 if (GetVerboseLevel()>1) {
56 G4cout << "G4UnknownDecay constructor " << " Name:" << processName << G4endl;
57 }
58#endif
59 pParticleChange = &fParticleChangeForDecay;
60}
61
63{
64}
65
67{
68 if(aParticleType.GetParticleName()=="unknown") return true;
69 return false;
70}
71
73{
74 return 0.0;
75}
76
78{
79 return;
80}
81
83{
84 // The DecayIt() method returns by pointer a particle-change object.
85 // Units are expressed in GEANT4 internal units.
86
87 // Initialize ParticleChange
88 // all members of G4VParticleChange are set to equal to
89 // corresponding member in G4Track
90 fParticleChangeForDecay.Initialize(aTrack);
91
92 // get particle
93 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
94
95 //check if thePreAssignedDecayProducts exists
96 const G4DecayProducts* o_products = (aParticle->GetPreAssignedDecayProducts());
97 G4bool isPreAssigned = (o_products != nullptr);
98 G4DecayProducts* products = nullptr;
99
100 if (!isPreAssigned ){
101 fParticleChangeForDecay.SetNumberOfSecondaries(0);
102 // Kill the parent particle
103 fParticleChangeForDecay.ProposeTrackStatus( fStopAndKill ) ;
104 fParticleChangeForDecay.ProposeLocalEnergyDeposit(0.0);
105
107 return &fParticleChangeForDecay ;
108 }
109
110 // copy decay products
111 products = new G4DecayProducts(*o_products);
112
113 // get parent particle information ...................................
114 G4double ParentEnergy = aParticle->GetTotalEnergy();
115 G4double ParentMass = aParticle->GetMass();
116 if (ParentEnergy < ParentMass) {
117 ParentEnergy = ParentMass;
118#ifdef G4VERBOSE
119 if (GetVerboseLevel()>1) {
120 G4cout << "G4UnknownDecay::DoIt : Total Energy is less than its mass" << G4endl;
121 G4cout << " Particle: " << aParticle->GetDefinition()->GetParticleName();
122 G4cout << " Energy:" << ParentEnergy/MeV << "[MeV]";
123 G4cout << " Mass:" << ParentMass/MeV << "[MeV]";
124 G4cout << G4endl;
125 }
126#endif
127 }
128 G4ThreeVector ParentDirection(aParticle->GetMomentumDirection());
129
130 G4double energyDeposit = 0.0;
131 G4double finalGlobalTime = aTrack.GetGlobalTime();
132 //boost all decay products to laboratory frame
133 //if the particle has traveled
134 if(aParticle->GetPreAssignedDecayProperTime()>=0.) {
135 products->Boost( ParentEnergy, ParentDirection);
136 }
137
138 //add products in fParticleChangeForDecay
139 G4int numberOfSecondaries = products->entries();
140 fParticleChangeForDecay.SetNumberOfSecondaries(numberOfSecondaries);
141#ifdef G4VERBOSE
142 if (GetVerboseLevel()>1) {
143 G4cout << "G4UnknownDecay::DoIt : Decay vertex :";
144 G4cout << " Time: " << finalGlobalTime/ns << "[ns]";
145 G4cout << " X:" << (aTrack.GetPosition()).x() /cm << "[cm]";
146 G4cout << " Y:" << (aTrack.GetPosition()).y() /cm << "[cm]";
147 G4cout << " Z:" << (aTrack.GetPosition()).z() /cm << "[cm]";
148 G4cout << G4endl;
149 G4cout << "G4UnknownDecay::DoIt : decay products in Lab. Frame" << G4endl;
150 products->DumpInfo();
151 }
152#endif
153 G4int index;
154 G4ThreeVector currentPosition;
155 const G4TouchableHandle thand = aTrack.GetTouchableHandle();
156 for (index=0; index < numberOfSecondaries; index++){
157 // get current position of the track
158 currentPosition = aTrack.GetPosition();
159 // create a new track object
160 G4Track* secondary = new G4Track( products->PopProducts(),
161 finalGlobalTime ,
162 currentPosition );
163 // switch on good for tracking flag
164 secondary->SetGoodForTrackingFlag();
165 secondary->SetTouchableHandle(thand);
166 // add the secondary track in the List
167 fParticleChangeForDecay.AddSecondary(secondary);
168 }
169 delete products;
170
171 // Kill the parent particle
172 fParticleChangeForDecay.ProposeTrackStatus( fStopAndKill ) ;
173 fParticleChangeForDecay.ProposeLocalEnergyDeposit(energyDeposit);
174 fParticleChangeForDecay.ProposeGlobalTime( finalGlobalTime );
175 // reset NumberOfInteractionLengthLeft
177
178 return &fParticleChangeForDecay ;
179}
180
181void G4UnknownDecay::ProcessDescription(std::ostream& outFile) const
182{
183 outFile << GetProcessName()
184 << ": Decay of 'unknown' particles. \n"
185 << "kinematics of daughters are dertermined "
186 << "by PreAssignedDecayProducts. \n";
187}
188
189
190
191
192
@ DECAY_Unknown
G4ForceCondition
@ fDecay
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void DumpInfo() const
G4int entries() const
G4DynamicParticle * PopProducts()
void Boost(G4double totalEnergy, const G4ThreeVector &momentumDirection)
G4double GetMass() const
const G4ThreeVector & GetMomentumDirection() const
const G4DecayProducts * GetPreAssignedDecayProducts() const
G4ParticleDefinition * GetDefinition() const
G4double GetTotalEnergy() const
G4double GetPreAssignedDecayProperTime() const
void Initialize(const G4Track &) final
const G4String & GetParticleName() const
Definition: G4Step.hh:62
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
void SetGoodForTrackingFlag(G4bool value=true)
virtual void BuildPhysicsTable(const G4ParticleDefinition &) override
G4UnknownDecay(const G4String &processName="UnknownDecay")
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override
virtual G4bool IsApplicable(const G4ParticleDefinition &) override
virtual ~G4UnknownDecay()
virtual G4VParticleChange * DecayIt(const G4Track &aTrack, const G4Step &aStep)
virtual void ProcessDescription(std::ostream &outFile) const override
void ProposeTrackStatus(G4TrackStatus status)
void AddSecondary(G4Track *aSecondary)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4int GetVerboseLevel() const
Definition: G4VProcess.hh:422
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:428
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:325
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
#define ns(x)
Definition: xmltok.c:1649