G4UnknownDecay Class Reference

#include <G4UnknownDecay.hh>

Inheritance diagram for G4UnknownDecay:

G4VDiscreteProcess G4VProcess

Public Member Functions

 G4UnknownDecay (const G4String &processName="UnknownDecay")
virtual ~G4UnknownDecay ()
virtual G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
virtual G4bool IsApplicable (const G4ParticleDefinition &)
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
void SetVerboseLevel (G4int value)
G4int GetVerboseLevel () const

Protected Member Functions

virtual G4VParticleChangeDecayIt (const G4Track &aTrack, const G4Step &aStep)
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)

Detailed Description

Definition at line 44 of file G4UnknownDecay.hh.


Constructor & Destructor Documentation

G4UnknownDecay::G4UnknownDecay ( const G4String processName = "UnknownDecay"  ) 

Definition at line 47 of file G4UnknownDecay.cc.

References DECAY_Unknown, G4cout, G4endl, GetVerboseLevel(), G4VProcess::pParticleChange, and G4VProcess::SetProcessSubType().

00048                                :G4VDiscreteProcess(processName, fDecay),
00049                                 verboseLevel(1),
00050                                 HighestValue(20.0)
00051 {
00052   // set Process Sub Type
00053   SetProcessSubType(static_cast<int>(DECAY_Unknown));
00054 
00055 #ifdef G4VERBOSE
00056   if (GetVerboseLevel()>1) {
00057     G4cout << "G4UnknownDecay  constructor " << "  Name:" << processName << G4endl;
00058   }
00059 #endif
00060   pParticleChange = &fParticleChangeForDecay;
00061 }

G4UnknownDecay::~G4UnknownDecay (  )  [virtual]

Definition at line 63 of file G4UnknownDecay.cc.

00064 {
00065 }


Member Function Documentation

void G4UnknownDecay::BuildPhysicsTable ( const G4ParticleDefinition  )  [virtual]

Reimplemented from G4VProcess.

Definition at line 78 of file G4UnknownDecay.cc.

00079 {
00080   return;
00081 }

G4VParticleChange * G4UnknownDecay::DecayIt ( const G4Track aTrack,
const G4Step aStep 
) [protected, virtual]

Definition at line 83 of file G4UnknownDecay.cc.

References G4VParticleChange::AddSecondary(), G4DecayProducts::Boost(), G4VProcess::ClearNumberOfInteractionLengthLeft(), G4DecayProducts::DumpInfo(), G4DecayProducts::entries(), fStopAndKill, G4cout, G4endl, G4DynamicParticle::GetDefinition(), G4Track::GetDynamicParticle(), G4Track::GetGlobalTime(), G4DynamicParticle::GetMass(), G4DynamicParticle::GetMomentumDirection(), G4ParticleDefinition::GetParticleName(), G4Track::GetPosition(), G4DynamicParticle::GetPreAssignedDecayProducts(), G4DynamicParticle::GetPreAssignedDecayProperTime(), G4DynamicParticle::GetTotalEnergy(), G4Track::GetTouchableHandle(), GetVerboseLevel(), G4ParticleChangeForDecay::Initialize(), ns, G4DecayProducts::PopProducts(), G4ParticleChangeForDecay::ProposeGlobalTime(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4VParticleChange::ProposeTrackStatus(), G4Track::SetGoodForTrackingFlag(), G4VParticleChange::SetNumberOfSecondaries(), and G4Track::SetTouchableHandle().

Referenced by PostStepDoIt().

00084 {
00085   // The DecayIt() method returns by pointer a particle-change object.
00086   // Units are expressed in GEANT4 internal units.
00087 
00088   //   Initialize ParticleChange
00089   //     all members of G4VParticleChange are set to equal to 
00090   //     corresponding member in G4Track
00091   fParticleChangeForDecay.Initialize(aTrack);
00092 
00093   // get particle 
00094   const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
00095 
00096   //check if thePreAssignedDecayProducts exists
00097   const G4DecayProducts* o_products = (aParticle->GetPreAssignedDecayProducts());
00098   G4bool isPreAssigned = (o_products != 0);   
00099   G4DecayProducts* products = 0;
00100 
00101   if (!isPreAssigned ){
00102     fParticleChangeForDecay.SetNumberOfSecondaries(0);
00103     // Kill the parent particle
00104     fParticleChangeForDecay.ProposeTrackStatus( fStopAndKill ) ;
00105     fParticleChangeForDecay.ProposeLocalEnergyDeposit(0.0); 
00106     
00107     ClearNumberOfInteractionLengthLeft();
00108     return &fParticleChangeForDecay ;
00109   }
00110 
00111   // copy decay products 
00112   products = new G4DecayProducts(*o_products); 
00113   
00114   // get parent particle information ...................................
00115   G4double   ParentEnergy  = aParticle->GetTotalEnergy();
00116   G4double   ParentMass    = aParticle->GetMass();
00117   if (ParentEnergy < ParentMass) {
00118     ParentEnergy = ParentMass;
00119 #ifdef G4VERBOSE
00120     if (GetVerboseLevel()>1) {
00121       G4cout << "G4UnknownDecay::DoIt  : Total Energy is less than its mass" << G4endl;
00122       G4cout << " Particle: " << aParticle->GetDefinition()->GetParticleName();
00123       G4cout << " Energy:"    << ParentEnergy/MeV << "[MeV]";
00124       G4cout << " Mass:"    << ParentMass/MeV << "[MeV]";
00125       G4cout << G4endl;
00126     }
00127 #endif
00128   }
00129   G4ThreeVector ParentDirection(aParticle->GetMomentumDirection());
00130 
00131   G4double energyDeposit = 0.0;
00132   G4double finalGlobalTime = aTrack.GetGlobalTime();
00133   //boost all decay products to laboratory frame
00134   //if the particle has traveled 
00135   if(aParticle->GetPreAssignedDecayProperTime()>0.) {
00136     products->Boost( ParentEnergy, ParentDirection);
00137   }
00138 
00139   //add products in fParticleChangeForDecay
00140   G4int numberOfSecondaries = products->entries();
00141   fParticleChangeForDecay.SetNumberOfSecondaries(numberOfSecondaries);
00142 #ifdef G4VERBOSE
00143   if (GetVerboseLevel()>1) {
00144     G4cout << "G4UnknownDecay::DoIt  : Decay vertex :";
00145     G4cout << " Time: " << finalGlobalTime/ns << "[ns]";
00146     G4cout << " X:" << (aTrack.GetPosition()).x() /cm << "[cm]";
00147     G4cout << " Y:" << (aTrack.GetPosition()).y() /cm << "[cm]";
00148     G4cout << " Z:" << (aTrack.GetPosition()).z() /cm << "[cm]";
00149     G4cout << G4endl;
00150     G4cout << "G4UnknownDecay::DoIt  : decay products in Lab. Frame" << G4endl;
00151     products->DumpInfo();
00152   }
00153 #endif
00154   G4int index;
00155   G4ThreeVector currentPosition;
00156   const G4TouchableHandle thand = aTrack.GetTouchableHandle();
00157   for (index=0; index < numberOfSecondaries; index++)
00158   {
00159      // get current position of the track
00160      currentPosition = aTrack.GetPosition();
00161      // create a new track object
00162      G4Track* secondary = new G4Track( products->PopProducts(),
00163                                       finalGlobalTime ,
00164                                       currentPosition );
00165      // switch on good for tracking flag
00166      secondary->SetGoodForTrackingFlag();
00167      secondary->SetTouchableHandle(thand);
00168      // add the secondary track in the List
00169      fParticleChangeForDecay.AddSecondary(secondary);
00170   }
00171   delete products;
00172 
00173   // Kill the parent particle
00174   fParticleChangeForDecay.ProposeTrackStatus( fStopAndKill ) ;
00175   fParticleChangeForDecay.ProposeLocalEnergyDeposit(energyDeposit); 
00176   fParticleChangeForDecay.ProposeGlobalTime( finalGlobalTime );
00177   // reset NumberOfInteractionLengthLeft
00178   ClearNumberOfInteractionLengthLeft();
00179 
00180   return &fParticleChangeForDecay ;
00181 } 

G4double G4UnknownDecay::GetMeanFreePath ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
) [protected, virtual]

Implements G4VDiscreteProcess.

Definition at line 73 of file G4UnknownDecay.cc.

References DBL_MIN.

00074 {
00075    return  DBL_MIN;
00076 }

G4int G4UnknownDecay::GetVerboseLevel (  )  const [inline]

Reimplemented from G4VProcess.

Definition at line 154 of file G4UnknownDecay.hh.

Referenced by DecayIt(), and G4UnknownDecay().

00154 { return verboseLevel; }

G4bool G4UnknownDecay::IsApplicable ( const G4ParticleDefinition  )  [virtual]

Reimplemented from G4VProcess.

Definition at line 67 of file G4UnknownDecay.cc.

References G4ParticleDefinition::GetParticleName().

00068 {
00069   if(aParticleType.GetParticleName()=="unknown") return true;
00070   return false;
00071 }

G4VParticleChange * G4UnknownDecay::PostStepDoIt ( const G4Track aTrack,
const G4Step aStep 
) [inline, virtual]

Reimplemented from G4VDiscreteProcess.

Definition at line 157 of file G4UnknownDecay.hh.

References DecayIt().

00161 {
00162   return DecayIt(aTrack, aStep);
00163 }

G4double G4UnknownDecay::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
) [inline, virtual]

Reimplemented from G4VDiscreteProcess.

Definition at line 126 of file G4UnknownDecay.hh.

References DBL_MIN, G4Track::GetDynamicParticle(), G4DynamicParticle::GetPreAssignedDecayProperTime(), G4Track::GetProperTime(), and NotForced.

00131 {
00132   // pre-assigned UnknownDecay time
00133   G4double pTime = track.GetDynamicParticle()->GetPreAssignedDecayProperTime();
00134 
00135   if (pTime < 0.) pTime = DBL_MIN;
00136 
00137   // condition is set to "Not Forced"
00138   *condition = NotForced;
00139   
00140   // reminder proper time
00141   G4double remainder = pTime - track.GetProperTime();
00142   if (remainder <= 0.0) remainder = DBL_MIN;
00143   
00144   // use pre-assigned Decay time to determine PIL
00145   //return GetMeanFreePath(track, previousStepSize, condition);
00146   return remainder*CLHEP::c_light;
00147 
00148 }

void G4UnknownDecay::SetVerboseLevel ( G4int  value  )  [inline]

Reimplemented from G4VProcess.

Definition at line 151 of file G4UnknownDecay.hh.

00151 { verboseLevel = value; }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:53:37 2013 for Geant4 by  doxygen 1.4.7