Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ParticleChangeForLoss Class Reference

#include <G4ParticleChangeForLoss.hh>

Inheritance diagram for G4ParticleChangeForLoss:
Collaboration diagram for G4ParticleChangeForLoss:

Public Member Functions

 G4ParticleChangeForLoss ()
 
virtual ~G4ParticleChangeForLoss ()
 
G4StepUpdateStepForAlongStep (G4Step *Step)
 
G4StepUpdateStepForPostStep (G4Step *Step)
 
void InitializeForAlongStep (const G4Track &)
 
void InitializeForPostStep (const G4Track &)
 
G4double GetProposedCharge () const
 
void SetProposedCharge (G4double theCharge)
 
G4double GetCharge () const
 
void ProposeCharge (G4double finalCharge)
 
G4double GetProposedKineticEnergy () const
 
void SetProposedKineticEnergy (G4double proposedKinEnergy)
 
const G4ThreeVectorGetProposedMomentumDirection () const
 
void SetProposedMomentumDirection (const G4ThreeVector &dir)
 
const G4ThreeVectorGetMomentumDirection () const
 
void ProposeMomentumDirection (G4double Px, G4double Py, G4double Pz)
 
void ProposeMomentumDirection (const G4ThreeVector &Pfinal)
 
const G4ThreeVectorGetProposedPolarization () const
 
void ProposePolarization (const G4ThreeVector &dir)
 
void ProposePolarization (G4double Px, G4double Py, G4double Pz)
 
const G4TrackGetCurrentTrack () const
 
void SetLowEnergyLimit (G4double elimit)
 
virtual void DumpInfo () const
 
virtual G4bool CheckIt (const G4Track &)
 
- Public Member Functions inherited from G4VParticleChange
 G4VParticleChange ()
 
virtual ~G4VParticleChange ()
 
G4bool operator== (const G4VParticleChange &right) const
 
G4bool operator!= (const G4VParticleChange &right) const
 
virtual G4StepUpdateStepForAtRest (G4Step *Step)
 
virtual void Initialize (const G4Track &)
 
G4double GetTrueStepLength () const
 
void ProposeTrueStepLength (G4double truePathLength)
 
G4double GetLocalEnergyDeposit () const
 
void ProposeLocalEnergyDeposit (G4double anEnergyPart)
 
G4double GetNonIonizingEnergyDeposit () const
 
void ProposeNonIonizingEnergyDeposit (G4double anEnergyPart)
 
G4TrackStatus GetTrackStatus () const
 
void ProposeTrackStatus (G4TrackStatus status)
 
G4SteppingControl GetSteppingControl () const
 
void ProposeSteppingControl (G4SteppingControl StepControlFlag)
 
G4bool GetFirstStepInVolume () const
 
G4bool GetLastStepInVolume () const
 
void ProposeFirstStepInVolume (G4bool flag)
 
void ProposeLastStepInVolume (G4bool flag)
 
void Clear ()
 
void SetNumberOfSecondaries (G4int totSecondaries)
 
G4int GetNumberOfSecondaries () const
 
G4TrackGetSecondary (G4int anIndex) const
 
void AddSecondary (G4Track *aSecondary)
 
G4double GetWeight () const
 
G4double GetParentWeight () const
 
void ProposeWeight (G4double finalWeight)
 
void ProposeParentWeight (G4double finalWeight)
 
void SetSecondaryWeightByProcess (G4bool)
 
G4bool IsSecondaryWeightSetByProcess () const
 
void SetParentWeightByProcess (G4bool)
 
G4bool IsParentWeightSetByProcess () const
 
void SetVerboseLevel (G4int vLevel)
 
G4int GetVerboseLevel () const
 
void ClearDebugFlag ()
 
void SetDebugFlag ()
 
G4bool GetDebugFlag () const
 

Protected Member Functions

 G4ParticleChangeForLoss (const G4ParticleChangeForLoss &right)
 
G4ParticleChangeForLossoperator= (const G4ParticleChangeForLoss &right)
 
- Protected Member Functions inherited from G4VParticleChange
 G4VParticleChange (const G4VParticleChange &right)
 
G4VParticleChangeoperator= (const G4VParticleChange &right)
 
G4StepUpdateStepInfo (G4Step *Step)
 
void InitializeTrueStepLength (const G4Track &)
 
void InitializeLocalEnergyDeposit (const G4Track &)
 
void InitializeSteppingControl (const G4Track &)
 
void InitializeParentWeight (const G4Track &)
 
void InitializeParentGlobalTime (const G4Track &)
 
void InitializeStatusChange (const G4Track &)
 
void InitializeSecondaries (const G4Track &)
 
void InitializeStepInVolumeFlags (const G4Track &)
 
G4bool CheckSecondary (G4Track &)
 
G4double GetAccuracyForWarning () const
 
G4double GetAccuracyForException () const
 

Additional Inherited Members

- Protected Attributes inherited from G4VParticleChange
G4TrackFastVectortheListOfSecondaries
 
G4int theNumberOfSecondaries
 
G4int theSizeOftheListOfSecondaries
 
G4TrackStatus theStatusChange
 
G4SteppingControl theSteppingControlFlag
 
G4double theLocalEnergyDeposit
 
G4double theNonIonizingEnergyDeposit
 
G4double theTrueStepLength
 
G4bool theFirstStepInVolume
 
G4bool theLastStepInVolume
 
G4double theParentWeight
 
G4bool isParentWeightProposed
 
G4bool fSetSecondaryWeightByProcess
 
G4double theParentGlobalTime
 
G4int verboseLevel
 
G4bool debugFlag
 
- Static Protected Attributes inherited from G4VParticleChange
static const G4double accuracyForWarning = 1.0e-9
 
static const G4double accuracyForException = 0.001
 

Detailed Description

Definition at line 60 of file G4ParticleChangeForLoss.hh.

Constructor & Destructor Documentation

G4ParticleChangeForLoss::G4ParticleChangeForLoss ( )

Definition at line 52 of file G4ParticleChangeForLoss.cc.

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 }
G4GLOB_DLL std::ostream G4cout
G4SteppingControl theSteppingControlFlag
static constexpr double eV
Definition: G4SIunits.hh:215
#define G4endl
Definition: G4ios.hh:61
G4ParticleChangeForLoss::~G4ParticleChangeForLoss ( )
virtual

Definition at line 65 of file G4ParticleChangeForLoss.cc.

66 {
67 #ifdef G4VERBOSE
68  if (verboseLevel>2) {
69  G4cout << "G4ParticleChangeForLoss::~G4ParticleChangeForLoss() " << G4endl;
70  }
71 #endif
72 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
G4ParticleChangeForLoss::G4ParticleChangeForLoss ( const G4ParticleChangeForLoss right)
protected

Definition at line 75 of file G4ParticleChangeForLoss.cc.

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 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Member Function Documentation

G4bool G4ParticleChangeForLoss::CheckIt ( const G4Track aTrack)
virtual

Reimplemented from G4VParticleChange.

Definition at line 160 of file G4ParticleChangeForLoss.cc.

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 }
G4ParticleDefinition * GetDefinition() const
virtual void DumpInfo() const
double x() const
const G4ThreeVector & GetPosition() const
const G4String & GetParticleName() const
double z() const
static const G4double accuracyForException
virtual G4bool CheckIt(const G4Track &)
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
static constexpr double m
Definition: G4SIunits.hh:129
bool G4bool
Definition: G4Types.hh:79
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
double y() const
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
static const G4double accuracyForWarning

Here is the call graph for this function:

void G4ParticleChangeForLoss::DumpInfo ( ) const
virtual

Reimplemented from G4VParticleChange.

Definition at line 136 of file G4ParticleChangeForLoss.cc.

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 }
double x() const
virtual void DumpInfo() const
int G4int
Definition: G4Types.hh:78
double z() const
G4GLOB_DLL std::ostream G4cout
static constexpr double eplus
Definition: G4SIunits.hh:199
double y() const
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4ParticleChangeForLoss::GetCharge ( ) const
inline

Definition at line 160 of file G4ParticleChangeForLoss.hh.

161 {
162  return currentCharge;
163 }
const G4Track * G4ParticleChangeForLoss::GetCurrentTrack ( ) const
inline

Definition at line 207 of file G4ParticleChangeForLoss.hh.

208 {
209  return currentTrack;
210 }

Here is the caller graph for this function:

const G4ThreeVector & G4ParticleChangeForLoss::GetMomentumDirection ( ) const
inline

Definition at line 182 of file G4ParticleChangeForLoss.hh.

183 {
184  return proposedMomentumDirection;
185 }
G4double G4ParticleChangeForLoss::GetProposedCharge ( ) const
inline

Definition at line 155 of file G4ParticleChangeForLoss.hh.

156 {
157  return currentCharge;
158 }
G4double G4ParticleChangeForLoss::GetProposedKineticEnergy ( ) const
inline

Definition at line 145 of file G4ParticleChangeForLoss.hh.

146 {
147  return proposedKinEnergy;
148 }

Here is the caller graph for this function:

const G4ThreeVector & G4ParticleChangeForLoss::GetProposedMomentumDirection ( ) const
inline

Definition at line 176 of file G4ParticleChangeForLoss.hh.

177 {
178  return proposedMomentumDirection;
179 }

Here is the caller graph for this function:

const G4ThreeVector & G4ParticleChangeForLoss::GetProposedPolarization ( ) const
inline

Definition at line 213 of file G4ParticleChangeForLoss.hh.

214 {
215  return proposedPolarization;
216 }
void G4ParticleChangeForLoss::InitializeForAlongStep ( const G4Track track)
inline

Definition at line 232 of file G4ParticleChangeForLoss.hh.

233 {
235  theLocalEnergyDeposit = 0.0;
237  InitializeSecondaries(track);
238  theParentWeight = track.GetWeight();
239  // isParentWeightProposed = false;
240  proposedKinEnergy = track.GetKineticEnergy();
241  currentCharge = track.GetDynamicParticle()->GetCharge();
242 }
const G4DynamicParticle * GetDynamicParticle() const
G4TrackStatus GetTrackStatus() const
G4double GetKineticEnergy() const
G4double GetCharge() const
void InitializeSecondaries(const G4Track &)
G4double GetWeight() const
G4TrackStatus theStatusChange
G4double theNonIonizingEnergyDeposit

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ParticleChangeForLoss::InitializeForPostStep ( const G4Track track)
inline

Definition at line 244 of file G4ParticleChangeForLoss.hh.

245 {
247  theLocalEnergyDeposit = 0.0;
249  InitializeSecondaries(track);
250  theParentWeight = track.GetWeight();
251  // isParentWeightProposed = false;
252  proposedKinEnergy = track.GetKineticEnergy();
253  currentCharge = track.GetDynamicParticle()->GetCharge();
254  proposedMomentumDirection = track.GetMomentumDirection();
255  proposedPolarization = track.GetPolarization();
256  currentTrack = &track;
257 }
const G4ThreeVector & GetPolarization() const
const G4DynamicParticle * GetDynamicParticle() const
G4TrackStatus GetTrackStatus() const
G4double GetKineticEnergy() const
G4double GetCharge() const
void InitializeSecondaries(const G4Track &)
const G4ThreeVector & GetMomentumDirection() const
G4double GetWeight() const
G4TrackStatus theStatusChange
G4double theNonIonizingEnergyDeposit

Here is the call graph for this function:

Here is the caller graph for this function:

G4ParticleChangeForLoss & G4ParticleChangeForLoss::operator= ( const G4ParticleChangeForLoss right)
protected

Definition at line 89 of file G4ParticleChangeForLoss.cc.

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;
112  theNumberOfSecondaries = right.theNumberOfSecondaries;
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 }
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:76
G4TrackFastVector * theListOfSecondaries
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4SteppingControl theSteppingControlFlag
G4FastVector< G4Track, G4TrackFastVectorSize > G4TrackFastVector
#define G4endl
Definition: G4ios.hh:61
G4TrackStatus theStatusChange

Here is the call graph for this function:

void G4ParticleChangeForLoss::ProposeCharge ( G4double  finalCharge)
inline

Definition at line 170 of file G4ParticleChangeForLoss.hh.

171 {
172  currentCharge = theCharge;
173 }
void G4ParticleChangeForLoss::ProposeMomentumDirection ( G4double  Px,
G4double  Py,
G4double  Pz 
)
inline

Definition at line 200 of file G4ParticleChangeForLoss.hh.

201 {
202  proposedMomentumDirection.setX(Px);
203  proposedMomentumDirection.setY(Py);
204  proposedMomentumDirection.setZ(Pz);
205 }
void setY(double)
void setZ(double)
void setX(double)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ParticleChangeForLoss::ProposeMomentumDirection ( const G4ThreeVector Pfinal)
inline

Definition at line 188 of file G4ParticleChangeForLoss.hh.

189 {
190  proposedMomentumDirection = dir;
191 }
void G4ParticleChangeForLoss::ProposePolarization ( const G4ThreeVector dir)
inline

Definition at line 219 of file G4ParticleChangeForLoss.hh.

220 {
221  proposedPolarization = dir;
222 }

Here is the caller graph for this function:

void G4ParticleChangeForLoss::ProposePolarization ( G4double  Px,
G4double  Py,
G4double  Pz 
)
inline

Definition at line 225 of file G4ParticleChangeForLoss.hh.

226 {
227  proposedPolarization.setX(Px);
228  proposedPolarization.setY(Py);
229  proposedPolarization.setZ(Pz);
230 }
void setY(double)
void setZ(double)
void setX(double)

Here is the call graph for this function:

void G4ParticleChangeForLoss::SetLowEnergyLimit ( G4double  elimit)
inline

Definition at line 273 of file G4ParticleChangeForLoss.hh.

274 {
275  lowEnergyLimit = elimit;
276 }
void G4ParticleChangeForLoss::SetProposedCharge ( G4double  theCharge)
inline

Definition at line 165 of file G4ParticleChangeForLoss.hh.

166 {
167  currentCharge = theCharge;
168 }

Here is the caller graph for this function:

void G4ParticleChangeForLoss::SetProposedKineticEnergy ( G4double  proposedKinEnergy)
inline

Definition at line 150 of file G4ParticleChangeForLoss.hh.

151 {
152  proposedKinEnergy = energy;
153 }
G4double energy(const ThreeVector &p, const G4double m)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ParticleChangeForLoss::SetProposedMomentumDirection ( const G4ThreeVector dir)
inline

Definition at line 194 of file G4ParticleChangeForLoss.hh.

195 {
196  proposedMomentumDirection = dir;
197 }

Here is the caller graph for this function:

G4Step * G4ParticleChangeForLoss::UpdateStepForAlongStep ( G4Step Step)
virtual

Reimplemented from G4VParticleChange.

Definition at line 210 of file G4ParticleChangeForLoss.cc.

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 
239  pStep->AddTotalEnergyDeposit( theLocalEnergyDeposit );
240  pStep->AddNonIonizingEnergyDeposit( theNonIonizingEnergyDeposit );
241  return pStep;
242 }
void SetWeight(G4double aValue)
G4double CalculateVelocity() const
Definition: G4Track.cc:222
void SetVelocity(G4double v)
void SetCharge(G4double value)
void SetKineticEnergy(const G4double aValue)
G4double GetKineticEnergy() const
double G4double
Definition: G4Types.hh:76
void SetKineticEnergy(const G4double aValue)
G4double theNonIonizingEnergyDeposit

Here is the call graph for this function:

G4Step * G4ParticleChangeForLoss::UpdateStepForPostStep ( G4Step Step)
virtual

Reimplemented from G4VParticleChange.

Definition at line 244 of file G4ParticleChangeForLoss.cc.

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 
264  pStep->AddTotalEnergyDeposit( theLocalEnergyDeposit );
265  pStep->AddNonIonizingEnergyDeposit( theNonIonizingEnergyDeposit );
266  return pStep;
267 }
void SetWeight(G4double aValue)
void SetMomentumDirection(const G4ThreeVector &aValue)
void SetPolarization(const G4ThreeVector &aValue)
G4double CalculateVelocity() const
Definition: G4Track.cc:222
void SetVelocity(G4double v)
void SetCharge(G4double value)
void SetKineticEnergy(const G4double aValue)
void SetKineticEnergy(const G4double aValue)
G4double theNonIonizingEnergyDeposit

Here is the call graph for this function:


The documentation for this class was generated from the following files: