Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions
G4ParticleChangeForMSC Class Reference

#include <G4ParticleChangeForMSC.hh>

Inheritance diagram for G4ParticleChangeForMSC:
G4VParticleChange

Public Member Functions

 G4ParticleChangeForMSC ()
 
virtual ~G4ParticleChangeForMSC ()
 
virtual G4StepUpdateStepForAlongStep (G4Step *Step)
 
virtual G4StepUpdateStepForPostStep (G4Step *Step)
 
virtual void Initialize (const G4Track &)
 
void ProposeMomentumDirection (const G4ThreeVector &Pfinal)
 
void ProposeMomentumDirection (G4double Px, G4double Py, G4double Pz)
 
const G4ThreeVectorGetMomentumDirection () const
 
const G4ThreeVectorGetProposedMomentumDirection () const
 
void SetProposedMomentumDirection (const G4ThreeVector &Pfinal)
 
const G4ThreeVectorGetPosition () const
 
void ProposePosition (const G4ThreeVector &finalPosition)
 
const G4ThreeVectorGetProposedPosition () const
 
void SetProposedPosition (const G4ThreeVector &finalPosition)
 
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)
 
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

 G4ParticleChangeForMSC (const G4ParticleChangeForMSC &right)
 
G4ParticleChangeForMSCoperator= (const G4ParticleChangeForMSC &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 54 of file G4ParticleChangeForMSC.hh.

Constructor & Destructor Documentation

G4ParticleChangeForMSC::G4ParticleChangeForMSC ( )

Definition at line 47 of file G4ParticleChangeForMSC.cc.

References G4cout, G4endl, and G4VParticleChange::verboseLevel.

49 {
50 #ifdef G4VERBOSE
51  if (verboseLevel>2) {
52  G4cout << "G4ParticleChangeForMSC::G4ParticleChangeForMSC() " << G4endl;
53  }
54 #endif
55 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
G4ParticleChangeForMSC::~G4ParticleChangeForMSC ( )
virtual

Definition at line 57 of file G4ParticleChangeForMSC.cc.

References G4cout, G4endl, and G4VParticleChange::verboseLevel.

58 {
59 #ifdef G4VERBOSE
60  if (verboseLevel>2) {
61  G4cout << "G4ParticleChangeForMSC::~G4ParticleChangeForMSC() " << G4endl;
62  }
63 #endif
64 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
G4ParticleChangeForMSC::G4ParticleChangeForMSC ( const G4ParticleChangeForMSC right)
protected

Definition at line 67 of file G4ParticleChangeForMSC.cc.

References G4cout, G4endl, and G4VParticleChange::verboseLevel.

68  : G4VParticleChange(right)
69 {
70  if (verboseLevel>1) {
71  G4cout << "G4ParticleChangeForMSC:: copy constructor is called " << G4endl;
72  }
73  theMomentumDirection = right.theMomentumDirection;
74  thePosition = right.thePosition;
75 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Member Function Documentation

G4bool G4ParticleChangeForMSC::CheckIt ( const G4Track aTrack)
virtual

Reimplemented from G4VParticleChange.

Definition at line 169 of file G4ParticleChangeForMSC.cc.

References G4VParticleChange::accuracyForException, G4VParticleChange::accuracyForWarning, G4VParticleChange::CheckIt(), DumpInfo(), EventMustBeAborted, G4cout, G4endl, G4Exception(), G4Track::GetDefinition(), G4Track::GetKineticEnergy(), G4ParticleDefinition::GetParticleName(), G4Track::GetPosition(), python.hepunit::m, CLHEP::Hep3Vector::mag(), CLHEP::Hep3Vector::mag2(), python.hepunit::MeV, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

170 {
171  G4bool itsOK = true;
172  G4bool exitWithError = false;
173 
174  G4double accuracy;
175 
176  // check
177 
178  // MomentumDirection should be unit vector
179  accuracy = std::fabs(theMomentumDirection.mag2()-1.0);
180  if (accuracy > accuracyForWarning) {
181  itsOK = false;
182  exitWithError = (accuracy > accuracyForException);
183 #ifdef G4VERBOSE
184  G4cout << " G4ParticleChangeForMSC::CheckIt : ";
185  G4cout << "the Momentum Change is not unit vector !!"
186  << " Difference: " << accuracy << G4endl;
187  G4cout << aTrack.GetDefinition()->GetParticleName()
188  << " E=" << aTrack.GetKineticEnergy()/MeV
189  << " pos=" << aTrack.GetPosition().x()/m
190  << ", " << aTrack.GetPosition().y()/m
191  << ", " << aTrack.GetPosition().z()/m
192  <<G4endl;
193 #endif
194  }
195 
196  // dump out information of this particle change
197 #ifdef G4VERBOSE
198  if (!itsOK) DumpInfo();
199 #endif
200 
201  // Exit with error
202  if (exitWithError) {
203  G4Exception("G4ParticleChangeForMSC::CheckIt",
204  "300",
206  "momentum direction was illegal");
207  }
208  //correction
209  if (!itsOK) {
210  G4double vmag = theMomentumDirection.mag();
211  theMomentumDirection = (1./vmag)*theMomentumDirection;
212  }
213 
214  itsOK = (itsOK) && G4VParticleChange::CheckIt(aTrack);
215  return itsOK;
216 }
G4ParticleDefinition * GetDefinition() 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
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
double mag2() const
virtual void DumpInfo() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static const G4double accuracyForWarning
double mag() const
void G4ParticleChangeForMSC::DumpInfo ( ) const
virtual

Reimplemented from G4VParticleChange.

Definition at line 141 of file G4ParticleChangeForMSC.cc.

References G4VParticleChange::DumpInfo(), G4cout, G4endl, python.hepunit::mm, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by CheckIt().

142 {
143 // use base-class DumpInfo
145 
146  G4int oldprc = G4cout.precision(3);
147  G4cout << " Position - x (mm) : "
148  << std::setw(20) << thePosition.x()/mm
149  << G4endl;
150  G4cout << " Position - y (mm) : "
151  << std::setw(20) << thePosition.y()/mm
152  << G4endl;
153  G4cout << " Position - z (mm) : "
154  << std::setw(20) << thePosition.z()/mm
155  << G4endl;
156  G4cout << " Momentum Direct - x : "
157  << std::setw(20) << theMomentumDirection.x()
158  << G4endl;
159  G4cout << " Momentum Direct - y : "
160  << std::setw(20) << theMomentumDirection.y()
161  << G4endl;
162  G4cout << " Momentum Direct - z : "
163  << std::setw(20) << theMomentumDirection.z()
164  << G4endl;
165  G4cout.precision(oldprc);
166 }
double x() const
virtual void DumpInfo() const
int G4int
Definition: G4Types.hh:78
double z() const
G4GLOB_DLL std::ostream G4cout
double y() const
#define G4endl
Definition: G4ios.hh:61
const G4ThreeVector* G4ParticleChangeForMSC::GetMomentumDirection ( ) const
const G4ThreeVector* G4ParticleChangeForMSC::GetPosition ( ) const
const G4ThreeVector* G4ParticleChangeForMSC::GetProposedMomentumDirection ( ) const
const G4ThreeVector* G4ParticleChangeForMSC::GetProposedPosition ( ) const
virtual void G4ParticleChangeForMSC::Initialize ( const G4Track )
virtual

Reimplemented from G4VParticleChange.

Referenced by G4VMultipleScattering::PostStepDoIt().

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

Definition at line 78 of file G4ParticleChangeForMSC.cc.

References G4cout, G4endl, G4VParticleChange::theListOfSecondaries, G4VParticleChange::theLocalEnergyDeposit, G4VParticleChange::theNumberOfSecondaries, G4VParticleChange::theSizeOftheListOfSecondaries, G4VParticleChange::theStatusChange, G4VParticleChange::theSteppingControlFlag, G4VParticleChange::theTrueStepLength, and G4VParticleChange::verboseLevel.

80 {
81  if (verboseLevel>1) {
82  G4cout << "G4ParticleChangeForMSC:: assignment operator is called " << G4endl;
83  }
84  if (this != &right)
85  {
93 
94  theMomentumDirection = right.theMomentumDirection;
95  thePosition = right.thePosition;
96  }
97  return *this;
98 }
G4TrackFastVector * theListOfSecondaries
G4GLOB_DLL std::ostream G4cout
G4SteppingControl theSteppingControlFlag
#define G4endl
Definition: G4ios.hh:61
G4TrackStatus theStatusChange
void G4ParticleChangeForMSC::ProposeMomentumDirection ( const G4ThreeVector Pfinal)
void G4ParticleChangeForMSC::ProposeMomentumDirection ( G4double  Px,
G4double  Py,
G4double  Pz 
)
void G4ParticleChangeForMSC::ProposePosition ( const G4ThreeVector finalPosition)
void G4ParticleChangeForMSC::SetProposedMomentumDirection ( const G4ThreeVector Pfinal)
void G4ParticleChangeForMSC::SetProposedPosition ( const G4ThreeVector finalPosition)
G4Step * G4ParticleChangeForMSC::UpdateStepForAlongStep ( G4Step Step)
virtual

Reimplemented from G4VParticleChange.

Definition at line 104 of file G4ParticleChangeForMSC.cc.

References G4Step::GetPostStepPoint(), G4Step::GetTrack(), G4Track::GetTrackStatus(), G4StepPoint::SetMomentumDirection(), G4StepPoint::SetPosition(), G4Step::SetStepLength(), G4VParticleChange::theStatusChange, and G4VParticleChange::theTrueStepLength.

105 {
106  // Update the G4Step specific attributes
107  pStep->SetStepLength(theTrueStepLength);
108  theStatusChange = pStep->GetTrack()->GetTrackStatus();
109 
110  // Multiple scattering calculates the final state of the particle
111  G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
112 
113  // update momentum direction
114  pPostStepPoint->SetMomentumDirection(theMomentumDirection);
115 
116  // update position
117  pPostStepPoint->SetPosition( thePosition );
118  return pStep;
119 }
void SetPosition(const G4ThreeVector &aValue)
void SetMomentumDirection(const G4ThreeVector &aValue)
G4TrackStatus theStatusChange
G4Step * G4ParticleChangeForMSC::UpdateStepForPostStep ( G4Step Step)
virtual

Reimplemented from G4VParticleChange.

Definition at line 121 of file G4ParticleChangeForMSC.cc.

References G4Step::GetPostStepPoint(), G4StepPoint::SetMomentumDirection(), and G4StepPoint::SetPosition().

122 {
123  // Multiple scattering calculates the final state of the particle
124  G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
125 
126  // update momentum direction
127  pPostStepPoint->SetMomentumDirection(theMomentumDirection);
128 
129  // update position
130  pPostStepPoint->SetPosition( thePosition );
131 
132  // Update the G4Step specific attributes
133  return pStep;
134 }
void SetPosition(const G4ThreeVector &aValue)
void SetMomentumDirection(const G4ThreeVector &aValue)

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