// // ******************************************************************** // * License and Disclaimer * // * * // * The Geant4 software is copyright of the Copyright Holders of * // * the Geant4 Collaboration. It is provided under the terms and * // * conditions of the Geant4 Software License, included in the file * // * LICENSE and available at http://cern.ch/geant4/license . These * // * include a list of copyright holders. * // * * // * Neither the authors of this software system, nor their employing * // * institutes,nor the agencies providing financial support for this * // * work make any representation or warranty, express or implied, * // * regarding this software system or assume any liability for its * // * use. Please see the license in the file LICENSE and URL above * // * for the full disclaimer and the limitation of liability. * // * * // * This code implementation is the result of the scientific and * // * technical work of the GEANT4 collaboration. * // * By using, copying, modifying or distributing the software (or * // * any work based on the software) you agree to acknowledge its * // * use in resulting scientific publications, and indicate your * // * acceptance of all terms of the Geant4 Software license. * // ******************************************************************** // // // $Id: G4VParticleChange.icc,v 1.16 2007/03/25 22:54:52 kurasige Exp $ // GEANT4 tag $Name: geant4-09-02-ref-02 $ // // remove obsolete methods of SetXXX 19 Sep, 04 H.Kurashige //---------------------------------------------------------------- // Virtual methods for updating G4Step // inline G4Step* G4VParticleChange::UpdateStepInfo(G4Step* pStep) { // Update the G4Step specific attributes pStep->SetStepLength( theTrueStepLength ); pStep->AddTotalEnergyDeposit( theLocalEnergyDeposit ); pStep->AddNonIonizingEnergyDeposit( theNonIonizingEnergyDeposit ); pStep->SetControlFlag( theSteppingControlFlag ); if (theFirstStepInVolume) {pStep->SetFirstStepFlag();} else {pStep->ClearFirstStepFlag();} if (theLastStepInVolume) {pStep->SetLastStepFlag();} else {pStep->ClearLastStepFlag();} return pStep; } inline G4Step* G4VParticleChange::UpdateStepForAtRest(G4Step* Step) { if (!fSetParentWeightByProcess){ Step->GetPostStepPoint()->SetWeight( theParentWeight ); } return UpdateStepInfo(Step); } inline G4Step* G4VParticleChange::UpdateStepForAlongStep(G4Step* Step) { if (!fSetParentWeightByProcess){ G4double newWeight= theParentWeight/(Step->GetPreStepPoint()->GetWeight())*(Step->GetPostStepPoint()->GetWeight()); Step->GetPostStepPoint()->SetWeight( newWeight ); } return UpdateStepInfo(Step); } inline G4Step* G4VParticleChange::UpdateStepForPostStep(G4Step* Step) { if (!fSetParentWeightByProcess){ Step->GetPostStepPoint()->SetWeight( theParentWeight ); } return UpdateStepInfo(Step); } //---------------------------------------------------------------- // Set/Get inline functions // inline G4Track* G4VParticleChange::GetSecondary(G4int anIndex) const { return (*theListOfSecondaries)[anIndex]; } inline G4int G4VParticleChange::GetNumberOfSecondaries() const { return theNumberOfSecondaries; } inline void G4VParticleChange::ProposeTrackStatus(G4TrackStatus aStatus) { theStatusChange = aStatus; } inline G4TrackStatus G4VParticleChange::GetTrackStatus() const { return theStatusChange; } inline G4SteppingControl G4VParticleChange::GetSteppingControl() const { return theSteppingControlFlag; } inline void G4VParticleChange::ProposeSteppingControl(G4SteppingControl StepControlFlag) { theSteppingControlFlag = StepControlFlag; } inline G4bool G4VParticleChange::GetFirstStepInVolume() const { return theFirstStepInVolume; } inline G4bool G4VParticleChange::GetLastStepInVolume() const { return theLastStepInVolume; } inline void G4VParticleChange::ProposeFirstStepInVolume(G4bool flag) { theFirstStepInVolume = flag; } inline void G4VParticleChange::ProposeLastStepInVolume(G4bool flag) { theLastStepInVolume = flag; } //---------------------------------------------------------------- // Set/Get inline functions // inline G4double G4VParticleChange::GetLocalEnergyDeposit() const { return theLocalEnergyDeposit; } inline void G4VParticleChange::ProposeLocalEnergyDeposit(G4double anEnergyPart) { theLocalEnergyDeposit = anEnergyPart; } inline G4double G4VParticleChange::GetNonIonizingEnergyDeposit() const { return theNonIonizingEnergyDeposit; } inline void G4VParticleChange::ProposeNonIonizingEnergyDeposit(G4double anEnergyPart) { theNonIonizingEnergyDeposit = anEnergyPart; } inline G4double G4VParticleChange::GetTrueStepLength() const { return theTrueStepLength; } inline void G4VParticleChange::ProposeTrueStepLength(G4double aLength) { theTrueStepLength = aLength; } inline void G4VParticleChange::SetVerboseLevel(G4int vLevel) { verboseLevel = vLevel; } inline G4int G4VParticleChange::GetVerboseLevel() const { return verboseLevel; } inline G4double G4VParticleChange::GetParentWeight() const { return theParentWeight; } //---------------------------------------------------------------- // inline functions for Initialization // inline void G4VParticleChange::InitializeLocalEnergyDeposit(const G4Track&) { // clear theLocalEnergyDeposited theLocalEnergyDeposit = 0.0; theNonIonizingEnergyDeposit = 0.0; } inline void G4VParticleChange::InitializeSteppingControl(const G4Track& ) { // SteppingControlFlag theSteppingControlFlag = NormalCondition; } inline void G4VParticleChange::Clear() { theNumberOfSecondaries = 0; theFirstStepInVolume = false; theLastStepInVolume = false; } //---------------------------------------------------------------- // functions for Initialization // inline void G4VParticleChange::InitializeStatusChange(const G4Track& track) { // set TrackStatus equal to the parent track's one theStatusChange = track.GetTrackStatus(); } inline void G4VParticleChange::InitializeParentWeight(const G4Track& track) { // set the parent track's weight theParentWeight = track.GetWeight(); } inline void G4VParticleChange::InitializeTrueStepLength(const G4Track& track) { // Reset theTrueStepLength // !! Caution !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! theTrueStepLength = track.GetStep()->GetStepLength(); // !! TrueStepLength should be copied from G4Step not G4Track // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! } //---------------------------------------------------------------- // methods for initialize inline void G4VParticleChange::InitializeStepInVolumeFlags(const G4Track& track) { const G4Step* aStep = track.GetStep(); theFirstStepInVolume = aStep-> IsFirstStepInVolume(); theLastStepInVolume = aStep-> IsLastStepInVolume(); } inline void G4VParticleChange::InitializeSecondaries(const G4Track&) { // clear secondaries if (theNumberOfSecondaries>0) { #ifdef G4VERBOSE if (verboseLevel>0) { G4cerr << "G4VParticleChange::Initialize() Warning "; G4cerr << "theListOfSecondaries is not empty " << G4endl; G4cerr << "All objects in theListOfSecondaries are destroyed!" << G4endl; } #endif for (G4int index= 0; index0) { #ifdef G4VERBOSE if (verboseLevel>0) { G4cerr << "G4VParticleChange::SetNumberOfSecondaries() Warning "; G4cerr << "theListOfSecondaries is not empty "; } #endif for (G4int index= 0; indexInitialize(totSecondaries); } inline void G4VParticleChange::Initialize(const G4Track& track) { InitializeStatusChange(track); InitializeLocalEnergyDeposit(track); InitializeSteppingControl(track); InitializeTrueStepLength(track); InitializeSecondaries(track); InitializeParentWeight(track); InitializeStepInVolumeFlags(track); } inline void G4VParticleChange::ClearDebugFlag() { debugFlag = false; } inline void G4VParticleChange::SetDebugFlag() { debugFlag = true; } inline G4bool G4VParticleChange::GetDebugFlag() const { return debugFlag; } inline void G4VParticleChange::SetSecondaryWeightByProcess(G4bool flag) { fSetSecondaryWeightByProcess = flag; } inline G4bool G4VParticleChange::IsSecondaryWeightSetByProcess() const { return fSetSecondaryWeightByProcess; } inline void G4VParticleChange::ProposeParentWeight(G4double w) { theParentWeight = w; } inline void G4VParticleChange::SetParentWeightByProcess(G4bool flag) { fSetParentWeightByProcess = flag; } inline G4bool G4VParticleChange::IsParentWeightSetByProcess() const { return fSetParentWeightByProcess; }