// // ******************************************************************** // * 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: G4VRestDiscreteProcess.hh,v 1.9 2007/11/15 04:10:18 kurasige Exp $ // GEANT4 tag $Name: geant4-09-04-beta-01 $ // // // ------------------------------------------------------------ // GEANT 4 class header file // // // Class Description // Abstract class which defines the public behavior of // rest + discrete physics interactions. // // ------------------------------------------------------------ // New Physics scheme 8 Mar. 1997 H.Kurahige // ------------------------------------------------------------ // modified for new ParticleChange 12 Mar. 1998 H.Kurashige // Fixed a bug in PostStepGetPhysicalInteractionLength // 15 Apr. 2002 H.Kurashige // #ifndef G4VRestDiscreteProcess_h #define G4VRestDiscreteProcess_h 1 #include "globals.hh" #include "G4ios.hh" #include "G4VProcess.hh" class G4VRestDiscreteProcess : public G4VProcess { // Abstract class which defines the public behavior of // rest + discrete physics interactions. public: G4VRestDiscreteProcess(const G4String& , G4ProcessType aType = fNotDefined ); G4VRestDiscreteProcess(G4VRestDiscreteProcess &); virtual ~G4VRestDiscreteProcess(); public :// with description virtual G4double PostStepGetPhysicalInteractionLength( const G4Track& track, G4double previousStepSize, G4ForceCondition* condition ); virtual G4VParticleChange* PostStepDoIt( const G4Track& , const G4Step& ); virtual G4double AtRestGetPhysicalInteractionLength( const G4Track& , G4ForceCondition* ); virtual G4VParticleChange* AtRestDoIt( const G4Track& , const G4Step& ); // no operation in AlongStepDoIt virtual G4double AlongStepGetPhysicalInteractionLength( const G4Track&, G4double , G4double , G4double& , G4GPILSelection* ){ return -1.0; }; // no operation in AlongStepDoIt virtual G4VParticleChange* AlongStepDoIt( const G4Track& , const G4Step& ) {return 0;}; protected:// with description virtual G4double GetMeanFreePath(const G4Track& aTrack, G4double previousStepSize, G4ForceCondition* condition )=0; // Calculates from the macroscopic cross section a mean // free path, the value is returned in units of distance. virtual G4double GetMeanLifeTime(const G4Track& aTrack,G4ForceCondition* condition)=0; // Calculates the mean life-time (i.e. for decays) of the // particle at rest due to the occurence of the given process, // or converts the probability of interaction (i.e. for // annihilation) into the life-time of the particle for the // occurence of the given process. private: // hide default constructor and assignment operator as private G4VRestDiscreteProcess(); G4VRestDiscreteProcess & operator=(const G4VRestDiscreteProcess &right); }; // ----------------------------------------- // inlined function members implementation // ----------------------------------------- #include "Randomize.hh" #include "G4Step.hh" #include "G4Track.hh" #include "G4MaterialTable.hh" #include "G4VParticleChange.hh" inline G4double G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength( const G4Track& track, G4double previousStepSize, G4ForceCondition* condition ) { if ( (previousStepSize < 0.0) || (theNumberOfInteractionLengthLeft<=0.0)) { // beggining of tracking (or just after DoIt of this process) ResetNumberOfInteractionLengthLeft(); } else if ( previousStepSize > 0.0) { // subtract NumberOfInteractionLengthLeft SubtractNumberOfInteractionLengthLeft(previousStepSize); } else { // zero step // DO NOTHING } // condition is set to "Not Forced" *condition = NotForced; // get mean free path currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition); G4double value; if (currentInteractionLength 1){ G4cout << "G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength "; G4cout << "[ " << GetProcessName() << "]" <DumpInfo(); G4cout << " in Material " << track.GetMaterial()->GetName() <2)){ G4cout << "G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength "; G4cout << "[ " << GetProcessName() << "]" <DumpInfo(); G4cout << " in Material " << track.GetMaterial()->GetName() <