Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4KaonMinusAbsorptionAtRest Class Reference

#include <G4KaonMinusAbsorptionAtRest.hh>

+ Inheritance diagram for G4KaonMinusAbsorptionAtRest:

Public Member Functions

 G4KaonMinusAbsorptionAtRest (const G4String &processName="KaonMinusAbsorptionAtRest", G4ProcessType aType=fHadronic)
 
 ~G4KaonMinusAbsorptionAtRest ()
 
G4bool IsApplicable (const G4ParticleDefinition &particle)
 
void PreparePhysicsTable (const G4ParticleDefinition &)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4VParticleChangeAtRestDoIt (const G4Track &aTrack, const G4Step &aStep)
 
- Public Member Functions inherited from G4VRestProcess
 G4VRestProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VRestProcess (G4VRestProcess &)
 
virtual ~G4VRestProcess ()
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double, G4ForceCondition *)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

Protected Member Functions

G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *)
 
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 54 of file G4KaonMinusAbsorptionAtRest.hh.

Constructor & Destructor Documentation

◆ G4KaonMinusAbsorptionAtRest()

G4KaonMinusAbsorptionAtRest::G4KaonMinusAbsorptionAtRest ( const G4String processName = "KaonMinusAbsorptionAtRest",
G4ProcessType  aType = fHadronic 
)

Definition at line 48 of file G4KaonMinusAbsorptionAtRest.cc.

49 :
50 G4VRestProcess (processName, aType)
51{
52 G4HadronicDeprecate("G4KaonMinusAbsorptionAtRest");
53 if (verboseLevel>0) {
54 G4cout << GetProcessName() << " is created "<< G4endl;
55 }
57
58 // see Cohn et al, PLB27(1968) 527;
59 // Davis et al, PLB1(1967) 434;
60
61 pionAbsorptionRate = 0.07;
62
63 // see VanderVelde-Wilquet et al, Nuov.Cim.39A(1978)538;
64 // see VanderVelde-Wilquet et al, Nuov.Cim.38A(1977)178;
65 // see VanderVelde-Wilquet et al, Nucl.Phys.A241(1975)511;
66 // primary production rates ( for absorption on Carbon)
67 // .. other elements are extrapolated by the halo factor.
68
69 rateLambdaZeroPiZero = 0.052;
70 rateSigmaMinusPiPlus = 0.199;
71 rateSigmaPlusPiMinus = 0.446;
72 rateSigmaZeroPiZero = 0.303;
73 rateLambdaZeroPiMinus = 0.568;
74 rateSigmaZeroPiMinus = 0.216;
75 rateSigmaMinusPiZero = 0.216;
76
77 // for sigma- p -> lambda n
78 // sigma+ n -> lambda p
79 // sigma- n -> lambda
80 // all values compatible with 0.55 same literature as above.
81
82 sigmaPlusLambdaConversionRate = 0.55;
83 sigmaMinusLambdaConversionRate = 0.55;
84 sigmaZeroLambdaConversionRate = 0.55;
85
87}
#define G4HadronicDeprecate(name)
@ fHadronAtRest
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void RegisterExtraProcess(G4VProcess *)
static G4HadronicProcessStore * Instance()
G4int verboseLevel
Definition: G4VProcess.hh:368
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379

◆ ~G4KaonMinusAbsorptionAtRest()

G4KaonMinusAbsorptionAtRest::~G4KaonMinusAbsorptionAtRest ( )

Definition at line 90 of file G4KaonMinusAbsorptionAtRest.cc.

Member Function Documentation

◆ AtRestDoIt()

G4VParticleChange * G4KaonMinusAbsorptionAtRest::AtRestDoIt ( const G4Track aTrack,
const G4Step aStep 
)
virtual

Reimplemented from G4VRestProcess.

Definition at line 105 of file G4KaonMinusAbsorptionAtRest.cc.

107{
108 stoppedHadron = track.GetDynamicParticle();
109
110 // Check applicability
111
112 if (!IsApplicable(*(stoppedHadron->GetDefinition())))
113 {
114 G4cerr <<"G4KaonMinusAbsorptionAtRest:ERROR, particle must be a Kaon!" <<G4endl;
115 return 0;
116 }
117
118 G4Material* material;
119 material = track.GetMaterial();
120 nucleus = 0;
121 do
122 {
123 // Select the nucleus, get nucleon
124 nucleus = new G4Nucleus(material);
125 if (nucleus->GetA_asInt() < 1.5)
126 {
127 delete nucleus;
128 nucleus = 0;
129 }
130 } while(nucleus == 0);
131
132 G4double Z = nucleus->GetZ_asInt();
133 G4double A = nucleus->GetA_asInt();
134
135 // Do the interaction with the nucleon
136 G4DynamicParticleVector* absorptionProducts = KaonNucleonReaction();
137
138 //A.R. 26-Jul-2012 Coverity fix
139 if ( ! absorptionProducts ) {
140 G4Exception("G4KaonMinusAbsorptionAtRest::AtRestDoIt()", "HAD_STOP_0001",
141 FatalException, "NULL absorptionProducts");
142 return 0;
143 }
144
145 // Secondary interactions
146
147 G4DynamicParticle* thePion;
148 unsigned int i;
149 for(i = 0; i < absorptionProducts->size(); i++)
150 {
151 thePion = (*absorptionProducts)[i];
152 if (thePion->GetDefinition() == G4PionMinus::PionMinus()
153 || thePion->GetDefinition() == G4PionPlus::PionPlus()
154 || thePion->GetDefinition() == G4PionZero::PionZero())
155 {
156 if (AbsorbPionByNucleus(thePion))
157 {
158 absorptionProducts->erase(absorptionProducts->begin()+i);
159 i--;
160 delete thePion;
161 if (verboseLevel > 1)
162 G4cout << "G4KaonMinusAbsorption::AtRestDoIt: Pion absorbed in Nucleus"
163 << G4endl;
164 }
165 }
166 }
167
168 G4DynamicParticle* theSigma;
169 G4DynamicParticle* theLambda;
170 for (i = 0; i < absorptionProducts->size(); i++)
171 {
172 theSigma = (*absorptionProducts)[i];
173 if (theSigma->GetDefinition() == G4SigmaMinus::SigmaMinus()
174 || theSigma->GetDefinition() == G4SigmaPlus::SigmaPlus()
175 || theSigma->GetDefinition() == G4SigmaZero::SigmaZero())
176 {
177 theLambda = SigmaLambdaConversion(theSigma);
178 if (theLambda != 0){
179 absorptionProducts->erase(absorptionProducts->begin()+i);
180 i--;
181 delete theSigma;
182 absorptionProducts->push_back(theLambda);
183
184 if (verboseLevel > 1)
185 G4cout << "G4KaonMinusAbsorption::AtRestDoIt: SigmaLambdaConversion Done"
186 << G4endl;
187 }
188 }
189 }
190
191 // Nucleus deexcitation
192
193 G4double productEnergy = 0.;
194 G4ThreeVector pProducts(0.,0.,0.);
195
196 unsigned int nAbsorptionProducts = 0;
197 if (absorptionProducts != 0) nAbsorptionProducts = absorptionProducts->size();
198
199 for ( i = 0; i<nAbsorptionProducts; i++)
200 {
201 pProducts += (*absorptionProducts)[i]->GetMomentum();
202 productEnergy += (*absorptionProducts)[i]->GetKineticEnergy();
203 }
204
205 G4double newZ = nucleus->GetZ_asInt();
206 G4double newA = nucleus->GetA_asInt();
207
208 G4double bDiff = G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(A),static_cast<G4int>(Z)) -
209 G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(newA), static_cast<G4int>(newZ));
210
211 G4StopDeexcitationAlgorithm* nucleusAlgorithm = new G4StopTheoDeexcitation();
212 G4StopDeexcitation stopDeexcitation(nucleusAlgorithm);
213
214 nucleus->AddExcitationEnergy(bDiff);
215
216 // returns excitation energy for the moment ..
217 G4double energyDeposit = nucleus->GetEnergyDeposit();
218 if (verboseLevel>0)
219 {
220 G4cout << " -- KaonAtRest -- excitation = "
221 << energyDeposit
222 << ", pNucleus = "
223 << pProducts
224 << ", A: "
225 << A
226 << ", "
227 << newA
228 << ", Z: "
229 << Z
230 << ", "
231 << newZ
232 << G4endl;
233 }
234
235 if (energyDeposit < 0.)
236 G4Exception("G4KaonMinusAbsorptionAtRest::AtRestDoIt()", "HAD_STOP_0002",
237 FatalException, "Excitation energy < 0");
238 delete nucleus;
239
240 G4ReactionProductVector* fragmentationProducts = stopDeexcitation.DoBreakUp(newA,newZ,energyDeposit,pProducts);
241
242 unsigned int nFragmentationProducts = 0;
243 if (fragmentationProducts != 0) nFragmentationProducts = fragmentationProducts->size();
244
245 //Initialize ParticleChange
247 aParticleChange.SetNumberOfSecondaries(G4int(nAbsorptionProducts+nFragmentationProducts) );
248
249 // update List of alive particles. put energy deposit at the right place ...
250 for (i = 0; i < nAbsorptionProducts; i++)
251 {aParticleChange.AddSecondary((*absorptionProducts)[i]); }
252 if (absorptionProducts != 0) delete absorptionProducts;
253
254// for (i = 0; i < nFragmentationProducts; i++)
255// { aParticleChange.AddSecondary(fragmentationProducts->at(i)); }
256 for(i=0; i<nFragmentationProducts; i++)
257 {
258 G4DynamicParticle * aNew =
259 new G4DynamicParticle((*fragmentationProducts)[i]->GetDefinition(),
260 (*fragmentationProducts)[i]->GetTotalEnergy(),
261 (*fragmentationProducts)[i]->GetMomentum());
262 G4double newTime = aParticleChange.GetGlobalTime((*fragmentationProducts)[i]->GetFormationTime());
263 aParticleChange.AddSecondary(aNew, newTime);
264 delete (*fragmentationProducts)[i];
265 }
266 if (fragmentationProducts != 0) delete fragmentationProducts;
267
268 // finally ...
269 aParticleChange.ProposeTrackStatus(fStopAndKill); // Kill the incident Kaon
270 return &aParticleChange;
271}
std::vector< G4DynamicParticle * > G4DynamicParticleVector
@ FatalException
std::vector< G4ReactionProduct * > G4ReactionProductVector
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4DLLIMPORT std::ostream G4cerr
G4ParticleDefinition * GetDefinition() const
G4bool IsApplicable(const G4ParticleDefinition &particle)
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:576
static G4double GetBindingEnergy(const G4int A, const G4int Z)
void AddExcitationEnergy(G4double anEnergy)
Definition: G4Nucleus.cc:435
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4double GetEnergyDeposit()
Definition: G4Nucleus.hh:184
void AddSecondary(G4Track *aSecondary)
virtual void Initialize(const G4Track &)
G4double GetGlobalTime(G4double timeDelay=0.0) const
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static G4PionZero * PionZero()
Definition: G4PionZero.cc:104
static G4SigmaMinus * SigmaMinus()
static G4SigmaPlus * SigmaPlus()
Definition: G4SigmaPlus.cc:108
static G4SigmaZero * SigmaZero()
Definition: G4SigmaZero.cc:99
void ProposeTrackStatus(G4TrackStatus status)
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

◆ BuildPhysicsTable()

void G4KaonMinusAbsorptionAtRest::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Definition at line 100 of file G4KaonMinusAbsorptionAtRest.cc.

101{
103}
void PrintInfo(const G4ParticleDefinition *)

◆ GetMeanLifeTime()

G4double G4KaonMinusAbsorptionAtRest::GetMeanLifeTime ( const G4Track aTrack,
G4ForceCondition  
)
inlineprotectedvirtual

Implements G4VRestProcess.

Definition at line 81 of file G4KaonMinusAbsorptionAtRest.hh.

83 {
84 G4double result = 0;
85 if(aTrack.GetMaterial()->GetNumberOfElements() == 1)
86 if(aTrack.GetMaterial()->GetZ()<1.5) result = DBL_MAX;
87 return result;
88 }
G4double GetZ() const
Definition: G4Material.cc:604
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
G4Material * GetMaterial() const
#define DBL_MAX
Definition: templates.hh:83

◆ IsApplicable()

G4bool G4KaonMinusAbsorptionAtRest::IsApplicable ( const G4ParticleDefinition particle)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 68 of file G4KaonMinusAbsorptionAtRest.hh.

68 {
69 return( particle == *(G4KaonMinus::KaonMinus()) );
70 }
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113

Referenced by AtRestDoIt().

◆ PreparePhysicsTable()

void G4KaonMinusAbsorptionAtRest::PreparePhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Definition at line 95 of file G4KaonMinusAbsorptionAtRest.cc.

96{
98}
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)

The documentation for this class was generated from the following files: