Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LatticeLogical.hh
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/// \file materials/include/G4LatticeLogical.hh
27/// \brief Definition of the G4LatticeLogical class
28//
29// 20131114 Add verbosity for diagnostic output
30// 20131115 Expose maximum array dimensions for use by LatticeReader
31
32#ifndef G4LatticeLogical_h
33#define G4LatticeLogical_h
34
35#include "G4ThreeVector.hh"
36#include "globals.hh"
37
38#include <iosfwd>
39
41{
42 public:
43 enum
44 {
45 MAXRES = 322
46 }; // Maximum map resolution (bins)
47
48 public:
50 virtual ~G4LatticeLogical() = default;
51
52 void SetVerboseLevel(G4int vb) { verboseLevel = vb; }
53
56
57 // Get group velocity magnitude for input polarization and wavevector
58 virtual G4double MapKtoV(G4int, const G4ThreeVector&) const;
59
60 // Get group velocity direction (unit vector) for input polarization and K
61 virtual G4ThreeVector MapKtoVDir(G4int, const G4ThreeVector&) const;
62
63 // Dump structure in format compatible with reading back
64 void Dump(std::ostream& os) const;
65 void DumpMap(std::ostream& os, G4int pol, const G4String& name) const;
66 void Dump_NMap(std::ostream& os, G4int pol, const G4String& name) const;
67
69 {
70 fBeta = Beta;
71 fGamma = Gamma;
72 fLambda = Lambda;
73 fMu = Mu;
74 }
75
76 void SetScatteringConstant(G4double b) { fB = b; }
77 void SetAnhDecConstant(G4double a) { fA = a; }
78 void SetLDOS(G4double LDOS) { fLDOS = LDOS; }
79 void SetSTDOS(G4double STDOS) { fSTDOS = STDOS; }
80 void SetFTDOS(G4double FTDOS) { fFTDOS = FTDOS; }
81
82 G4double GetBeta() const { return fBeta; }
83 G4double GetGamma() const { return fGamma; }
84 G4double GetLambda() const { return fLambda; }
85 G4double GetMu() const { return fMu; }
86 G4double GetScatteringConstant() const { return fB; }
87 G4double GetAnhDecConstant() const { return fA; }
88 G4double GetLDOS() const { return fLDOS; }
89 G4double GetSTDOS() const { return fSTDOS; }
90 G4double GetFTDOS() const { return fFTDOS; }
91
92 private:
93 G4int verboseLevel{0}; // Enable diagnostic output
94
95 G4double fMap[3][MAXRES][MAXRES]; // map for group velocity scalars
96 G4ThreeVector fN_map[3][MAXRES][MAXRES]; // map for direction vectors
97
98 G4int fVresTheta{0}; // velocity map theta resolution (inclination)
99 G4int fVresPhi{0}; // velocity map phi resolution (azimuth)
100 G4int fDresTheta{0}; // direction map theta resn
101 G4int fDresPhi{0}; // direction map phi resn
102
103 G4double fA{0}; // Scaling constant for Anh.Dec. mean free path
104 G4double fB{0}; // Scaling constant for Iso.Scat. mean free path
105 G4double fLDOS{0}; // Density of states for L-phonons
106 G4double fSTDOS{0}; // Density of states for ST-phonons
107 G4double fFTDOS{0}; // Density of states for FT-phonons
108 G4double fBeta{0}, fGamma{0}, fLambda{0}, fMu{0}; // dynamical constants for material
109};
110
111// Write lattice structure to output stream
112inline std::ostream& operator<<(std::ostream& os, const G4LatticeLogical& lattice)
113{
114 lattice.Dump(os);
115 return os;
116}
117
118#endif
std::ostream & operator<<(std::ostream &os, const G4LatticeLogical &lattice)
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void SetLDOS(G4double LDOS)
G4double GetLDOS() const
void SetAnhDecConstant(G4double a)
virtual G4double MapKtoV(G4int, const G4ThreeVector &) const
G4double GetMu() const
void SetFTDOS(G4double FTDOS)
G4double GetLambda() const
G4double GetGamma() const
G4double GetBeta() const
G4double GetScatteringConstant() const
void DumpMap(std::ostream &os, G4int pol, const G4String &name) const
G4double GetSTDOS() const
G4bool LoadMap(G4int, G4int, G4int, G4String)
G4double GetAnhDecConstant() const
void SetSTDOS(G4double STDOS)
virtual ~G4LatticeLogical()=default
void SetDynamicalConstants(G4double Beta, G4double Gamma, G4double Lambda, G4double Mu)
void SetScatteringConstant(G4double b)
void SetVerboseLevel(G4int vb)
void Dump_NMap(std::ostream &os, G4int pol, const G4String &name) const
virtual G4ThreeVector MapKtoVDir(G4int, const G4ThreeVector &) const
G4double GetFTDOS() const
void Dump(std::ostream &os) const
G4bool Load_NMap(G4int, G4int, G4int, G4String)