Geant4  10.03
G4PhononDownconversion.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 //
28 //
29 // $Id: G4PhononDownconversion.cc 79778 2014-03-13 17:24:53Z gcosmo $
30 //
31 // 20131111 Add verbose output for MFP calculation
32 // 20131115 Initialize data buffers in ctor
33 
35 #include "G4LatticePhysical.hh"
36 #include "G4PhononLong.hh"
37 #include "G4PhononPolarization.hh"
38 #include "G4PhononTrackMap.hh"
39 #include "G4PhononTransFast.hh"
40 #include "G4PhononTransSlow.hh"
41 #include "G4PhysicalConstants.hh"
42 #include "G4RandomDirection.hh"
43 #include "G4Step.hh"
44 #include "G4SystemOfUnits.hh"
45 #include "G4VParticleChange.hh"
46 #include "Randomize.hh"
47 #include <cmath>
48 
49 
50 
52  : G4VPhononProcess(aName), fBeta(0.), fGamma(0.), fLambda(0.), fMu(0.) {;}
53 
55 
56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
57 
59  G4double /*previousStepSize*/,
61  //Determines mean free path for longitudinal phonons to split
63  G4double Eoverh = aTrack.GetKineticEnergy()/h_Planck;
64 
65  //Calculate mean free path for anh. decay
66  G4double mfp = aTrack.GetVelocity()/(Eoverh*Eoverh*Eoverh*Eoverh*Eoverh*A);
67 
68  if (verboseLevel > 1)
69  G4cout << "G4PhononDownconversion::GetMeanFreePath = " << mfp << G4endl;
70 
71  *condition = NotForced;
72  return mfp;
73 }
74 
75 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
76 
77 
79  const G4Step&) {
81 
82  //Obtain dynamical constants from this volume's lattice
86  fMu=theLattice->GetMu();
87 
88  //Destroy the parent phonon and create the daughter phonons.
89  //74% chance that daughter phonons are both transverse
90  //26% Transverse and Longitudinal
91  if (G4UniformRand()>0.740) MakeLTSecondaries(aTrack);
92  else MakeTTSecondaries(aTrack);
93 
96 
97  return &aParticleChange;
98 }
99 
100 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
101 
103  //Only L-phonons decay
104  return (&aPD==G4PhononLong::PhononDefinition());
105 }
106 
107 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
108 
109 //probability density of energy distribution of L'-phonon in L->L'+T process
110 
112  //d=delta= ratio of group velocities vl/vt and x is the fraction of energy in the longitudinal mode, i.e. x=EL'/EL
113  return (1/(x*x))*(1-x*x)*(1-x*x)*((1+x)*(1+x)-d*d*((1-x)*(1-x)))*(1+x*x-d*d*(1-x)*(1-x))*(1+x*x-d*d*(1-x)*(1-x));
114 }
115 
116 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
117 
118 //probability density of energy distribution of T-phonon in L->T+T process
120  //dynamic constants from Tamura, PRL31, 1985
121  G4double A = 0.5*(1-d*d)*(fBeta+fLambda+(1+d*d)*(fGamma+fMu));
122  G4double B = fBeta+fLambda+2*d*d*(fGamma+fMu);
123  G4double C = fBeta + fLambda + 2*(fGamma+fMu);
124  G4double D = (1-d*d)*(2*fBeta+4*fGamma+fLambda+3*fMu);
125 
126  return (A+B*d*x-B*x*x)*(A+B*d*x-B*x*x)+(C*x*(d-x)-D/(d-x)*(x-d-(1-d*d)/(4*x)))*(C*x*(d-x)-D/(d-x)*(x-d-(1-d*d)/(4*x)));
127 }
128 
129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
130 
131 
133  //change in L'-phonon propagation direction after decay
134 
135  return std::acos((1+(x*x)-((d*d)*(1-x)*(1-x)))/(2*x));
136 }
137 
138 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
139 
140 
142  //change in T-phonon propagation direction after decay (L->L+T process)
143 
144  return std::acos((1-x*x+d*d*(1-x)*(1-x))/(2*d*(1-x)));
145 }
146 
147 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
148 
149 
151  //change in T-phonon propagation direction after decay (L->T+T process)
152 
153  return std::acos((1-d*d*(1-x)*(1-x)+d*d*x*x)/(2*d*x));
154 }
155 
156 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
157 
158 
159 //Generate daughter phonons from L->T+T process
160 
162  //d is the velocity ratio vL/vT
163  G4double d=1.6338;
164  G4double upperBound=(1+(1/d))/2;
165  G4double lowerBound=(1-(1/d))/2;
166 
167  //Use MC method to generate point from distribution:
168  //if a random point on the energy-probability plane is
169  //smaller that the curve of the probability density,
170  //then accept that point.
171  //x=fraction of parent phonon energy in first T phonon
172  G4double x = G4UniformRand()*(upperBound-lowerBound) + lowerBound;
173  G4double p = 1.5*G4UniformRand();
174  while(p >= GetTTDecayProb(d, x*d)) {
175  x = G4UniformRand()*(upperBound-lowerBound) + lowerBound;
176  p = 1.5*G4UniformRand();
177  }
178 
179  //using energy fraction x to calculate daughter phonon directions
180  G4double theta1=MakeTTDeviation(d, x);
181  G4double theta2=MakeTTDeviation(d, 1-x);
182  G4ThreeVector dir1=trackKmap->GetK(aTrack);
183  G4ThreeVector dir2=dir1;
184 
185  // FIXME: These extra randoms change timing and causting outputs of example!
186  G4ThreeVector ran = G4RandomDirection(); // FIXME: Drop this line
187 
189  dir1 = dir1.rotate(dir1.orthogonal(),theta1).rotate(dir1, ph);
190  dir2 = dir2.rotate(dir2.orthogonal(),-theta2).rotate(dir2,ph);
191 
192  G4double E=aTrack.GetKineticEnergy();
193  G4double Esec1 = x*E, Esec2 = E-Esec1;
194 
195  // Make FT or ST phonon (0. means no longitudinal)
196  G4int polarization1 = ChoosePolarization(0., theLattice->GetSTDOS(),
197  theLattice->GetFTDOS());
198 
199  // Make FT or ST phonon (0. means no longitudinal)
200  G4int polarization2 = ChoosePolarization(0., theLattice->GetSTDOS(),
201  theLattice->GetFTDOS());
202 
203  // Construct the secondaries and set their wavevectors
204  G4Track* sec1 = CreateSecondary(polarization1, dir1, Esec1);
205  G4Track* sec2 = CreateSecondary(polarization2, dir2, Esec2);
206 
210 }
211 
212 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
213 
214 
215 //Generate daughter phonons from L->L'+T process
216 
218  //d is the velocity ratio vL/v
219  G4double d=1.6338;
220  G4double upperBound=1;
221  G4double lowerBound=(d-1)/(d+1);
222 
223  //Use MC method to generate point from distribution:
224  //if a random point on the energy-probability plane is
225  //smaller that the curve of the probability density,
226  //then accept that point.
227  //x=fraction of parent phonon energy in L phonon
228  G4double x = G4UniformRand()*(upperBound-lowerBound) + lowerBound;
229  G4double p = 4.0*G4UniformRand();
230  while(p >= GetLTDecayProb(d, x)) {
231  x = G4UniformRand()*(upperBound-lowerBound) + lowerBound;
232  p = 4.0*G4UniformRand(); //4.0 is about the max in the PDF
233  }
234 
235  //using energy fraction x to calculate daughter phonon directions
236  G4double thetaL=MakeLDeviation(d, x);
237  G4double thetaT=MakeTDeviation(d, x); // FIXME: Should be 1-x?
238  G4ThreeVector dir1=trackKmap->GetK(aTrack);
239  G4ThreeVector dir2=dir1;
240 
242  dir1 = dir1.rotate(dir1.orthogonal(),thetaL).rotate(dir1, ph);
243  dir2 = dir2.rotate(dir2.orthogonal(),-thetaT).rotate(dir2,ph);
244 
245  G4double E=aTrack.GetKineticEnergy();
246  G4double Esec1 = x*E, Esec2 = E-Esec1;
247 
248  // First secondary is longitudnal
249  G4int polarization1 = G4PhononPolarization::Long;
250 
251  // Make FT or ST phonon (0. means no longitudinal)
252  G4int polarization2 = ChoosePolarization(0., theLattice->GetSTDOS(),
253  theLattice->GetFTDOS());
254 
255  // Construct the secondaries and set their wavevectors
256  G4Track* sec1 = CreateSecondary(polarization1, dir1, Esec1);
257  G4Track* sec2 = CreateSecondary(polarization2, dir2, Esec2);
258 
262 }
263 
264 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
265 
G4PhononDownconversion(const G4String &processName="phononDownconversion")
G4double condition(const G4ErrorSymMatrix &m)
G4double GetGamma() const
const G4LatticePhysical * theLattice
G4double GetMu() const
void MakeTTSecondaries(const G4Track &)
G4double MakeTDeviation(G4double, G4double) const
G4double GetAnhDecConstant() const
G4int verboseLevel
Definition: G4VProcess.hh:368
void MakeLTSecondaries(const G4Track &)
G4double GetTTDecayProb(G4double, G4double) const
CLHEP::Hep3Vector G4ThreeVector
G4double GetVelocity() const
virtual G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *)
G4double GetLambda() const
G4double GetSTDOS() const
G4ThreeVector G4RandomDirection()
Definition of the G4PhononPolarization enum.
Definition of the G4PhononDownconversion class.
double C(double temp)
double B(double temperature)
int G4int
Definition: G4Types.hh:78
G4double GetFTDOS() const
G4double GetLTDecayProb(G4double, G4double) const
G4PhononTrackMap * trackKmap
static constexpr double twopi
Definition: G4SIunits.hh:76
G4double GetKineticEnergy() const
Definition of the G4PhononTrackMap base class.
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
const G4ThreeVector & GetK(const G4Track *track) const
bool G4bool
Definition: G4Types.hh:79
static G4PhononLong * PhononDefinition()
Definition: G4PhononLong.cc:74
Definition: G4Step.hh:76
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
virtual G4bool IsApplicable(const G4ParticleDefinition &)
virtual void Initialize(const G4Track &)
virtual G4int ChoosePolarization(G4double Ldos, G4double STdos, G4double FTdos) const
void SetNumberOfSecondaries(G4int totSecondaries)
void ProposeEnergy(G4double finalEnergy)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
double D(double temp)
Definition of the G4LatticePhysical class.
void AddSecondary(G4Track *aSecondary)
virtual G4Track * CreateSecondary(G4int polarization, const G4ThreeVector &K, G4double energy) const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4ForceCondition
G4double MakeTTDeviation(G4double, G4double) const
G4double GetBeta() const
G4double MakeLDeviation(G4double, G4double) const