Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4CrystalExtension.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//
26
27//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
28
29// 21-04-16, created by E.Bagli
30
31//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32
33#include "G4CrystalExtension.hh"
34#include "G4AtomicFormFactor.hh"
35
38 , fMaterial(mat)
39 , theUnitCell(nullptr)
40{}
41
42//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
43
45{}
46
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48
50ComputeStructureFactor(G4double kScatteringVector,
51 G4int h,
52 G4int k,
53 G4int l){
54 //SF == Structure Factor
55 //AFF == Atomic Form Factor
56 //GFS == Geometrical Structure Factor
57 G4complex SF = G4complex(0.,0.);
58
59 for(auto & anElement: *(fMaterial->GetElementVector())){
60 G4double AFF = G4AtomicFormFactor::GetManager()->Get(kScatteringVector,anElement->GetZ());
61
62 G4complex GFS = G4complex(0.,0.);
63
64 for(const auto& anAtomPos : GetAtomBase(anElement)->GetPos())
65 {
66 G4double aDouble = h * anAtomPos.x()
67 + k * anAtomPos.y()
68 + l * anAtomPos.z();
69 GFS += G4complex(std::cos(CLHEP::twopi * aDouble),
70 std::sin(CLHEP::twopi * aDouble));
71 }
72
73
74 SF += G4complex(AFF * GFS.real(),AFF * GFS.imag());
75 }
76 return SF;
77}
78
79//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
80
83 G4int k,
84 G4int l){
85 //GFS == Geometrical Structure Form Factor
86 G4complex GFS = G4complex(0.,0.);
87
88 for(auto & anElement: *(fMaterial->GetElementVector())){
89 for(const auto& anAtomPos : GetAtomBase(anElement)->GetPos())
90 {
91 G4double aDouble =
92 h * anAtomPos.x() + k * anAtomPos.y() + l * anAtomPos.z();
93 GFS += G4complex(std::cos(CLHEP::twopi * aDouble),
94 std::sin(CLHEP::twopi * aDouble));
95 }
96 }
97 return GFS;
98}
99
100//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
101
102void G4CrystalExtension::SetElReduced(const ReducedElasticity& mat) {
103 for (size_t i=0; i<6; ++i) {
104 for (size_t j=0; j<6; ++j) {
105 fElReduced[i][j] = mat[i][j];
106 }
107 }
108}
109
110//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
111
113 if(p > 0 && p < 7 && q > 0 && q < 7)
114 {
115 fElReduced[p - 1][q - 1] = value;
116 }
117}
118
119//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
120
122 if((theCrystalAtomBaseMap.count(anElement)<1)){
123 G4String astring = "Atom base for element " + anElement->GetName()
124 + " is not registered." ;
125 G4Exception ("G4CrystalExtension::GetAtomBase()", "cry001", JustWarning,astring);
126
127 AddAtomBase(anElement, new G4CrystalAtomBase());
128 }
129 return theCrystalAtomBaseMap[anElement];
130}
131
132//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
133
134G4bool G4CrystalExtension::GetAtomPos(const G4Element* anEl, std::vector<G4ThreeVector>& vecout){
135 std::vector<G4ThreeVector> pos;
136 for(auto & asinglepos: GetAtomBase(anEl)->GetPos()){
137 pos.clear();
138 theUnitCell->FillAtomicPos(asinglepos,pos);
139 vecout.insert(std::end(vecout), std::begin(pos), std::end(pos));
140 }
141 return true;
142}
143
144//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
145
146G4bool G4CrystalExtension::GetAtomPos(std::vector<G4ThreeVector>& vecout){
147 std::vector<G4ThreeVector> pos;
148 vecout.clear();
149 for(auto & anElement: *(fMaterial->GetElementVector())){
150 pos.clear();
151 GetAtomPos(anElement,pos);
152 vecout.insert(std::end(vecout), std::begin(pos), std::end(pos));
153 }
154 return true;
155}
156
157//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
158
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
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 Get(G4double kScatteringVector, G4int Z, G4int charge=0)
static G4AtomicFormFactor * GetManager()
void SetElReduced(const ReducedElasticity &mat)
G4complex ComputeStructureFactorGeometrical(G4int h, G4int k, G4int l)
G4CrystalExtension(G4Material *, const G4String &name="crystal")
G4CrystalAtomBase * GetAtomBase(const G4Element *anElement)
void AddAtomBase(const G4Element *anElement, G4CrystalAtomBase *aBase)
G4complex ComputeStructureFactor(G4double kScatteringVector, G4int h, G4int k, G4int l)
G4bool GetAtomPos(const G4Element *anEl, std::vector< G4ThreeVector > &vecout)
ReducedElasticity fElReduced
~G4CrystalExtension() override
void SetCpq(G4int p, G4int q, G4double value)
G4bool FillAtomicPos(G4ThreeVector &pos, std::vector< G4ThreeVector > &vecout)
const G4String & GetName() const
Definition: G4Element.hh:127
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185