Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:59:26

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 // G4VParticleChange inline methods implementation
0027 //
0028 // Author: Hisaya Kurashige, 23 March 1998
0029 // --------------------------------------------------------------------
0030 
0031 inline G4Track* G4VParticleChange::GetSecondary(G4int anIndex) const
0032 {
0033   return theListOfSecondaries[anIndex];
0034 }
0035 
0036 inline G4int G4VParticleChange::GetNumberOfSecondaries() const
0037 {
0038   return theNumberOfSecondaries;
0039 }
0040 
0041 inline void G4VParticleChange::ProposeTrackStatus(G4TrackStatus aStatus)
0042 {
0043   theStatusChange = aStatus;
0044 }
0045 
0046 inline const G4Track* G4VParticleChange::GetCurrentTrack() const
0047 {
0048   return theCurrentTrack;
0049 }
0050 
0051 inline G4TrackStatus G4VParticleChange::GetTrackStatus() const
0052 {
0053   return theStatusChange;
0054 }
0055 
0056 inline G4SteppingControl G4VParticleChange::GetSteppingControl() const
0057 {
0058   return theSteppingControlFlag;
0059 }
0060 
0061 inline void
0062 G4VParticleChange::ProposeSteppingControl(G4SteppingControl StepControlFlag)
0063 {
0064   theSteppingControlFlag = StepControlFlag;
0065 }
0066 
0067 inline G4bool G4VParticleChange::GetFirstStepInVolume() const
0068 {
0069   return theFirstStepInVolume;
0070 }
0071 
0072 inline G4bool G4VParticleChange::GetLastStepInVolume() const
0073 {
0074   return theLastStepInVolume;
0075 }
0076 
0077 inline void G4VParticleChange::ProposeFirstStepInVolume(G4bool flag)
0078 {
0079   theFirstStepInVolume = flag;
0080 }
0081 
0082 inline void G4VParticleChange::ProposeLastStepInVolume(G4bool flag)
0083 {
0084   theLastStepInVolume = flag;
0085 }
0086 
0087 inline G4double G4VParticleChange::GetLocalEnergyDeposit() const
0088 {
0089   return theLocalEnergyDeposit;
0090 }
0091 
0092 inline void G4VParticleChange::ProposeLocalEnergyDeposit(G4double anEnergyPart)
0093 {
0094   theLocalEnergyDeposit = anEnergyPart;
0095 }
0096 
0097 inline G4double G4VParticleChange::GetNonIonizingEnergyDeposit() const
0098 {
0099   return theNonIonizingEnergyDeposit;
0100 }
0101 
0102 inline void
0103 G4VParticleChange::ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
0104 {
0105   theNonIonizingEnergyDeposit = anEnergyPart;
0106 }
0107 
0108 inline G4double G4VParticleChange::GetTrueStepLength() const
0109 {
0110   return theTrueStepLength;
0111 }
0112 
0113 inline void G4VParticleChange::ProposeTrueStepLength(G4double aLength)
0114 {
0115   theTrueStepLength = aLength;
0116 }
0117 
0118 inline void G4VParticleChange::SetVerboseLevel(G4int vLevel)
0119 {
0120   verboseLevel = vLevel;
0121 }
0122 
0123 inline G4int G4VParticleChange::GetVerboseLevel() const
0124 {
0125   return verboseLevel;
0126 }
0127 
0128 inline G4double G4VParticleChange::GetParentWeight() const
0129 {
0130   return theParentWeight;
0131 }
0132 
0133 inline G4double G4VParticleChange::GetWeight() const
0134 {
0135   return theParentWeight;
0136 }
0137 
0138 // ----------------------------------------------------------------
0139 // inline functions for Initialization
0140 //
0141 inline void G4VParticleChange::InitializeLocalEnergyDeposit()
0142 {
0143   // clear theLocalEnergyDeposited
0144   theLocalEnergyDeposit       = 0.0;
0145   theNonIonizingEnergyDeposit = 0.0;
0146 }
0147 
0148 inline void G4VParticleChange::InitializeSteppingControl()
0149 {
0150   // SteppingControlFlag
0151   theSteppingControlFlag = NormalCondition;
0152 }
0153 
0154 inline void G4VParticleChange::Clear()
0155 {
0156   theNumberOfSecondaries = 0;
0157   theFirstStepInVolume   = false;
0158   theLastStepInVolume    = false;
0159 }
0160 
0161 inline void G4VParticleChange::InitializeStatusChange(const G4Track& track)
0162 {
0163   // set TrackStatus equal to the parent track's one
0164   theStatusChange = track.GetTrackStatus();
0165   theCurrentTrack = &track;
0166 }
0167 
0168 inline void G4VParticleChange::InitializeParentWeight(const G4Track& track)
0169 {
0170   // set the parent track's weight
0171   theParentWeight        = track.GetWeight();
0172   isParentWeightProposed = false;
0173 }
0174 
0175 inline void G4VParticleChange::InitializeFromStep(const G4Step* step)
0176 {
0177   // set the parent track's global time at the pre-step point
0178   theParentGlobalTime = step->GetPreStepPoint()->GetGlobalTime();
0179 
0180   theTrueStepLength = step->GetStepLength();
0181   theFirstStepInVolume = step->IsFirstStepInVolume();
0182   theLastStepInVolume = step->IsLastStepInVolume();
0183 }
0184 
0185 // ----------------------------------------------------------------
0186 // methods for handling secondaries
0187 //
0188 inline void G4VParticleChange::InitializeSecondaries()
0189 {
0190   theNumberOfSecondaries = 0;
0191 }
0192 
0193 inline void G4VParticleChange::SetNumberOfSecondaries(G4int totSecondaries)
0194 {
0195   if(totSecondaries > theSizeOftheListOfSecondaries)
0196   {
0197     theListOfSecondaries.resize(totSecondaries, nullptr);
0198     theSizeOftheListOfSecondaries = totSecondaries;
0199   }
0200   theNumberOfSecondaries = 0;
0201 }
0202 
0203 // ----------------------------------------------------------------
0204 // total 
0205 //
0206 inline void G4VParticleChange::Initialize(const G4Track& track)
0207 {
0208   InitializeStatusChange(track);
0209   InitializeLocalEnergyDeposit();
0210   InitializeSteppingControl();
0211   InitializeSecondaries();
0212   InitializeParentWeight(track);
0213   InitializeFromStep(track.GetStep());
0214 }
0215 
0216 inline void G4VParticleChange::ClearDebugFlag()
0217 {
0218   debugFlag = false;
0219 }
0220 
0221 inline void G4VParticleChange::SetDebugFlag()
0222 {
0223   debugFlag = true;
0224 }
0225 
0226 inline G4bool G4VParticleChange::GetDebugFlag() const
0227 {
0228   return debugFlag;
0229 }
0230 
0231 inline void G4VParticleChange::SetSecondaryWeightByProcess(G4bool flag)
0232 {
0233   fSetSecondaryWeightByProcess = flag;
0234 }
0235 
0236 inline G4bool G4VParticleChange::IsSecondaryWeightSetByProcess() const
0237 {
0238   return fSetSecondaryWeightByProcess;
0239 }
0240 
0241 inline void G4VParticleChange::ProposeWeight(G4double w)
0242 {
0243   theParentWeight        = w;
0244   isParentWeightProposed = true;
0245 }
0246 
0247 inline void G4VParticleChange::ProposeParentWeight(G4double w)
0248 {
0249   ProposeWeight(w);
0250 }
0251 
0252 inline G4double G4VParticleChange::ComputeBeta(G4double ekin)
0253 {
0254   G4double mass = theCurrentTrack->GetDefinition()->GetPDGMass();
0255   return std::sqrt(ekin*(ekin + 2*mass))/(ekin + mass);
0256 }