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