Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ContinuousGainOfEnergy.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// Class: G4ContinuousGainOfEnergy
29// Author: L. Desorgher
30// Organisation: SpaceIT GmbH
31// Contract: ESA contract 21435/08/NL/AT
32// Customer: ESA/ESTEC
33/////////////////////////////////////////////////////////////////////////////////
34//
35// CHANGE HISTORY
36// --------------
37// ChangeHistory:
38// -10 May 2007 creation by L. Desorgher
39// -February-March 2009 Update for protons by L.Desorgher
40// -July August 2009 Update for ion by L.Desorgher
41//
42//-------------------------------------------------------------
43// Documentation:
44// Continuous process acting on adjoint particles to compute the continuous gain of energy of charged particles when they are tracked back!
45//
46//
47#ifndef G4ContinuousGainOfEnergy_h
48#define G4ContinuousGainOfEnergy_h 1
49
51#include "globals.hh"
52#include "G4Material.hh"
54#include "G4Track.hh"
55#include "G4UnitsTable.hh"
56#include "G4ParticleChange.hh"
59
60
61class G4Step;
63class G4VEmModel;
65
66
67
69{
70public:
71
72 G4ContinuousGainOfEnergy(const G4String& name = "EnergyGain",
74
76
77
78protected:
79
80
81 //------------------------------------------------------------------------
82 // Methods with standard implementation; may be overwritten if needed
83 //------------------------------------------------------------------------
84protected:
85
86
87 virtual G4double GetContinuousStepLimit(const G4Track& track,
88 G4double previousStepSize,
89 G4double currentMinimumStep,
90 G4double& currentSafety);
91
92
93 //------------------------------------------------------------------------
94 // Generic methods common to all processes
95 //------------------------------------------------------------------------
96public:
97
98
99
101
103
104
106
107
108 void SetLossFluctuations(G4bool val);
109 inline void SetIsIntegral(G4bool val){is_integral= val;}
110
111 inline void SetDirectEnergyLossProcess(G4VEnergyLossProcess* aProcess){theDirectEnergyLossProcess=aProcess;};
112
114
115protected:
116
117
118
119
120private:
121
122 void DefineMaterial(const G4MaterialCutsCouple* couple);
123 void SetDynamicMassCharge(const G4Track& track, G4double energy);
124
125
126 // hide assignment operator
127
129 G4ContinuousGainOfEnergy & operator=(const G4ContinuousGainOfEnergy &right);
130
131
132private:
133
134 const G4Material* currentMaterial;
135 const G4MaterialCutsCouple* currentCouple;
136 size_t currentMaterialIndex;
137 size_t currentCoupleIndex;
138 G4double currentTcut;
139 G4double currentCutInRange;
140 G4double preStepKinEnergy;
141
142
143
144 G4double linLossLimit;
145 G4bool lossFluctuationFlag;
146 G4bool lossFluctuationArePossible;
147
148 G4VEnergyLossProcess* theDirectEnergyLossProcess;
149 G4ParticleDefinition* theDirectPartDef;
150
151
152 G4bool is_integral;
153
154 //adding for Ions
155 //----------------
156 G4bool IsIon;
157 G4double massRatio;
158 G4double chargeSqRatio;
159 G4VEmModel* currentModel;
160 G4double preStepChargeSqRatio;
161 G4double preStepScaledKinEnergy;
162 G4double preStepRange;
163
164
165
166
167
168};
169
170///////////////////////////////////////////////////////
171//
172inline void G4ContinuousGainOfEnergy::DefineMaterial(
173 const G4MaterialCutsCouple* couple)
174{
175 if(couple != currentCouple) {
176 currentCouple = couple;
177 currentMaterial = couple->GetMaterial();
178 currentCoupleIndex = couple->GetIndex();
179 currentMaterialIndex = currentMaterial->GetIndex();
180
181 size_t idx=1;
182 const std::vector<G4double>* aVec = G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(idx);
183 currentTcut=(*aVec)[currentCoupleIndex];
184 currentCutInRange = couple->GetProductionCuts()->GetProductionCut(theDirectPartDef->GetParticleName());
185 //G4cout<<"Define Material"<<G4endl;
186 //if(!meanFreePath) ResetNumberOfInteractionLengthLeft();
187 }
188}
189
190#endif
G4ProcessType
@ fElectromagnetic
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
void SetDirectEnergyLossProcess(G4VEnergyLossProcess *aProcess)
virtual G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
void BuildPhysicsTable(const G4ParticleDefinition &)
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
void SetDirectParticle(G4ParticleDefinition *p)
void PreparePhysicsTable(const G4ParticleDefinition &)
const G4Material * GetMaterial() const
G4ProductionCuts * GetProductionCuts() const
size_t GetIndex() const
Definition: G4Material.hh:258
const G4String & GetParticleName() const
const std::vector< G4double > * GetEnergyCutsVector(std::size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4double GetProductionCut(G4int index) const
Definition: G4Step.hh:62