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
G4UCNMultiScattering.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// UCN Multiple Scattering Class Implementation
29///////////////////////////////////////////////////////////////////////
30//
31// File: G4UCNMultiScattering.cc
32// Description: G4VDiscreteProcess -- MultiScattering of UCNs
33// Version: 1.0
34// Created: 2014-05-12
35// Author: Peter Gumplinger
36// adopted from Geant4UCN by Peter Fierlinger (7.9.04) and
37// Marcin Kuzniak (21.4.06)
38// Calculate "elastic scattering" inside materials
39// Updated:
40//
41// mail: gum@triumf.ca
42//
43///////////////////////////////////////////////////////////////////////
44
46
48
49#include "G4SystemOfUnits.hh"
51
52/////////////////////////
53// Class Implementation
54/////////////////////////
55
56 //////////////
57 // Operators
58 //////////////
59
60 // G4UCNMultiScattering::operator=(const G4UCNMultiScattering &right)
61 // {
62 // }
63
64 /////////////////
65 // Constructors
66 /////////////////
67
69 G4ProcessType type)
70 : G4VDiscreteProcess(processName, type)
71{
72 if (verboseLevel>0) G4cout << GetProcessName() << " is created " << G4endl;
73
75}
76
77// G4UCNMultiScattering::G4UCNMultiScattering(const G4UCNMultiScattering &right)
78// {
79// }
80
81 ////////////////
82 // Destructors
83 ////////////////
84
86
87 ////////////
88 // Methods
89 ////////////
90
91// PostStepDoIt
92// ------------
93
96{
98
99 if ( verboseLevel > 0 ) G4cout << "UCNMULTISCATTER at: "
100 << aTrack.GetProperTime()/s << "s, "
101 << aTrack.GetGlobalTime()/s << "s. "
102 << ", after track length " << aTrack.GetTrackLength()/cm << "cm, "
103 << "in volume "
105 << G4endl;
106
107 G4ThreeVector scattered = Scatter();
108
110
111 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
112}
113
114// GetMeanFreePath
115// ---------------
116
118 G4double ,
120{
121 G4double AttenuationLength = DBL_MAX;
122
123 const G4Material* aMaterial = aTrack.GetMaterial();
124 G4MaterialPropertiesTable* aMaterialPropertiesTable =
125 aMaterial->GetMaterialPropertiesTable();
126
127 G4double crossect = 0.0;
128 if (aMaterialPropertiesTable) {
129 crossect = aMaterialPropertiesTable->GetConstProperty("SCATCS");
130// if(crossect == 0.0)
131// G4cout << "No UCN MultiScattering length specified" << G4endl;
132 }
133// else G4cout << "No UCN MultiScattering length specified" << G4endl;
134
135 if (crossect) {
136
137 // Calculate a UCN MultiScattering length for this cross section
138
139 G4double density = aMaterial->GetTotNbOfAtomsPerVolume();
140
141 crossect *= barn;
142
143 // attenuation length in mm
144 AttenuationLength = 1./density/crossect;
145 }
146
147 return AttenuationLength;
148}
149
151{
152
153 G4ThreeVector final(0.,0.,1.);
154
155 // Make a simple uniform distribution in 4 pi
156 // apply scattering, calculate angle phi, theta
157
158 G4double theta = std::acos(2*G4UniformRand()-1);
159 G4double phi = G4UniformRand() * 2 * pi;
160
161 final.rotateY(theta);
162 final.rotateZ(phi);
163 final = final.unit();
164
165 return final;
166}
G4ForceCondition
G4ProcessType
double G4double
Definition: G4Types.hh:83
@ fUCNMultiScattering
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
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
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4VPhysicalVolume * GetPhysicalVolume() const
Definition: G4Step.hh:62
G4StepPoint * GetPostStepPoint() const
G4double GetTrackLength() const
G4double GetGlobalTime() const
G4double GetProperTime() const
G4Material * GetMaterial() const
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
G4UCNMultiScattering(const G4String &processName="UCNMultiScattering", G4ProcessType type=fUCN)
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *condition)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
const G4String & GetName() const
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