Geant4 Cross Reference |
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 #include "G4CrystalExtension.hh" 29 30 #include "G4AtomicFormFactor.hh" 31 32 G4CrystalExtension::G4CrystalExtension(G4Material* mat, const G4String& name) 33 : G4VMaterialExtension(name), fMaterial(mat) 34 {} 35 36 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 37 38 G4complex G4CrystalExtension::ComputeStructureFactor( 39 G4double kScatteringVector, G4int h, G4int k, G4int l) 40 { 41 // SF == Structure Factor 42 // AFF == Atomic Form Factor 43 // GFS == Geometrical Structure Factor 44 G4complex SF = G4complex(0., 0.); 45 46 for (auto& anElement : *(fMaterial->GetElementVector())) { 47 G4double AFF = G4AtomicFormFactor::GetManager()->Get(kScatteringVector, anElement->GetZ()); 48 49 G4complex GFS = G4complex(0., 0.); 50 51 for (const auto& anAtomPos : GetAtomBase(anElement)->GetPos()) { 52 G4double aDouble = h * anAtomPos.x() + k * anAtomPos.y() + l * anAtomPos.z(); 53 GFS += G4complex(std::cos(CLHEP::twopi * aDouble), std::sin(CLHEP::twopi * aDouble)); 54 } 55 56 SF += G4complex(AFF * GFS.real(), AFF * GFS.imag()); 57 } 58 return SF; 59 } 60 61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 62 63 G4complex G4CrystalExtension::ComputeStructureFactorGeometrical(G4int h, G4int k, G4int l) 64 { 65 // GFS == Geometrical Structure Form Factor 66 G4complex GFS = G4complex(0., 0.); 67 68 for (auto& anElement : *(fMaterial->GetElementVector())) { 69 for (const auto& anAtomPos : GetAtomBase(anElement)->GetPos()) { 70 G4double aDouble = h * anAtomPos.x() + k * anAtomPos.y() + l * anAtomPos.z(); 71 GFS += G4complex(std::cos(CLHEP::twopi * aDouble), std::sin(CLHEP::twopi * aDouble)); 72 } 73 } 74 return GFS; 75 } 76 77 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 78 79 void G4CrystalExtension::SetElReduced(const ReducedElasticity& mat) 80 { 81 for (size_t i = 0; i < 6; ++i) { 82 for (size_t j = 0; j < 6; ++j) { 83 fElReduced[i][j] = mat[i][j]; 84 } 85 } 86 } 87 88 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 89 90 void G4CrystalExtension::SetCpq(G4int p, G4int q, G4double value) 91 { 92 if (p > 0 && p < 7 && q > 0 && q < 7) { 93 fElReduced[p - 1][q - 1] = value; 94 } 95 } 96 97 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 98 99 G4CrystalAtomBase* G4CrystalExtension::GetAtomBase(const G4Element* anElement) 100 { 101 if ((theCrystalAtomBaseMap.count(anElement) < 1)) { 102 G4String astring = "Atom base for element " + anElement->GetName() + " is not registered."; 103 G4Exception("G4CrystalExtension::GetAtomBase()", "cry001", JustWarning, astring); 104 105 AddAtomBase(anElement, new G4CrystalAtomBase()); 106 } 107 return theCrystalAtomBaseMap[anElement]; 108 } 109 110 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 111 112 G4bool G4CrystalExtension::GetAtomPos(const G4Element* anEl, std::vector<G4ThreeVector>& vecout) 113 { 114 std::vector<G4ThreeVector> pos; 115 for (auto& asinglepos : GetAtomBase(anEl)->GetPos()) { 116 pos.clear(); 117 theUnitCell->FillAtomicPos(asinglepos, pos); 118 vecout.insert(std::end(vecout), std::begin(pos), std::end(pos)); 119 } 120 return true; 121 } 122 123 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 124 125 G4bool G4CrystalExtension::GetAtomPos(std::vector<G4ThreeVector>& vecout) 126 { 127 std::vector<G4ThreeVector> pos; 128 vecout.clear(); 129 for (auto& anElement : *(fMaterial->GetElementVector())) { 130 pos.clear(); 131 GetAtomPos(anElement, pos); 132 vecout.insert(std::end(vecout), std::begin(pos), std::end(pos)); 133 } 134 return true; 135 } 136