source: trunk/source/processes/hadronic/stopping/include/G4KaonMinusAbsorptionAtRest.hh@ 1036

Last change on this file since 1036 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

File size: 5.2 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// File name: G4KaonMinusAbsorptionAtRest.hh
27//
28// Author: Christian V"olcker (Christian.Volcker@cern.ch),
29//
30// Creation date: 10. November 1997
31//
32// -------------------------------------------------------------------
33
34#ifndef G4KaonMinusAbsorptionAtRest_h
35#define G4KaonMinusAbsorptionAtRest_h 1
36
37// Class Description:
38//
39// Process for nuclear absorption of K- at rest.
40// To be used in your physics list in case you need this physics.
41
42#include "globals.hh"
43#include "Randomize.hh"
44#include "G4VRestProcess.hh"
45#include "G4ParticleTypes.hh"
46#include "G4Nucleus.hh"
47#include "G4DynamicParticle.hh"
48#include "G4DynamicParticleVector.hh"
49#include "G4NucleiProperties.hh"
50#include "G4HadronicProcessType.hh"
51
52
53// *********************************************************
54class G4KaonMinusAbsorptionAtRest : public G4VRestProcess
55// *********************************************************
56{
57private:
58 // hide assignment operator as private
59 G4KaonMinusAbsorptionAtRest& operator=(const G4KaonMinusAbsorptionAtRest &right);
60 G4KaonMinusAbsorptionAtRest(const G4KaonMinusAbsorptionAtRest& );
61public:
62 G4KaonMinusAbsorptionAtRest(const G4String& processName ="KaonMinusAbsorptionAtRest",
63 G4ProcessType aType = fHadronic );
64 ~G4KaonMinusAbsorptionAtRest();
65
66 //override methods...
67public:
68 G4bool IsApplicable(const G4ParticleDefinition& particle) {
69 return( particle == *(G4KaonMinus::KaonMinus()) );
70 }
71
72 // the main method ...
73 G4VParticleChange* AtRestDoIt(const G4Track& aTrack, const G4Step& aStep);
74
75protected: // why?? might be private....
76 // zero mean lifetime
77 G4double GetMeanLifeTime(const G4Track& aTrack,
78 G4ForceCondition* )
79 {
80 G4double result = 0;
81 if(aTrack.GetMaterial()->GetNumberOfElements() == 1)
82 if(aTrack.GetMaterial()->GetZ()<1.5) result = DBL_MAX;
83 return result;
84 }
85
86private:
87 // returns proton or neutron with fermi-momentum
88 G4DynamicParticle GetAbsorbingNucleon();
89
90 // returns proton or neutron particle definition;
91 G4ParticleDefinition* SelectAbsorbingNucleon();
92
93 // provides the neutron halo factor for absorption on nucleus surface.
94 // in the G4Nucleus
95 G4double NeutronHaloFactor(G4double Z, G4double N);
96
97 // creates the reaction products
98 G4DynamicParticleVector* KaonNucleonReaction();
99
100 // secondary pion absorption in parent nucleus
101 // if TRUE, then add excitation energy to the Nucleus
102 G4bool AbsorbPionByNucleus(G4DynamicParticle* aPion);
103
104 // secondary Sigma-Lambda conversion
105 // if conversion Done, then add excitation energy to the Nucleus
106 G4DynamicParticle *SigmaLambdaConversion(G4DynamicParticle* aSigma);
107
108 // instance variables ...
109private:
110 // pointer to current stopped hadron
111 const G4DynamicParticle *stoppedHadron;
112
113 // pointer to current target nucleus
114 G4Nucleus* nucleus;
115
116 // some constant parameters
117
118 G4double pionAbsorptionRate;
119
120 // primary production rates ( for absorption on Carbon)
121
122 G4double rateLambdaZeroPiZero;
123 G4double rateSigmaMinusPiPlus;
124 G4double rateSigmaPlusPiMinus;
125 G4double rateSigmaZeroPiZero;
126
127 G4double rateLambdaZeroPiMinus;
128 G4double rateSigmaZeroPiMinus;
129 G4double rateSigmaMinusPiZero;
130
131
132 // Sigma Lambda Conversion rates
133 // for sigma- p -> lambda n
134 // sigma+ n -> lambda p
135 // sigma- n -> lambda
136
137 G4double sigmaPlusLambdaConversionRate;
138 G4double sigmaMinusLambdaConversionRate;
139 G4double sigmaZeroLambdaConversionRate;
140
141};
142
143#endif
144
Note: See TracBrowser for help on using the repository browser.