Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4StepPoint.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// $Id$
28//
29//
30//---------------------------------------------------------------
31//
32// G4StepPoint.hh
33//
34// Class Description:
35// This class represents information associated with the
36// each end of a Step like the space/time data of the
37// particle.
38//
39// Contact:
40// Questions and comments to this code should be sent to
41// Hisaya Kurashige
42//
43// ---------------------------------------------------------------
44// Added fpMaterial 16 FEb. 2000 H.Kurahige
45// Added fpMaterialCutsCouple 8 Oct. 2002 H.Kurahige
46// Added fMagneticMoment Mar 2007 H.Kurashige
47
48#ifndef G4StepPoint_h
49#define G4StepPoint_h 1
50
51#include <cmath> // Include from 'system'
53
54#include "globals.hh" // Include from 'global'
55#include "G4Allocator.hh" // Include from 'global'
56#include "G4ThreeVector.hh" // Include from 'geometry'
57#include "G4VPhysicalVolume.hh" // Include from 'geometry'
58class G4VProcess;
59#include "G4SteppingControl.hh"
60#include "G4StepStatus.hh" // Include from 'track'
61#include "G4TouchableHandle.hh" // Include from 'geometry'
62#include "G4Material.hh"
63#include "G4LogicalVolume.hh"
64
67/////////////////
69/////////////////
70{
71
72//--------
73 public:
74
75
76// Constructor/Destructor
78
80
81// Copy Counstructor and assignment operator
82 G4StepPoint(const G4StepPoint& );
84
85//--------
86
87 public: // with description
88
89// Get/Set functions
90 const G4ThreeVector& GetPosition() const;
91 void SetPosition(const G4ThreeVector& aValue);
92 void AddPosition(const G4ThreeVector& aValue);
93
95 void SetLocalTime(const G4double aValue);
96 void AddLocalTime(const G4double aValue);
97 // Time since the track is created.
98
100 void SetGlobalTime(const G4double aValue);
101 void AddGlobalTime(const G4double aValue);
102 // Time since the event in which the track belongs is created.
103
105 void SetProperTime(const G4double aValue);
106 void AddProperTime(const G4double aValue);
107 // Proper time of the particle.
108
112 // Direction of momentum (should be an unit vector)
113
115 // Total momentum of the track
116
117
119 // Total energy of the track
120
122 void SetKineticEnergy(const G4double aValue);
123 void AddKineticEnergy(const G4double aValue);
124 // Kinetic Energy of the track
125
128 //
129
131 // Velocity of the track in unit of c(light velocity)
132
134 // Gamma factor (1/sqrt[1-beta*beta]) of the track
135
137
141
144
147
150
152 void SetSafety(const G4double aValue);
153
155 void SetPolarization(const G4ThreeVector& aValue);
156 void AddPolarization(const G4ThreeVector& aValue);
157
159 void SetStepStatus(const G4StepStatus aValue);
160
162 // If the pointer is 0, this means the Step is defined
163 // by the user defined limit in the current volume.
165
166
168 void SetMass(G4double value);
169
171 void SetCharge(G4double value);
172
175
176 void SetWeight(G4double aValue);
178
179//---------
180 private:
181//---------
182
183// Member data
184 G4ThreeVector fPosition;
185 G4double fGlobalTime;
186 // Time since event is created
187 G4double fLocalTime;
188 // Time since track is created
189 G4double fProperTime;
190 // Time since track is created (in rest frame of particle)
191 G4ThreeVector fMomentumDirection;
192 G4double fKineticEnergy;
193 G4double fVelocity;
194 // Momentum,energy and velocity
195 G4TouchableHandle fpTouchable;
196 // Touchable Handle
197 G4Material* fpMaterial;
198 // Material of the volmue
199 const G4MaterialCutsCouple* fpMaterialCutsCouple;
200 // MaterialCutsCouple of the volmue
201 G4VSensitiveDetector* fpSensitiveDetector;
202 G4double fSafety;
203 G4ThreeVector fPolarization;
204 G4StepStatus fStepStatus;
205 // DoIt type which defined the current Step.
206 const G4VProcess* fpProcessDefinedStep;
207 // Process which defined the current Step.
208 G4double fMass;
209 // Dynamical mass of the particle
210 G4double fCharge;
211 // Dynamical Charge of the particle
212 G4double fMagneticMoment;
213 // Dynamical MagneticMoment of the particle
214 G4double fWeight;
215 // Track Weight
216};
217
218#include "G4StepPoint.icc"
219
220#endif
G4StepStatus
Definition: G4StepStatus.hh:51
double G4double
Definition: G4Types.hh:64
void AddPolarization(const G4ThreeVector &aValue)
G4double GetMass() const
void SetLocalTime(const G4double aValue)
G4double GetTotalEnergy() const
void SetMagneticMoment(G4double value)
G4double GetMagneticMoment() const
void SetSensitiveDetector(G4VSensitiveDetector *)
void SetKineticEnergy(const G4double aValue)
void SetWeight(G4double aValue)
void SetMaterial(G4Material *)
void AddKineticEnergy(const G4double aValue)
G4StepStatus GetStepStatus() const
void SetMass(G4double value)
G4double GetVelocity() const
void SetSafety(const G4double aValue)
void SetCharge(G4double value)
G4double GetProperTime() const
G4double GetBeta() const
G4StepPoint & operator=(const G4StepPoint &)
Definition: G4StepPoint.cc:85
const G4VTouchable * GetTouchable() const
void SetStepStatus(const G4StepStatus aValue)
G4double GetGlobalTime() const
void AddProperTime(const G4double aValue)
void SetVelocity(G4double v)
G4double GetCharge() const
G4double GetSafety() const
const G4VProcess * GetProcessDefinedStep() const
G4double GetGamma() const
void AddPosition(const G4ThreeVector &aValue)
void SetProcessDefinedStep(const G4VProcess *aValue)
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4Material * GetMaterial() const
G4ThreeVector GetMomentum() const
const G4ThreeVector & GetPosition() const
void SetProperTime(const G4double aValue)
const G4ThreeVector & GetMomentumDirection() const
const G4TouchableHandle & GetTouchableHandle() const
void SetMaterialCutsCouple(const G4MaterialCutsCouple *)
void AddGlobalTime(const G4double aValue)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4double GetLocalTime() const
G4VSensitiveDetector * GetSensitiveDetector() const
void SetGlobalTime(const G4double aValue)
G4VPhysicalVolume * GetPhysicalVolume() const
void SetPosition(const G4ThreeVector &aValue)
const G4ThreeVector & GetPolarization() const
G4double GetKineticEnergy() const
G4double GetWeight() const
void AddLocalTime(const G4double aValue)
void SetMomentumDirection(const G4ThreeVector &aValue)
void SetPolarization(const G4ThreeVector &aValue)
void AddMomentumDirection(const G4ThreeVector &aValue)