Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// ClassName: G4CrystalExtension
29//
30// Description: Contains crystal properties
31//
32// Class description:
33//
34// Extension of G4Material for the management of a crystal
35// structure. It has to be attached to a G4ExtendedMaterial
36// in order to instantiate a G4LogicalCrystalVolume.
37//
38// 21-04-16, created by E.Bagli
39
40#ifndef G4CrystalExtension_HH
41#define G4CrystalExtension_HH 1
42
43#include "G4AtomicBond.hh"
44#include "G4CrystalAtomBase.hh"
45#include "G4CrystalUnitCell.hh"
46#include "G4NistManager.hh"
48
49#include <vector>
50
52{
53 public:
54 // Elasticity and reduced elasticity tensors
55 using Elasticity = G4double[3][3][3][3];
57
58 public: // with description
59 // Constructor to create a material
60 G4CrystalExtension(G4Material*, const G4String& name = "crystal");
61
62 ~G4CrystalExtension() override = default;
63
64 void Print() const override { ; };
65
66 G4Material* GetMaterial() { return fMaterial; };
67 void SetMaterial(G4Material* aMat) { fMaterial = aMat; };
68
69 inline void SetUnitCell(G4CrystalUnitCell* aUC) { theUnitCell = aUC; }
70 inline G4CrystalUnitCell* GetUnitCell() const { return theUnitCell; }
71
72 const Elasticity& GetElasticity() const { return fElasticity; }
73 const ReducedElasticity& GetElReduced() const { return fElReduced; }
74
75 G4double GetCijkl(G4int i, G4int j, G4int k, G4int l) const { return fElasticity[i][j][k][l]; }
76
77 // Reduced elasticity tensor: C11-C66 interface for clarity
78 void SetElReduced(const ReducedElasticity& mat);
79
80 void SetCpq(G4int p, G4int q, G4double value);
81 G4double GetCpq(G4int p, G4int q) const { return fElReduced[p - 1][q - 1]; }
82
83 G4CrystalAtomBase* GetAtomBase(const G4Element* anElement);
84 void AddAtomBase(const G4Element* anElement, G4CrystalAtomBase* aBase)
85 {
86 theCrystalAtomBaseMap.insert(std::pair<const G4Element*, G4CrystalAtomBase*>(anElement, aBase));
87 }
88
90 {
91 return GetAtomBase(fMaterial->GetElement(anElIdx));
92 }
93
94 void AddAtomBase(G4int anElIdx, G4CrystalAtomBase* aLattice)
95 {
96 AddAtomBase(fMaterial->GetElement(anElIdx), aLattice);
97 }
98
99 // Get the position of all the atoms in the unit cell
100 // for a single element or all the elements
101 G4bool GetAtomPos(const G4Element* anEl, std::vector<G4ThreeVector>& vecout);
102 G4bool GetAtomPos(std::vector<G4ThreeVector>& vecout);
103
104 G4bool GetAtomPos(G4int anElIdx, std::vector<G4ThreeVector>& vecout)
105 {
106 GetAtomPos(fMaterial->GetElement(anElIdx), vecout);
107 return true;
108 }
109
110 // Structure factor calculations
111 // Eq. 46, Chapter 2 , Introduction to solid state physics, C. Kittel
112 G4complex ComputeStructureFactor(G4double kScatteringVector, G4int h, G4int k, G4int l);
114
115 void AddAtomicBond(G4AtomicBond* aBond) { theAtomicBondVector.push_back(aBond); };
116 G4AtomicBond* GetAtomicBond(G4int idx) { return theAtomicBondVector[idx]; };
117 std::vector<G4AtomicBond*> GetAtomicBondVector() { return theAtomicBondVector; };
118
119 protected:
120 Elasticity fElasticity; // Full 4D elasticity tensor
121 ReducedElasticity fElReduced; // Reduced 2D elasticity tensor
122
123 private:
124 G4Material* fMaterial;
125
126 // Crystal cell description, i.e. space group
127 // and cell parameters
128 G4CrystalUnitCell* theUnitCell{nullptr};
129
130 // Map of atom positions for each element
131 // The coordinate system is the unit cell
132 std::map<const G4Element*, G4CrystalAtomBase*> theCrystalAtomBaseMap;
133
134 // Bond between atoms. Each bond is mapped with two Elements
135 // and the number of the atoms in the corresponding G4CrystalBaseAtomPos
136 std::vector<G4AtomicBond*> theAtomicBondVector;
137};
138
139#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)
~G4CrystalExtension() override=default
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
G4double[3][3][3][3] Elasticity
G4Material * GetMaterial()
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)
ReducedElasticity fElReduced
G4double[6][6] ReducedElasticity
void AddAtomBase(G4int anElIdx, G4CrystalAtomBase *aLattice)
void SetMaterial(G4Material *aMat)
void Print() const override
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