Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ParticleChangeForLoss.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 // $Id: G4ParticleChangeForLoss.cc 68795 2013-04-05 13:24:46Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------
31 // GEANT 4 class implementation file
32 //
33 // ------------------------------------------------------------
34 // Implemented for the new scheme 23 Mar. 1998 H.Kurahige
35 // --------------------------------------------------------------
36 //
37 // Modified:
38 // 16.01.04 V.Ivanchenko update for model variant of energy loss
39 // 15.04.05 V.Ivanchenko inline update methods
40 // 28.08.06 V.Ivanchenko Add access to current track and polarizaion
41 //
42 // ------------------------------------------------------------
43 //
45 #include "G4SystemOfUnits.hh"
46 #include "G4Track.hh"
47 #include "G4Step.hh"
48 #include "G4DynamicParticle.hh"
49 #include "G4ExceptionSeverity.hh"
50 //#include "G4VelocityTable.hh"
51 
53  : G4VParticleChange(), currentTrack(0), proposedKinEnergy(0.),
54  lowEnergyLimit(1.0*eV), currentCharge(0.)
55 {
57  debugFlag = false;
58 #ifdef G4VERBOSE
59  if (verboseLevel>2) {
60  G4cout << "G4ParticleChangeForLoss::G4ParticleChangeForLoss() " << G4endl;
61  }
62 #endif
63 }
64 
66 {
67 #ifdef G4VERBOSE
68  if (verboseLevel>2) {
69  G4cout << "G4ParticleChangeForLoss::~G4ParticleChangeForLoss() " << G4endl;
70  }
71 #endif
72 }
73 
76  : G4VParticleChange(right)
77 {
78  if (verboseLevel>1) {
79  G4cout << "G4ParticleChangeForLoss:: copy constructor is called " << G4endl;
80  }
81  currentTrack = right.currentTrack;
82  proposedKinEnergy = right.proposedKinEnergy;
83  lowEnergyLimit = right.lowEnergyLimit;
84  currentCharge = right.currentCharge;
85  proposedMomentumDirection = right.proposedMomentumDirection;
86 }
87 
88 // assignment operator
91 {
92 #ifdef G4VERBOSE
93  if (verboseLevel>1) {
94  G4cout << "G4ParticleChangeForLoss:: assignment operator is called " << G4endl;
95  }
96 #endif
97 
98  if (this != &right) {
99  if (theNumberOfSecondaries>0) {
100 #ifdef G4VERBOSE
101  if (verboseLevel>0) {
102  G4cout << "G4ParticleChangeForLoss: assignment operator Warning ";
103  G4cout << "theListOfSecondaries is not empty ";
104  }
105 #endif
106  for (G4int index= 0; index<theNumberOfSecondaries; index++){
107  if ( (*theListOfSecondaries)[index] ) delete (*theListOfSecondaries)[index] ;
108  }
109  }
110  delete theListOfSecondaries;
113  for (G4int index = 0; index<theNumberOfSecondaries; index++){
114  G4Track* newTrack = new G4Track(*((*right.theListOfSecondaries)[index] ));
115  theListOfSecondaries->SetElement(index, newTrack); }
116 
123 
124  currentTrack = right.currentTrack;
125  proposedKinEnergy = right.proposedKinEnergy;
126  currentCharge = right.currentCharge;
127  proposedMomentumDirection = right.proposedMomentumDirection;
128  }
129  return *this;
130 }
131 
132 //----------------------------------------------------------------
133 // methods for printing messages
134 //
135 
137 {
138 // use base-class DumpInfo
140 
141  G4int oldprc = G4cout.precision(3);
142  G4cout << " Charge (eplus) : "
143  << std::setw(20) << currentCharge/eplus
144  << G4endl;
145  G4cout << " Kinetic Energy (MeV): "
146  << std::setw(20) << proposedKinEnergy/MeV
147  << G4endl;
148  G4cout << " Momentum Direct - x : "
149  << std::setw(20) << proposedMomentumDirection.x()
150  << G4endl;
151  G4cout << " Momentum Direct - y : "
152  << std::setw(20) << proposedMomentumDirection.y()
153  << G4endl;
154  G4cout << " Momentum Direct - z : "
155  << std::setw(20) << proposedMomentumDirection.z()
156  << G4endl;
157  G4cout.precision(oldprc);
158 }
159 
161 {
162  G4bool itsOK = true;
163  G4bool exitWithError = false;
164 
165  G4double accuracy;
166 
167  // Energy should not be lager than initial value
168  accuracy = ( proposedKinEnergy - aTrack.GetKineticEnergy())/MeV;
169  if (accuracy > accuracyForWarning) {
170  itsOK = false;
171  exitWithError = (accuracy > accuracyForException);
172 #ifdef G4VERBOSE
173  G4cout << "G4ParticleChangeForLoss::CheckIt: ";
174  G4cout << "KinEnergy become larger than the initial value!"
175  << " Difference: " << accuracy << "[MeV] " <<G4endl;
176  G4cout << aTrack.GetDefinition()->GetParticleName()
177  << " E=" << aTrack.GetKineticEnergy()/MeV
178  << " pos=" << aTrack.GetPosition().x()/m
179  << ", " << aTrack.GetPosition().y()/m
180  << ", " << aTrack.GetPosition().z()/m
181  <<G4endl;
182 #endif
183  }
184 
185  // dump out information of this particle change
186 #ifdef G4VERBOSE
187  if (!itsOK) DumpInfo();
188 #endif
189 
190  // Exit with error
191  if (exitWithError) {
192  G4Exception("G4ParticleChangeForLoss::CheckIt",
193  "TRACK004", EventMustBeAborted,
194  "energy was illegal");
195  }
196 
197  //correction
198  if (!itsOK) {
199  proposedKinEnergy = aTrack.GetKineticEnergy();
200  }
201 
202  itsOK = (itsOK) && G4VParticleChange::CheckIt(aTrack);
203  return itsOK;
204 }
205 
206 //----------------------------------------------------------------
207 // methods for updating G4Step
208 //
209 
211 {
212  G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
213  G4StepPoint* pPreStepPoint = pStep->GetPreStepPoint();
214  G4Track* pTrack = pStep->GetTrack();
215 
216  // accumulate change of the kinetic energy
217  G4double preKinEnergy = pPreStepPoint->GetKineticEnergy();
218  G4double kinEnergy = pPostStepPoint->GetKineticEnergy() +
219  (proposedKinEnergy - preKinEnergy);
220 
221  // update kinetic energy and charge
222  if (kinEnergy < lowEnergyLimit) {
223  theLocalEnergyDeposit += kinEnergy;
224  kinEnergy = 0.0;
225  pPostStepPoint->SetVelocity(0.0);
226  } else {
227  pPostStepPoint->SetCharge( currentCharge );
228  // calculate velocity
229  pTrack->SetKineticEnergy(kinEnergy);
230  pPostStepPoint->SetVelocity(pTrack->CalculateVelocity());
231  pTrack->SetKineticEnergy(preKinEnergy);
232  }
233  pPostStepPoint->SetKineticEnergy( kinEnergy );
234 
236  pPostStepPoint->SetWeight( theParentWeight );
237  }
238 
241  return pStep;
242 }
243 
245 {
246  G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
247  G4Track* pTrack = pStep->GetTrack();
248 
249  pPostStepPoint->SetCharge( currentCharge );
250  pPostStepPoint->SetMomentumDirection( proposedMomentumDirection );
251  pPostStepPoint->SetKineticEnergy( proposedKinEnergy );
252  pTrack->SetKineticEnergy( proposedKinEnergy );
253  if(proposedKinEnergy > 0.0) {
254  pPostStepPoint->SetVelocity(pTrack->CalculateVelocity());
255  } else {
256  pPostStepPoint->SetVelocity(0.0);
257  }
258  pPostStepPoint->SetPolarization( proposedPolarization );
259 
261  pPostStepPoint->SetWeight( theParentWeight );
262  }
263 
266  return pStep;
267 }
268 
G4ParticleDefinition * GetDefinition() const
virtual void DumpInfo() const
virtual G4bool CheckIt(const G4Track &)
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:76
void AddNonIonizingEnergyDeposit(G4double value)
G4Step * UpdateStepForAlongStep(G4Step *Step)
double x() const
const G4ThreeVector & GetPosition() const
G4TrackFastVector * theListOfSecondaries
void SetWeight(G4double aValue)
virtual void DumpInfo() const
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
double z() const
static const G4double accuracyForException
G4Step * UpdateStepForPostStep(G4Step *Step)
void SetMomentumDirection(const G4ThreeVector &aValue)
G4StepPoint * GetPreStepPoint() const
G4ParticleChangeForLoss & operator=(const G4ParticleChangeForLoss &right)
virtual G4bool CheckIt(const G4Track &)
G4double GetKineticEnergy() const
void SetPolarization(const G4ThreeVector &aValue)
G4GLOB_DLL std::ostream G4cout
static constexpr double m
Definition: G4SIunits.hh:129
bool G4bool
Definition: G4Types.hh:79
G4SteppingControl theSteppingControlFlag
static constexpr double eplus
Definition: G4SIunits.hh:199
static constexpr double eV
Definition: G4SIunits.hh:215
Definition: G4Step.hh:76
G4double CalculateVelocity() const
Definition: G4Track.cc:222
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void SetVelocity(G4double v)
G4FastVector< G4Track, G4TrackFastVectorSize > G4TrackFastVector
void SetCharge(G4double value)
G4StepPoint * GetPostStepPoint() const
void AddTotalEnergyDeposit(G4double value)
double y() const
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4TrackStatus theStatusChange
void SetKineticEnergy(const G4double aValue)
G4double GetKineticEnergy() const
double G4double
Definition: G4Types.hh:76
static const G4double accuracyForWarning
void SetKineticEnergy(const G4double aValue)
G4Track * GetTrack() const
G4double theNonIonizingEnergyDeposit