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