Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4UCNAbsorption.cc
Go to the documentation of this file.
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 // $Id: G4UCNAbsorption.cc 69576 2013-05-08 13:48:13Z gcosmo $
27 //
29 // UCN Absorption Class Implementation
31 //
32 // File: G4UCNAbsorption.cc
33 // Description: Discrete Process -- Absorption of Ultra Cold Neutrons
34 // Version: 1.0
35 // Created: 2014-05-12
36 // Author: Peter Gumplinger
37 // adopted from Geant4UCN by Peter Fierlinger (7.9.04) and
38 // Marcin Kuzniak (21.4.06)
39 // 1/v energy dependent absorption cross section
40 // inside materials
41 // Updated:
42 //
43 // mail: gum@triumf.ca
44 //
46 
47 #include "G4UCNProcessSubType.hh"
48 
49 #include "G4UCNAbsorption.hh"
50 
51 //#include "G4Nucleus.hh"
52 //#include "G4ReactionProduct.hh"
53 //#include "G4NucleiPropertiesTable.hh"
54 
55 #include "G4SystemOfUnits.hh"
56 #include "G4PhysicalConstants.hh"
57 
59 // Class Implementation
61 
63  // Operators
65 
66 // G4UCNAbsorption::operator=(const G4UCNAbsorption &right)
67 // {
68 // }
69 
71  // Constructors
73 
75  : G4VDiscreteProcess(processName, type)
76 {
77  if (verboseLevel>0) G4cout << GetProcessName() << " is created " << G4endl;
78 
80 }
81 
82 // G4UCNAbsorption::G4UCNAbsorption(const G4UCNAbsorpton &right)
83 // {
84 // }
85 
87  // Destructors
89 
91 
93  // Methods
95 
96 // PostStepDoIt
97 // -------------
98 
100 G4UCNAbsorption::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
101 {
102  aParticleChange.Initialize(aTrack);
103 
105 
106  if ( verboseLevel > 0 ) G4cout << "UCNABSORPTION at: "
107  << aTrack.GetProperTime()/s << "s, "
108  << aTrack.GetGlobalTime()/s << "s. "
109  << ", after track length " << aTrack.GetTrackLength()/cm << "cm, "
110  << "in volume "
112  << G4endl;
113 
114  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
115 }
116 
117 // GetMeanFreePath
118 // ---------------
119 
121  G4double ,
123 {
124  G4double AttenuationLength = DBL_MAX;
125 
126  const G4Material* aMaterial = aTrack.GetMaterial();
127  G4MaterialPropertiesTable* aMaterialPropertiesTable =
128  aMaterial->GetMaterialPropertiesTable();
129 
130  G4double losscs = 0.0;
131  if (aMaterialPropertiesTable) {
132  losscs = aMaterialPropertiesTable->GetConstProperty("ABSCS");
133 // if (losscs == 0.0)
134 // G4cout << "No UCN Absorption length specified" << G4endl;
135  }
136 // else G4cout << "No UCN Absorption length specified" << G4endl;
137 
138  if (losscs) {
139 
140  // Calculate a UCN absorption length for this cross section
141 
142  // *** Thermal boost ***
143 
144  // Prepare neutron
145 
146  //G4double theA = aMaterial->GetElement(0)->GetN();
147  //G4double theZ = aMaterial->GetElement(0)->GetZ();
148 
149  //G4ReactionProduct
150  // theNeutron(const_cast<G4ParticleDefinition *>(aTrack.GetDefinition()));
151  //theNeutron.SetMomentum(aTrack.GetMomentum());
152  //theNeutron.SetKineticEnergy(aTrack.GetKineticEnergy());
153  //G4ThreeVector neuVelo = theNeutron.GetMomentum()/
154  // aTrack.GetDefinition()->GetPDGMass());
155 
156  // Prepare properly biased thermal nucleus
157 
158  //G4double theA = aMaterial->GetElement(0)->GetN();
159  //G4double theZ = aMaterial->GetElement(0)->GetZ();
160 
161  //G4double eps = 0.0001;
162 
163  //G4double eleMass =
164  // G4NucleiPropertiesTable::
165  // GetNuclearMass(static_cast<G4int>(theZ+eps),
166  // static_cast<G4int>(theA+eps)))
167  // / G4Neutron::Neutron()->GetPDGMass();
168 
169  //G4Nucleus aNuc;
170 
171  //G4ReactionProduct aThermalNuc =
172  // aNuc.GetBiasedThermalNucleus(eleMass,
173  // neuVelo,
174  // aMaterial->GetTemperature());
175 
176  // Boost to rest system and return
177 
178  //G4ReactionProduct boosted;
179  //boosted.Lorentz(theNeutron, aThermalNuc);
180 
181  //G4double vel = sqrt(2*boosted.GetKineticEnergy()/
182  // neutron_mass_c2*c_squared);
183 
184  G4double density = aMaterial->GetTotNbOfAtomsPerVolume();
185 
186  // Calculate cross section for a constant loss
187 
188  G4double vel = aTrack.GetVelocity();
189 
190  //G4cout << aTrack.GetVelocity()/meter*second << " "
191  // << vel/meter*second << "meters/second" << G4endl;
192 
193  // Input data is normally taken from the website:
194  // http://rrdjazz.nist.gov/resources/n-lengths/list.html
195  // and coresponds to 2200 m/s fast neutrons
196 
197  G4double crossect = losscs*barn*2200.*meter/second/vel;
198 
199  // In principle, if one asks for the MaterialProperty incoherent cross
200  // section, one could put the formula for inelastic up scattering here
201  // and add the cross section to the absorption
202 
203  // sigma inelastic = ... ignatovic, p. 174.
204 
205  // attenuation length in mm
206  AttenuationLength = 1./density/crossect;
207 
208  if (verboseLevel>0) G4cout << "UCNABSORPTION with"
209  << " AttenuationLength: " << AttenuationLength/m << "m"
210  << " CrossSection: " << crossect/barn << "barn" << G4endl;
211  }
212 
213  return AttenuationLength;
214 }
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
G4int verboseLevel
Definition: G4VProcess.hh:368
virtual ~G4UCNAbsorption()
G4double GetProperTime() const
G4double GetVelocity() const
static constexpr double second
Definition: G4SIunits.hh:157
static constexpr double meter
Definition: G4SIunits.hh:82
const XML_Char * s
Definition: expat.h:262
G4GLOB_DLL std::ostream G4cout
static constexpr double m
Definition: G4SIunits.hh:129
G4VPhysicalVolume * GetPhysicalVolume() const
const G4String & GetName() const
static constexpr double cm
Definition: G4SIunits.hh:119
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
Definition: G4Step.hh:76
G4double GetGlobalTime() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4double GetTrackLength() const
G4Material * GetMaterial() const
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *condition)
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:209
virtual void Initialize(const G4Track &)
G4UCNAbsorption(const G4String &processName="UCNAbsorption", G4ProcessType type=fUCN)
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
G4StepPoint * GetPostStepPoint() const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4ForceCondition
static constexpr double barn
Definition: G4SIunits.hh:105
#define DBL_MAX
Definition: templates.hh:83
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4double GetConstProperty(const char *key) const
G4ProcessType