Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4LatticeLogical.hh
이 파일의 문서화 페이지로 가기
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: G4LatticeLogical.hh 76883 2013-11-18 12:50:08Z gcosmo $
30 //
31 // 20131114 Add verbosity for diagnostic output
32 // 20131115 Expose maximum array dimensions for use by LatticeReader
33 
34 #ifndef G4LatticeLogical_h
35 #define G4LatticeLogical_h
36 
37 #include "globals.hh"
38 #include "G4ThreeVector.hh"
39 #include <iosfwd>
40 
41 
43 public:
45  virtual ~G4LatticeLogical();
46 
47  void SetVerboseLevel(G4int vb) { verboseLevel = vb; }
48 
51 
52  // Get group velocity magnitude for input polarization and wavevector
53  virtual G4double MapKtoV(G4int, const G4ThreeVector& ) const;
54 
55  // Get group velocity direction (unit vector) for input polarization and K
56  virtual G4ThreeVector MapKtoVDir(G4int, const G4ThreeVector& ) const;
57 
58  // Dump structure in format compatible with reading back
59  void Dump(std::ostream& os) const;
60  void DumpMap(std::ostream& os, G4int pol, const G4String& name) const;
61  void Dump_NMap(std::ostream& os, G4int pol, const G4String& name) const;
62 
63 public:
65  G4double Lambda, G4double Mu) {
66  fBeta=Beta; fGamma=Gamma; fLambda=Lambda; fMu=Mu;
67  }
68 
71  void SetLDOS(G4double LDOS) { fLDOS=LDOS; }
72  void SetSTDOS(G4double STDOS) { fSTDOS=STDOS; }
73  void SetFTDOS(G4double FTDOS) { fFTDOS=FTDOS; }
74 
75  G4double GetBeta() const { return fBeta; }
76  G4double GetGamma() const { return fGamma; }
77  G4double GetLambda() const { return fLambda; }
78  G4double GetMu() const { return fMu; }
79  G4double GetScatteringConstant() const { return fB; }
80  G4double GetAnhDecConstant() const { return fA; }
81  G4double GetLDOS() const { return fLDOS; }
82  G4double GetSTDOS() const { return fSTDOS; }
83  G4double GetFTDOS() const { return fFTDOS; }
84 
85 public:
86  enum { MAXRES=322 }; // Maximum map resolution (bins)
87 
88 private:
89  G4int verboseLevel; // Enable diagnostic output
90 
91  G4double fMap[3][MAXRES][MAXRES]; // map for group velocity scalars
92  G4ThreeVector fN_map[3][MAXRES][MAXRES]; // map for direction vectors
93 
94  G4int fVresTheta; //velocity map theta resolution (inclination)
95  G4int fVresPhi; //velocity map phi resolution (azimuth)
96  G4int fDresTheta; //direction map theta resn
97  G4int fDresPhi; //direction map phi resn
98 
99  G4double fA; //Scaling constant for Anh.Dec. mean free path
100  G4double fB; //Scaling constant for Iso.Scat. mean free path
101  G4double fLDOS; //Density of states for L-phonons
102  G4double fSTDOS; //Density of states for ST-phonons
103  G4double fFTDOS; //Density of states for FT-phonons
104  G4double fBeta, fGamma, fLambda, fMu; //dynamical constants for material
105 };
106 
107 // Write lattice structure to output stream
108 
109 inline std::ostream&
110 operator<<(std::ostream& os, const G4LatticeLogical& lattice) {
111  lattice.Dump(os);
112  return os;
113 }
114 
115 #endif
G4double GetFTDOS() const
const XML_Char * name
Definition: expat.h:151
G4bool LoadMap(G4int, G4int, G4int, G4String)
G4double GetGamma() const
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
std::ostream & operator<<(std::ostream &, const BasicVector3D< float > &)
G4double GetLambda() const
G4double GetAnhDecConstant() const
virtual G4double MapKtoV(G4int, const G4ThreeVector &) const
void SetSTDOS(G4double STDOS)
void SetScatteringConstant(G4double b)
G4double GetLDOS() const
void Dump_NMap(std::ostream &os, G4int pol, const G4String &name) const
virtual G4ThreeVector MapKtoVDir(G4int, const G4ThreeVector &) const
G4double GetBeta() const
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4double GetSTDOS() const
G4bool Load_NMap(G4int, G4int, G4int, G4String)
void Dump(std::ostream &os) const
G4double GetScatteringConstant() const
void SetAnhDecConstant(G4double a)
void SetFTDOS(G4double FTDOS)
G4double GetMu() const
G4double fMap[3][MAXRES][MAXRES]
G4ThreeVector fN_map[3][MAXRES][MAXRES]
int G4int
Definition: G4Types.hh:78
void SetVerboseLevel(G4int vb)
void DumpMap(std::ostream &os, G4int pol, const G4String &name) const
void SetDynamicalConstants(G4double Beta, G4double Gamma, G4double Lambda, G4double Mu)
void SetLDOS(G4double LDOS)
virtual ~G4LatticeLogical()