source: trunk/source/processes/electromagnetic/adjoint/include/G4ContinuousGainOfEnergy.hh@ 1036

Last change on this file since 1036 was 966, checked in by garnier, 17 years ago

fichier ajoutes

File size: 5.9 KB
Line 
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// Module: G4ContinuousGainOfEnergy.hh
28// Author: L. Desorgher
29// Date: 10 May 2007
30// Organisation: SpaceIT GmbH
31// Customer: ESA/ESTEC
32/////////////////////////////////////////////////////////////////////////////////
33//
34// CHANGE HISTORY
35// --------------
36// ChangeHistory:
37// 10 May 2007 creation by L. Desorgher
38//
39//-------------------------------------------------------------
40// Documentation:
41// Continuous process acting on adjoint particles to compute the continuous gain of energy of charged particels whern they are tracked back!
42//
43#ifndef G4ContinuousGainOfEnergy_h
44#define G4ContinuousGainOfEnergy_h 1
45
46#include "G4VContinuousProcess.hh"
47#include "globals.hh"
48#include "G4Material.hh"
49#include "G4MaterialCutsCouple.hh"
50#include "G4Track.hh"
51#include "G4UnitsTable.hh"
52#include "G4ParticleChange.hh"
53#include "G4VEnergyLossProcess.hh"
54
55
56class G4Step;
57class G4ParticleDefinition;
58class G4VEmModel;
59class G4VEmFluctuationModel;
60
61
62
63class G4ContinuousGainOfEnergy : public G4VContinuousProcess
64{
65public:
66
67 G4ContinuousGainOfEnergy(const G4String& name = "EnergyGain",
68 G4ProcessType type = fElectromagnetic);
69
70 virtual ~G4ContinuousGainOfEnergy();
71
72
73protected:
74
75
76 //------------------------------------------------------------------------
77 // Methods with standard implementation; may be overwritten if needed
78 //------------------------------------------------------------------------
79protected:
80
81
82 virtual G4double GetContinuousStepLimit(const G4Track& track,
83 G4double previousStepSize,
84 G4double currentMinimumStep,
85 G4double& currentSafety);
86
87
88 //------------------------------------------------------------------------
89 // Generic methods common to all processes
90 //------------------------------------------------------------------------
91public:
92
93
94
95 void PreparePhysicsTable(const G4ParticleDefinition&);
96
97 void BuildPhysicsTable(const G4ParticleDefinition&);
98
99
100 G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&);
101
102
103 void SetLossFluctuations(G4bool val);
104 inline void SetIsIntegral(G4bool val){is_integral= val;}
105
106 inline void SetDirectEnergyLossProcess(G4VEnergyLossProcess* aProcess){theDirectEnergyLossProcess=aProcess;};
107
108 inline void SetDirectParticle(G4ParticleDefinition* p){theDirectPartDef=p;};
109
110protected:
111
112
113
114
115private:
116
117 void DefineMaterial(const G4MaterialCutsCouple* couple);
118
119 // hide assignment operator
120
121 G4ContinuousGainOfEnergy(G4ContinuousGainOfEnergy &);
122 G4ContinuousGainOfEnergy & operator=(const G4ContinuousGainOfEnergy &right);
123
124
125private:
126
127 const G4Material* currentMaterial;
128 const G4MaterialCutsCouple* currentCouple;
129 size_t currentMaterialIndex;
130 G4double currentTcut;
131 G4double preStepKinEnergy;
132
133
134 G4double linLossLimit;
135 G4bool lossFluctuationFlag;
136 G4bool lossFluctuationArePossible;
137
138 G4VEnergyLossProcess* theDirectEnergyLossProcess;
139 G4ParticleDefinition* theDirectPartDef;
140
141
142 G4bool is_integral;
143
144
145};
146
147///////////////////////////////////////////////////////
148//
149inline void G4ContinuousGainOfEnergy::DefineMaterial(
150 const G4MaterialCutsCouple* couple)
151{
152 if(couple != currentCouple) {
153 currentCouple = couple;
154 currentMaterial = couple->GetMaterial();
155 currentMaterialIndex = couple->GetIndex();
156 currentTcut = couple->GetProductionCuts()->GetProductionCut(theDirectPartDef->GetParticleName());
157 //G4cout<<"Define Material"<<std::endl;
158 //if(!meanFreePath) ResetNumberOfInteractionLengthLeft();
159 }
160}
161///////////////////////////////////////////////////////
162//
163inline G4double G4ContinuousGainOfEnergy::GetContinuousStepLimit(const G4Track& track,
164 G4double , G4double , G4double& )
165{
166 G4double x = DBL_MAX;
167 x=.1*mm;
168
169 //G4cout<<x<<std::endl;
170 DefineMaterial(track.GetMaterialCutsCouple());
171 preStepKinEnergy = track.GetKineticEnergy();
172 G4double maxE=1.2*preStepKinEnergy;
173 G4double r = theDirectEnergyLossProcess->GetRange(preStepKinEnergy, currentCouple);
174 G4double r1 = theDirectEnergyLossProcess->GetRange(maxE, currentCouple);
175 x=std::max(r1-r,.1);
176
177 return x;
178
179
180}
181#endif
Note: See TracBrowser for help on using the repository browser.