Geant4-11
G4CrystalExtension.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//
27//
28
29//---------------------------------------------------------------------------
30//
31// ClassName: G4CrystalExtension
32//
33// Description: Contains crystal properties
34//
35// Class description:
36//
37// Extension of G4Material for the management of a crystal
38// structure. It has to be attached to a G4ExtendedMaterial
39// in order to instantiate a G4LogicalCrystalVolume.
40//
41
42//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
43
44// 21-04-16, created by E.Bagli
45
46//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
47
48#ifndef G4CrystalExtension_HH
49#define G4CrystalExtension_HH 1
50
52#include "G4CrystalAtomBase.hh"
53#include "G4AtomicBond.hh"
54#include "G4CrystalUnitCell.hh"
55#include "G4NistManager.hh"
56#include <vector>
57
58//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59
61{
62public: // with description
63 //
64 // Constructor to create a material
65 //
66 G4CrystalExtension(G4Material* ,const G4String& name = "crystal");
67
69
70 void Print() const {;};
71private:
73
74public:
76 void SetMaterial(G4Material* aMat) {fMaterial = aMat;};
77
78
79 //
80 // Crystal cell description, i.e. space group
81 // and cell parameters
82 //
83
84private:
86public:
87 inline void SetUnitCell(G4CrystalUnitCell* aUC) {theUnitCell=aUC;}
89 const {return theUnitCell;}
90
91 //
92 // Elasticity and reduced elasticity tensors
93 //
94
95public:
96 typedef G4double Elasticity[3][3][3][3];
98
99protected:
100 Elasticity fElasticity; // Full 4D elasticity tensor
101 ReducedElasticity fElReduced; // Reduced 2D elasticity tensor
102
103public:
104 const Elasticity& GetElasticity() const { return fElasticity; }
105 const ReducedElasticity& GetElReduced() const { return fElReduced; }
106
107public:
109 return fElasticity[i][j][k][l];
110 }
111
112 // Reduced elasticity tensor: C11-C66 interface for clarity
113 void SetElReduced(const ReducedElasticity& mat);
114
115 void SetCpq(G4int p, G4int q, G4double value);
116 G4double GetCpq(G4int p, G4int q) const { return fElReduced[p-1][q-1]; }
117
118 //
119 // Map of atom positions for each element
120 // The coordinate system is the unit cell
121 //
122private:
123 std::map<const G4Element*,G4CrystalAtomBase*> theCrystalAtomBaseMap;
124
125public:
126 G4CrystalAtomBase* GetAtomBase(const G4Element* anElement);
127
128 void AddAtomBase(const G4Element* anElement,G4CrystalAtomBase* aBase){
129 theCrystalAtomBaseMap.insert(std::pair<const G4Element*,G4CrystalAtomBase*>(anElement,aBase));
130 }
131
132
134 return GetAtomBase(fMaterial->GetElement(anElIdx));
135 }
136
137 void AddAtomBase(G4int anElIdx,G4CrystalAtomBase* aLattice){
138 AddAtomBase(fMaterial->GetElement(anElIdx),aLattice);
139 }
140
141 //
142 // Get the position of all the atoms in the unit cell
143 // for a single element or all the elements
144 //
145 G4bool GetAtomPos(const G4Element* anEl, std::vector<G4ThreeVector>& vecout);
146 G4bool GetAtomPos(std::vector<G4ThreeVector>& vecout);
147
148 G4bool GetAtomPos(G4int anElIdx, std::vector<G4ThreeVector>& vecout){
149 GetAtomPos(fMaterial->GetElement(anElIdx),vecout);
150 return true;
151 }
152
153
154 //
155 // Structure factor calculations
156 // Eq. 46, Chapter 2 , Introduction to solid state physics, C. Kittel
157 //
158 G4complex ComputeStructureFactor(G4double kScatteringVector,
159 G4int h,
160 G4int k,
161 G4int l);
163 G4int k,
164 G4int l);
165
166 //
167 // Bond between atoms. Each bond is mapped with two Elements
168 // and the number of the atoms in the corresponding G4CrystalBaseAtomPos
169 //
170
171private:
172 std::vector<G4AtomicBond*> theAtomicBondVector;
173
174public:
175 void AddAtomicBond(G4AtomicBond* aBond) {theAtomicBondVector.push_back(aBond);};
177 std::vector<G4AtomicBond*> GetAtomicBondVector() {return theAtomicBondVector;};
178
179};
180
181//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
182
183#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
std::complex< G4double > G4complex
Definition: G4Types.hh:88
int G4int
Definition: G4Types.hh:85
G4double GetCijkl(G4int i, G4int j, G4int k, G4int l) const
void AddAtomicBond(G4AtomicBond *aBond)
void SetElReduced(const ReducedElasticity &mat)
G4double ReducedElasticity[6][6]
std::vector< G4AtomicBond * > GetAtomicBondVector()
G4complex ComputeStructureFactorGeometrical(G4int h, G4int k, G4int l)
G4CrystalExtension(G4Material *, const G4String &name="crystal")
G4CrystalAtomBase * GetAtomBase(const G4Element *anElement)
G4double GetCpq(G4int p, G4int q) const
void AddAtomBase(const G4Element *anElement, G4CrystalAtomBase *aBase)
const ReducedElasticity & GetElReduced() const
const Elasticity & GetElasticity() const
std::map< const G4Element *, G4CrystalAtomBase * > theCrystalAtomBaseMap
G4Material * GetMaterial()
G4double Elasticity[3][3][3][3]
std::vector< G4AtomicBond * > theAtomicBondVector
G4CrystalUnitCell * GetUnitCell() const
G4AtomicBond * GetAtomicBond(G4int idx)
G4complex ComputeStructureFactor(G4double kScatteringVector, G4int h, G4int k, G4int l)
G4bool GetAtomPos(const G4Element *anEl, std::vector< G4ThreeVector > &vecout)
G4CrystalUnitCell * theUnitCell
ReducedElasticity fElReduced
void AddAtomBase(G4int anElIdx, G4CrystalAtomBase *aLattice)
void SetMaterial(G4Material *aMat)
void SetUnitCell(G4CrystalUnitCell *aUC)
void SetCpq(G4int p, G4int q, G4double value)
G4bool GetAtomPos(G4int anElIdx, std::vector< G4ThreeVector > &vecout)
G4CrystalAtomBase * GetAtomBase(G4int anElIdx)
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:198
const char * name(G4int ptype)