#include <G4PionMinusAbsorptionAtRest.hh>
Inheritance diagram for G4PionMinusAbsorptionAtRest:
Public Member Functions | |
G4PionMinusAbsorptionAtRest (const G4String &processName="PionMinusAbsorptionAtRest", G4ProcessType aType=fHadronic) | |
~G4PionMinusAbsorptionAtRest () | |
G4bool | IsApplicable (const G4ParticleDefinition &) |
void | PreparePhysicsTable (const G4ParticleDefinition &) |
void | BuildPhysicsTable (const G4ParticleDefinition &) |
G4double | AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *) |
G4double | GetMeanLifeTime (const G4Track &, G4ForceCondition *) |
G4VParticleChange * | AtRestDoIt (const G4Track &, const G4Step &) |
G4int | GetNumberOfSecondaries () |
G4GHEKinematicsVector * | GetSecondaryKinematics () |
Definition at line 46 of file G4PionMinusAbsorptionAtRest.hh.
G4PionMinusAbsorptionAtRest::G4PionMinusAbsorptionAtRest | ( | const G4String & | processName = "PionMinusAbsorptionAtRest" , |
|
G4ProcessType | aType = fHadronic | |||
) |
Definition at line 46 of file G4PionMinusAbsorptionAtRest.cc.
References fHadronAtRest, G4cout, G4endl, G4HadronicDeprecate, G4VProcess::GetProcessName(), G4HadronicProcessStore::Instance(), MAX_SECONDARIES, G4HadronicProcessStore::RegisterExtraProcess(), G4VProcess::SetProcessSubType(), and G4VProcess::verboseLevel.
00047 : 00048 G4VRestProcess (processName, aType), // initialization 00049 massPionMinus(G4PionMinus::PionMinus()->GetPDGMass()/GeV), 00050 pdefGamma(G4Gamma::Gamma()), 00051 pdefPionZero(G4PionZero::PionZero()), 00052 pdefPionMinus(G4PionMinus::PionMinus()), 00053 pdefProton(G4Proton::Proton()), 00054 pdefNeutron(G4Neutron::Neutron()), 00055 pdefDeuteron(G4Deuteron::Deuteron()), 00056 pdefTriton(G4Triton::Triton()), 00057 pdefAlpha(G4Alpha::Alpha()) 00058 { 00059 G4HadronicDeprecate("G4PiMinusAbsorptionAtRest"); 00060 00061 if (verboseLevel>0) { 00062 G4cout << GetProcessName() << " is created "<< G4endl; 00063 } 00064 SetProcessSubType(fHadronAtRest); 00065 pv = new G4GHEKinematicsVector [MAX_SECONDARIES+1]; 00066 eve = new G4GHEKinematicsVector [MAX_SECONDARIES]; 00067 gkin = new G4GHEKinematicsVector [MAX_SECONDARIES]; 00068 00069 G4HadronicProcessStore::Instance()->RegisterExtraProcess(this); 00070 }
G4PionMinusAbsorptionAtRest::~G4PionMinusAbsorptionAtRest | ( | ) |
Definition at line 74 of file G4PionMinusAbsorptionAtRest.cc.
References G4HadronicProcessStore::DeRegisterExtraProcess(), and G4HadronicProcessStore::Instance().
00075 { 00076 G4HadronicProcessStore::Instance()->DeRegisterExtraProcess(this); 00077 delete [] pv; 00078 delete [] eve; 00079 delete [] gkin; 00080 }
G4VParticleChange * G4PionMinusAbsorptionAtRest::AtRestDoIt | ( | const G4Track & | , | |
const G4Step & | ||||
) | [virtual] |
Reimplemented from G4VRestProcess.
Definition at line 142 of file G4PionMinusAbsorptionAtRest.cc.
References G4ParticleChange::AddSecondary(), G4VProcess::aParticleChange, fStopAndKill, G4cout, G4endl, G4UniformRand, G4Material::GetAtomicNumDensityVector(), G4Material::GetElementVector(), G4Track::GetGlobalTime(), G4Track::GetMaterial(), G4Material::GetNumberOfElements(), G4Track::GetPosition(), G4GHEKinematicsVector::GetTOF(), G4Track::GetTouchableHandle(), G4ParticleChange::Initialize(), position, G4VParticleChange::ProposeLocalEnergyDeposit(), G4VParticleChange::ProposeTrackStatus(), G4VProcess::ResetNumberOfInteractionLengthLeft(), G4VParticleChange::SetNumberOfSecondaries(), G4Track::SetTouchableHandle(), and G4VProcess::verboseLevel.
00151 { 00152 00153 // Initialize ParticleChange 00154 // all members of G4VParticleChange are set to equal to 00155 // corresponding member in G4Track 00156 00157 aParticleChange.Initialize(track); 00158 00159 // Store some global quantities that depend on current material and particle 00160 00161 globalTime = track.GetGlobalTime()/s; 00162 G4Material * aMaterial = track.GetMaterial(); 00163 const G4int numberOfElements = aMaterial->GetNumberOfElements(); 00164 const G4ElementVector* theElementVector = aMaterial->GetElementVector(); 00165 00166 const G4double* theAtomicNumberDensity = aMaterial->GetAtomicNumDensityVector(); 00167 G4double normalization = 0; 00168 for ( G4int i1=0; i1 < numberOfElements; i1++ ) 00169 { 00170 normalization += theAtomicNumberDensity[i1] ; // change when nucleon specific 00171 // probabilities are included. 00172 } 00173 G4double runningSum= 0.; 00174 G4double random = G4UniformRand()*normalization; 00175 for ( G4int i2=0; i2 < numberOfElements; i2++ ) 00176 { 00177 runningSum += theAtomicNumberDensity[i2]; // change when nucleon specific 00178 // probabilities are included. 00179 if (random<=runningSum) 00180 { 00181 targetCharge = G4double((*theElementVector)[i2]->GetZ()); 00182 targetAtomicMass = (*theElementVector)[i2]->GetN(); 00183 } 00184 } 00185 if (random>runningSum) 00186 { 00187 targetCharge = G4double((*theElementVector)[numberOfElements-1]->GetZ()); 00188 targetAtomicMass = (*theElementVector)[numberOfElements-1]->GetN(); 00189 00190 } 00191 00192 if (verboseLevel>1) { 00193 G4cout << "G4PionMinusAbsorptionAtRest::AtRestDoIt is invoked " <<G4endl; 00194 } 00195 00196 G4ParticleMomentum momentum; 00197 G4float localtime; 00198 00199 G4ThreeVector position = track.GetPosition(); 00200 00201 GenerateSecondaries(); // Generate secondaries 00202 00203 aParticleChange.SetNumberOfSecondaries( ngkine ); 00204 00205 for ( G4int isec = 0; isec < ngkine; isec++ ) { 00206 G4DynamicParticle* aNewParticle = new G4DynamicParticle; 00207 aNewParticle->SetDefinition( gkin[isec].GetParticleDef() ); 00208 aNewParticle->SetMomentum( gkin[isec].GetMomentum() * GeV ); 00209 00210 localtime = globalTime + gkin[isec].GetTOF(); 00211 00212 G4Track* aNewTrack = new G4Track( aNewParticle, localtime*s, position ); 00213 aNewTrack->SetTouchableHandle(track.GetTouchableHandle()); 00214 aParticleChange.AddSecondary( aNewTrack ); 00215 00216 } 00217 00218 aParticleChange.ProposeLocalEnergyDeposit( 0.0*GeV ); 00219 00220 aParticleChange.ProposeTrackStatus(fStopAndKill); // Kill the incident PionMinus 00221 00222 // clear InteractionLengthLeft 00223 00224 ResetNumberOfInteractionLengthLeft(); 00225 00226 return &aParticleChange; 00227 00228 }
G4double G4PionMinusAbsorptionAtRest::AtRestGetPhysicalInteractionLength | ( | const G4Track & | , | |
G4ForceCondition * | ||||
) | [virtual] |
Reimplemented from G4VRestProcess.
Definition at line 116 of file G4PionMinusAbsorptionAtRest.cc.
References G4VProcess::currentInteractionLength, G4DynamicParticle::DumpInfo(), G4cout, G4endl, G4Track::GetDynamicParticle(), G4Track::GetMaterial(), GetMeanLifeTime(), G4Material::GetName(), G4VProcess::GetProcessName(), NotForced, ns, G4VProcess::ResetNumberOfInteractionLengthLeft(), G4VProcess::theNumberOfInteractionLengthLeft, and G4VProcess::verboseLevel.
00120 { 00121 // beggining of tracking 00122 ResetNumberOfInteractionLengthLeft(); 00123 00124 // condition is set to "Not Forced" 00125 *condition = NotForced; 00126 00127 // get mean life time 00128 currentInteractionLength = GetMeanLifeTime(track, condition); 00129 00130 if ((currentInteractionLength <0.0) || (verboseLevel>2)){ 00131 G4cout << "G4PionMinusAbsorptionAtRestProcess::AtRestGetPhysicalInteractionLength "; 00132 G4cout << "[ " << GetProcessName() << "]" <<G4endl; 00133 track.GetDynamicParticle()->DumpInfo(); 00134 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl; 00135 G4cout << "MeanLifeTime = " << currentInteractionLength/ns << "[ns]" <<G4endl; 00136 } 00137 00138 return theNumberOfInteractionLengthLeft * currentInteractionLength; 00139 00140 }
void G4PionMinusAbsorptionAtRest::BuildPhysicsTable | ( | const G4ParticleDefinition & | ) | [virtual] |
Reimplemented from G4VProcess.
Definition at line 87 of file G4PionMinusAbsorptionAtRest.cc.
References G4HadronicProcessStore::Instance(), and G4HadronicProcessStore::PrintInfo().
00088 { 00089 G4HadronicProcessStore::Instance()->PrintInfo(&p); 00090 }
G4double G4PionMinusAbsorptionAtRest::GetMeanLifeTime | ( | const G4Track & | , | |
G4ForceCondition * | ||||
) | [inline, virtual] |
Implements G4VRestProcess.
Definition at line 72 of file G4PionMinusAbsorptionAtRest.hh.
Referenced by AtRestGetPhysicalInteractionLength().
G4int G4PionMinusAbsorptionAtRest::GetNumberOfSecondaries | ( | ) |
G4GHEKinematicsVector * G4PionMinusAbsorptionAtRest::GetSecondaryKinematics | ( | ) |
G4bool G4PionMinusAbsorptionAtRest::IsApplicable | ( | const G4ParticleDefinition & | ) | [virtual] |
void G4PionMinusAbsorptionAtRest::PreparePhysicsTable | ( | const G4ParticleDefinition & | ) | [virtual] |
Reimplemented from G4VProcess.
Definition at line 82 of file G4PionMinusAbsorptionAtRest.cc.
References G4HadronicProcessStore::Instance(), and G4HadronicProcessStore::RegisterParticleForExtraProcess().
00083 { 00084 G4HadronicProcessStore::Instance()->RegisterParticleForExtraProcess(this, &p); 00085 }