Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4UCNMultiScattering.cc
이 파일의 문서화 페이지로 가기
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: G4UCNMultiScattering.cc 69576 2013-05-08 13:48:13Z gcosmo $
27 //
29 // UCN Multiple Scattering Class Implementation
31 //
32 // File: G4UCNMultiScattering.cc
33 // Description: G4VDiscreteProcess -- MultiScattering of UCNs
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 // Calculate "elastic scattering" inside materials
40 // Updated:
41 //
42 // mail: gum@triumf.ca
43 //
45 
46 #include "G4UCNProcessSubType.hh"
47 
48 #include "G4UCNMultiScattering.hh"
49 
50 #include "G4SystemOfUnits.hh"
51 #include "G4PhysicalConstants.hh"
52 
54 // Class Implementation
56 
58  // Operators
60 
61  // G4UCNMultiScattering::operator=(const G4UCNMultiScattering &right)
62  // {
63  // }
64 
66  // Constructors
68 
70  G4ProcessType type)
71  : G4VDiscreteProcess(processName, type)
72 {
73  if (verboseLevel>0) G4cout << GetProcessName() << " is created " << G4endl;
74 
76 }
77 
78 // G4UCNMultiScattering::G4UCNMultiScattering(const G4UCNMultiScattering &right)
79 // {
80 // }
81 
83  // Destructors
85 
87 
89  // Methods
91 
92 // PostStepDoIt
93 // ------------
94 
97 {
99 
100  if ( verboseLevel > 0 ) G4cout << "UCNMULTISCATTER at: "
101  << aTrack.GetProperTime()/s << "s, "
102  << aTrack.GetGlobalTime()/s << "s. "
103  << ", after track length " << aTrack.GetTrackLength()/cm << "cm, "
104  << "in volume "
106  << G4endl;
107 
108  G4ThreeVector scattered = Scatter();
109 
111 
112  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
113 }
114 
115 // GetMeanFreePath
116 // ---------------
117 
119  G4double ,
121 {
122  G4double AttenuationLength = DBL_MAX;
123 
124  const G4Material* aMaterial = aTrack.GetMaterial();
125  G4MaterialPropertiesTable* aMaterialPropertiesTable =
126  aMaterial->GetMaterialPropertiesTable();
127 
128  G4double crossect = 0.0;
129  if (aMaterialPropertiesTable) {
130  crossect = aMaterialPropertiesTable->GetConstProperty("SCATCS");
131 // if(crossect == 0.0)
132 // G4cout << "No UCN MultiScattering length specified" << G4endl;
133  }
134 // else G4cout << "No UCN MultiScattering length specified" << G4endl;
135 
136  if (crossect) {
137 
138  // Calculate a UCN MultiScattering length for this cross section
139 
140  G4double density = aMaterial->GetTotNbOfAtomsPerVolume();
141 
142  crossect *= barn;
143 
144  // attenuation length in mm
145  AttenuationLength = 1./density/crossect;
146  }
147 
148  return AttenuationLength;
149 }
150 
152 {
153 
154  G4ThreeVector final(0.,0.,1.);
155 
156  // Make a simple uniform distribution in 4 pi
157  // apply scattering, calculate angle phi, theta
158 
159  G4double theta = std::acos(2*G4UniformRand()-1);
160  G4double phi = G4UniformRand() * 2 * pi;
161 
162  final.rotateY(theta);
163  final.rotateZ(phi);
164  final = final.unit();
165 
166  return final;
167 }
G4double GetConstProperty(const char *key) const
#define G4endl
Definition: G4ios.hh:61
G4UCNMultiScattering(const G4String &processName="UCNMultiScattering", G4ProcessType type=fUCN)
G4VPhysicalVolume * GetPhysicalVolume() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4double GetProperTime() const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
G4double GetTrackLength() const
const XML_Char * s
Definition: expat.h:262
double G4double
Definition: G4Types.hh:76
virtual void Initialize(const G4Track &)
G4ProcessType
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *condition)
#define G4UniformRand()
Definition: Randomize.hh:53
G4double GetGlobalTime() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:411
Definition: G4Step.hh:76
G4StepPoint * GetPostStepPoint() const
G4Material * GetMaterial() const
G4int verboseLevel
Definition: G4VProcess.hh:371
static constexpr double barn
Definition: G4SIunits.hh:105
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:435
G4ForceCondition
static constexpr double cm
Definition: G4SIunits.hh:119
G4GLOB_DLL std::ostream G4cout
static constexpr double pi
Definition: G4SIunits.hh:75
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:210
#define DBL_MAX
Definition: templates.hh:83
const G4String & GetName() const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252