Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4UCNLoss.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//
27////////////////////////////////////////////////////////////////////////
28// Ultra Cold Neutron Loss Class Implementation
29////////////////////////////////////////////////////////////////////////
30//
31// File: G4UCNLoss.cc
32// Description: Discrete Process -- Loss of UCN
33// energy-independent loss cross section inside a material
34// Version: 1.0
35// Created: 2014-05-05
36// Author: Peter Gumplinger
37// Adopted from: UCN simple absorption, Peter Fierlinger 7.9.04
38// Marcin Kuzniaka 21.4.06
39// Updated:
40// mail: gum@triumf.ca
41//
42////////////////////////////////////////////////////////////////////////
43
45
46#include "G4UCNLoss.hh"
47
48#include "G4SystemOfUnits.hh"
50
51/////////////////////////
52// Class Implementation
53/////////////////////////
54
55 //////////////
56 // Operators
57 //////////////
58
59// G4UCNLoss::operator=(const G4UCNLoss &right)
60// {
61// }
62
63 /////////////////
64 // Constructors
65 /////////////////
66
68 : G4VDiscreteProcess(processName, type)
69{
70 if (verboseLevel>0) G4cout << GetProcessName() << " is created " << G4endl;
71
73}
74
75// G4UCNLoss::G4UCNLoss(const G4UCNLoss &right)
76// {
77// }
78
79 ////////////////
80 // Destructors
81 ////////////////
82
84
85 ////////////
86 // Methods
87 ////////////
88
89// PostStepDoIt
90// -------------
91
93G4UCNLoss::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
94{
96
98
99 if (verboseLevel>0) G4cout << "\n** UCN lost! **" << G4endl;
100
101 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
102}
103
104// GetMeanFreePath
105// ---------------
106
108 G4double ,
110{
111 G4double AttenuationLength = DBL_MAX;
112
113 const G4Material* aMaterial = aTrack.GetMaterial();
114 G4MaterialPropertiesTable* aMaterialPropertiesTable =
115 aMaterial->GetMaterialPropertiesTable();
116
117 G4double crossect = 0.0;
118 if (aMaterialPropertiesTable) {
119 crossect = aMaterialPropertiesTable->GetConstProperty("LOSSCS");
120// if (crossect == 0.0)
121// G4cout << "No Loss Cross Section specified" << G4endl;
122 }
123// else G4cout << "No Loss Cross Section specified" << G4endl;
124
125 if (crossect) {
126
127 // Calculate a UCN absorption length for this cross section
128
129 G4double density = aMaterial->GetTotNbOfAtomsPerVolume();
130
131 // Calculate cross section for a constant loss
132
133 crossect *= barn;
134
135 // Attenuation length
136
137 AttenuationLength = 1./density/crossect;
138 }
139
140 return AttenuationLength;
141}
G4ForceCondition
G4ProcessType
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
@ fUCNLoss
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetConstProperty(const G4String &key) const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:251
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
void Initialize(const G4Track &) override
Definition: G4Step.hh:62
G4Material * GetMaterial() const
G4UCNLoss(const G4String &processName="UCNLoss", G4ProcessType type=fUCN)
Definition: G4UCNLoss.cc:67
virtual ~G4UCNLoss()
Definition: G4UCNLoss.cc:83
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4UCNLoss.cc:93
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *condition)
Definition: G4UCNLoss.cc:107
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void ProposeTrackStatus(G4TrackStatus status)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:331
G4int verboseLevel
Definition: G4VProcess.hh:360
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
#define DBL_MAX
Definition: templates.hh:62