Geant4  10.01.p01
G4VPhononProcess.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: G4VPhononProcess.cc 76885 2013-11-18 12:55:15Z gcosmo $
30 //
31 // 20131111 Add verbosity to report creating secondaries
32 
33 #include "G4VPhononProcess.hh"
34 #include "G4DynamicParticle.hh"
35 #include "G4ExceptionSeverity.hh"
36 #include "G4LatticeManager.hh"
37 #include "G4LatticePhysical.hh"
38 #include "G4ParticleDefinition.hh"
39 #include "G4PhononLong.hh"
40 #include "G4PhononPolarization.hh"
41 #include "G4PhononTrackMap.hh"
42 #include "G4PhononTransFast.hh"
43 #include "G4PhononTransSlow.hh"
44 #include "G4ProcessType.hh"
45 #include "G4ThreeVector.hh"
46 #include "G4Track.hh"
47 
48 namespace {
49  const G4ThreeVector nullVec(0.,0.,0.); // For convenience below
50 }
51 
52 // Constructor and destructor
53 
55  : G4VDiscreteProcess(processName, fPhonon),
56  trackKmap(G4PhononTrackMap::GetInstance()), theLattice(0),
57  currentTrack(0) {}
58 
60 
61 
62 // Only applies to the known phonon polarization states
63 
65  return (&aPD==G4PhononLong::Definition() ||
68 }
69 
70 
71 // Initialize wave vectors for currently active track(s)
72 
74  G4VProcess::StartTracking(track); // Apply base class actions
75 
76  // FIXME: THE WAVEVECTOR SHOULD BE COMPUTED BY INVERTING THE K/V MAP
77  if (!trackKmap->Find(track))
78  trackKmap->SetK(track, track->GetMomentumDirection());
79 
80  currentTrack = track; // Save for use by EndTracking
81 
82  // Fetch lattice for current track once, use in subsequent steps
84  theLattice = LM->GetLattice(track->GetVolume());
85 }
86 
88  G4VProcess::EndTracking(); // Apply base class actions
90  currentTrack = 0;
91  theLattice = 0;
92 }
93 
94 
95 // For convenience, map phonon type to polarization code
96 
99 }
100 
101 
102 // Generate random polarization from density of states
103 
105  G4double FTdos) const {
106  G4double norm = Ldos + STdos + FTdos;
107  G4double cProbST = STdos/norm;
108  G4double cProbFT = FTdos/norm + cProbST;
109 
110  // NOTE: Order of selection done to match previous random sequences
111  G4double modeMixer = G4UniformRand();
112  if (modeMixer<cProbST) return G4PhononPolarization::TransSlow;
113  if (modeMixer<cProbFT) return G4PhononPolarization::TransFast;
115 }
116 
117 
118 // Create new secondary track from phonon configuration
119 
121  const G4ThreeVector& waveVec,
122  G4double energy) const {
123  if (verboseLevel>1) {
124  G4cout << GetProcessName() << " CreateSecondary pol " << polarization
125  << " K " << waveVec << " E " << energy << G4endl;
126  }
127 
128  G4ThreeVector vgroup = theLattice->MapKtoVDir(polarization, waveVec);
129  if (verboseLevel>1) G4cout << " MapKtoVDir returned " << vgroup << G4endl;
130 
131  vgroup = theLattice->RotateToGlobal(vgroup);
132  if (verboseLevel>1) G4cout << " RotateToGlobal returned " << vgroup << G4endl;
133 
134  if (verboseLevel && std::fabs(vgroup.mag()-1.) > 0.01) {
135  G4cout << "WARNING: " << GetProcessName() << " vgroup not a unit vector: "
136  << vgroup << G4endl;
137  }
138 
139  G4ParticleDefinition* thePhonon = G4PhononPolarization::Get(polarization);
140 
141  // Secondaries are created at the current track coordinates
142  G4Track* sec = new G4Track(new G4DynamicParticle(thePhonon, vgroup, energy),
145 
146  // Store wavevector in lookup table for future tracking
147  trackKmap->SetK(sec, theLattice->RotateToGlobal(waveVec));
148 
149  if (verboseLevel>1) {
150  G4cout << GetProcessName() << " secondary K rotated to "
151  << trackKmap->GetK(sec) << G4endl;
152  }
153 
154  sec->SetVelocity(theLattice->MapKtoV(polarization, waveVec));
155  sec->UseGivenVelocity(true);
156 
157  return sec;
158 }
static G4PhononLong * Definition()
Definition: G4PhononLong.cc:40
const G4LatticePhysical * theLattice
virtual G4bool IsApplicable(const G4ParticleDefinition &aPD)
static G4LatticeManager * GetLatticeManager()
void SetVelocity(G4double val)
G4int verboseLevel
Definition: G4VProcess.hh:368
CLHEP::Hep3Vector G4ThreeVector
const G4ThreeVector & GetPosition() const
void RemoveTrack(const G4Track *track)
G4double MapKtoV(G4int, G4ThreeVector) const
Definition of the G4PhononPolarization enum.
int G4int
Definition: G4Types.hh:78
virtual ~G4VPhononProcess()
G4PhononTrackMap * trackKmap
G4VPhononProcess(const G4String &processName)
G4bool Find(const G4Track *track) const
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:101
void SetK(const G4Track *track, const G4ThreeVector &K)
Definition of the G4PhononTrackMap base class.
Definition of the G4VPhononProcess base class.
#define G4UniformRand()
Definition: Randomize.hh:95
G4GLOB_DLL std::ostream G4cout
static G4PhononTransFast * Definition()
const G4ThreeVector & GetK(const G4Track *track) const
bool G4bool
Definition: G4Types.hh:79
G4ThreeVector MapKtoVDir(G4int, G4ThreeVector) const
const G4ParticleDefinition * GetParticleDefinition() const
G4int Get(const G4ParticleDefinition *aPD)
const G4Track * currentTrack
G4double GetGlobalTime() const
G4bool UseGivenVelocity() const
virtual G4int GetPolarization(const G4Track &track) const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
static G4PhononTransSlow * Definition()
virtual G4int ChoosePolarization(G4double Ldos, G4double STdos, G4double FTdos) const
const G4ThreeVector & GetMomentumDirection() const
G4double energy(const ThreeVector &p, const G4double m)
virtual void EndTracking()
virtual void StartTracking(G4Track *track)
G4VPhysicalVolume * GetVolume() const
Definition of the G4LatticePhysical class.
virtual G4Track * CreateSecondary(G4int polarization, const G4ThreeVector &K, G4double energy) const
G4ThreeVector RotateToGlobal(const G4ThreeVector &dir) const
#define G4endl
Definition: G4ios.hh:61
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)
double G4double
Definition: G4Types.hh:76
virtual void EndTracking()
Definition: G4VProcess.cc:113
G4LatticeLogical * GetLattice(G4Material *) const