Geant4-11
Public Types | Public Member Functions | Protected Attributes | Private Attributes
G4CrystalExtension Class Reference

#include <G4CrystalExtension.hh>

Inheritance diagram for G4CrystalExtension:
G4VMaterialExtension

Public Types

typedef G4double Elasticity[3][3][3][3]
 
typedef G4double ReducedElasticity[6][6]
 

Public Member Functions

void AddAtomBase (const G4Element *anElement, G4CrystalAtomBase *aBase)
 
void AddAtomBase (G4int anElIdx, G4CrystalAtomBase *aLattice)
 
void AddAtomicBond (G4AtomicBond *aBond)
 
G4complex ComputeStructureFactor (G4double kScatteringVector, G4int h, G4int k, G4int l)
 
G4complex ComputeStructureFactorGeometrical (G4int h, G4int k, G4int l)
 
 G4CrystalExtension (G4Material *, const G4String &name="crystal")
 
G4CrystalAtomBaseGetAtomBase (const G4Element *anElement)
 
G4CrystalAtomBaseGetAtomBase (G4int anElIdx)
 
G4AtomicBondGetAtomicBond (G4int idx)
 
std::vector< G4AtomicBond * > GetAtomicBondVector ()
 
G4bool GetAtomPos (const G4Element *anEl, std::vector< G4ThreeVector > &vecout)
 
G4bool GetAtomPos (G4int anElIdx, std::vector< G4ThreeVector > &vecout)
 
G4bool GetAtomPos (std::vector< G4ThreeVector > &vecout)
 
G4double GetCijkl (G4int i, G4int j, G4int k, G4int l) const
 
G4double GetCpq (G4int p, G4int q) const
 
const ElasticityGetElasticity () const
 
const ReducedElasticityGetElReduced () const
 
const std::size_t & GetHash () const
 
G4MaterialGetMaterial ()
 
const G4StringGetName () const
 
G4CrystalUnitCellGetUnitCell () const
 
void Print () const
 
void SetCpq (G4int p, G4int q, G4double value)
 
void SetElReduced (const ReducedElasticity &mat)
 
void SetMaterial (G4Material *aMat)
 
void SetUnitCell (G4CrystalUnitCell *aUC)
 
 ~G4CrystalExtension ()
 

Protected Attributes

Elasticity fElasticity
 
ReducedElasticity fElReduced
 
const std::size_t fHash
 
const G4StringfName
 

Private Attributes

G4MaterialfMaterial
 
std::vector< G4AtomicBond * > theAtomicBondVector
 
std::map< const G4Element *, G4CrystalAtomBase * > theCrystalAtomBaseMap
 
G4CrystalUnitCelltheUnitCell
 

Detailed Description

Definition at line 60 of file G4CrystalExtension.hh.

Member Typedef Documentation

◆ Elasticity

typedef G4double G4CrystalExtension::Elasticity[3][3][3][3]

Definition at line 96 of file G4CrystalExtension.hh.

◆ ReducedElasticity

typedef G4double G4CrystalExtension::ReducedElasticity[6][6]

Definition at line 97 of file G4CrystalExtension.hh.

Constructor & Destructor Documentation

◆ G4CrystalExtension()

G4CrystalExtension::G4CrystalExtension ( G4Material mat,
const G4String name = "crystal" 
)

Definition at line 39 of file G4CrystalExtension.cc.

39 :
41fMaterial(mat),
42theUnitCell(0){;}
G4CrystalUnitCell * theUnitCell
G4VMaterialExtension(const G4String &name)
const char * name(G4int ptype)

◆ ~G4CrystalExtension()

G4CrystalExtension::~G4CrystalExtension ( )

Definition at line 46 of file G4CrystalExtension.cc.

46{;}

Member Function Documentation

◆ AddAtomBase() [1/2]

void G4CrystalExtension::AddAtomBase ( const G4Element anElement,
G4CrystalAtomBase aBase 
)
inline

Definition at line 128 of file G4CrystalExtension.hh.

128 {
129 theCrystalAtomBaseMap.insert(std::pair<const G4Element*,G4CrystalAtomBase*>(anElement,aBase));
130 }
std::map< const G4Element *, G4CrystalAtomBase * > theCrystalAtomBaseMap

References theCrystalAtomBaseMap.

Referenced by AddAtomBase(), and GetAtomBase().

◆ AddAtomBase() [2/2]

void G4CrystalExtension::AddAtomBase ( G4int  anElIdx,
G4CrystalAtomBase aLattice 
)
inline

Definition at line 137 of file G4CrystalExtension.hh.

137 {
138 AddAtomBase(fMaterial->GetElement(anElIdx),aLattice);
139 }
void AddAtomBase(const G4Element *anElement, G4CrystalAtomBase *aBase)
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:198

References AddAtomBase(), fMaterial, and G4Material::GetElement().

◆ AddAtomicBond()

void G4CrystalExtension::AddAtomicBond ( G4AtomicBond aBond)
inline

Definition at line 175 of file G4CrystalExtension.hh.

175{theAtomicBondVector.push_back(aBond);};
std::vector< G4AtomicBond * > theAtomicBondVector

References theAtomicBondVector.

◆ ComputeStructureFactor()

G4complex G4CrystalExtension::ComputeStructureFactor ( G4double  kScatteringVector,
G4int  h,
G4int  k,
G4int  l 
)

Definition at line 50 of file G4CrystalExtension.cc.

54 {
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}
double G4double
Definition: G4Types.hh:83
std::complex< G4double > G4complex
Definition: G4Types.hh:88
G4double Get(G4double kScatteringVector, G4int Z, G4int charge=0)
static G4AtomicFormFactor * GetManager()
G4CrystalAtomBase * GetAtomBase(const G4Element *anElement)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:186
static constexpr double pi
Definition: SystemOfUnits.h:55

References fMaterial, G4AtomicFormFactor::Get(), GetAtomBase(), G4Material::GetElementVector(), G4AtomicFormFactor::GetManager(), and CLHEP::pi.

◆ ComputeStructureFactorGeometrical()

G4complex G4CrystalExtension::ComputeStructureFactorGeometrical ( G4int  h,
G4int  k,
G4int  l 
)

Definition at line 82 of file G4CrystalExtension.cc.

85 {
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}

References fMaterial, GetAtomBase(), G4Material::GetElementVector(), and CLHEP::pi.

◆ GetAtomBase() [1/2]

G4CrystalAtomBase * G4CrystalExtension::GetAtomBase ( const G4Element anElement)

Definition at line 120 of file G4CrystalExtension.cc.

120 {
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}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
const G4String & GetName() const
Definition: G4Element.hh:127

References AddAtomBase(), G4Exception(), G4Element::GetName(), JustWarning, and theCrystalAtomBaseMap.

Referenced by ComputeStructureFactor(), ComputeStructureFactorGeometrical(), GetAtomBase(), and GetAtomPos().

◆ GetAtomBase() [2/2]

G4CrystalAtomBase * G4CrystalExtension::GetAtomBase ( G4int  anElIdx)
inline

Definition at line 133 of file G4CrystalExtension.hh.

133 {
134 return GetAtomBase(fMaterial->GetElement(anElIdx));
135 }

References fMaterial, GetAtomBase(), and G4Material::GetElement().

◆ GetAtomicBond()

G4AtomicBond * G4CrystalExtension::GetAtomicBond ( G4int  idx)
inline

Definition at line 176 of file G4CrystalExtension.hh.

176{return theAtomicBondVector[idx];};

References theAtomicBondVector.

◆ GetAtomicBondVector()

std::vector< G4AtomicBond * > G4CrystalExtension::GetAtomicBondVector ( )
inline

Definition at line 177 of file G4CrystalExtension.hh.

177{return theAtomicBondVector;};

References theAtomicBondVector.

◆ GetAtomPos() [1/3]

G4bool G4CrystalExtension::GetAtomPos ( const G4Element anEl,
std::vector< G4ThreeVector > &  vecout 
)

Definition at line 133 of file G4CrystalExtension.cc.

133 {
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}
static const G4double pos
G4bool FillAtomicPos(G4ThreeVector &pos, std::vector< G4ThreeVector > &vecout)

References G4CrystalUnitCell::FillAtomicPos(), GetAtomBase(), pos, and theUnitCell.

Referenced by GetAtomPos().

◆ GetAtomPos() [2/3]

G4bool G4CrystalExtension::GetAtomPos ( G4int  anElIdx,
std::vector< G4ThreeVector > &  vecout 
)
inline

Definition at line 148 of file G4CrystalExtension.hh.

148 {
149 GetAtomPos(fMaterial->GetElement(anElIdx),vecout);
150 return true;
151 }
G4bool GetAtomPos(const G4Element *anEl, std::vector< G4ThreeVector > &vecout)

References fMaterial, GetAtomPos(), and G4Material::GetElement().

◆ GetAtomPos() [3/3]

G4bool G4CrystalExtension::GetAtomPos ( std::vector< G4ThreeVector > &  vecout)

Definition at line 145 of file G4CrystalExtension.cc.

145 {
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}

References fMaterial, GetAtomPos(), G4Material::GetElementVector(), and pos.

◆ GetCijkl()

G4double G4CrystalExtension::GetCijkl ( G4int  i,
G4int  j,
G4int  k,
G4int  l 
) const
inline

Definition at line 108 of file G4CrystalExtension.hh.

108 {
109 return fElasticity[i][j][k][l];
110 }

References fElasticity.

◆ GetCpq()

G4double G4CrystalExtension::GetCpq ( G4int  p,
G4int  q 
) const
inline

Definition at line 116 of file G4CrystalExtension.hh.

116{ return fElReduced[p-1][q-1]; }
ReducedElasticity fElReduced

References fElReduced.

◆ GetElasticity()

const Elasticity & G4CrystalExtension::GetElasticity ( ) const
inline

Definition at line 104 of file G4CrystalExtension.hh.

104{ return fElasticity; }

References fElasticity.

◆ GetElReduced()

const ReducedElasticity & G4CrystalExtension::GetElReduced ( ) const
inline

Definition at line 105 of file G4CrystalExtension.hh.

105{ return fElReduced; }

References fElReduced.

◆ GetHash()

const std::size_t & G4VMaterialExtension::GetHash ( ) const
inlineinherited

Definition at line 71 of file G4VMaterialExtension.hh.

71{ return fHash;}

References G4VMaterialExtension::fHash.

◆ GetMaterial()

G4Material * G4CrystalExtension::GetMaterial ( )
inline

Definition at line 75 of file G4CrystalExtension.hh.

75{return fMaterial;};

References fMaterial.

◆ GetName()

const G4String & G4VMaterialExtension::GetName ( ) const
inlineinherited

Definition at line 73 of file G4VMaterialExtension.hh.

73{ return fName; }

References G4VMaterialExtension::fName.

◆ GetUnitCell()

G4CrystalUnitCell * G4CrystalExtension::GetUnitCell ( ) const
inline

Definition at line 88 of file G4CrystalExtension.hh.

89 {return theUnitCell;}

References theUnitCell.

Referenced by G4LogicalCrystalVolume::GetBasis().

◆ Print()

void G4CrystalExtension::Print ( ) const
inlinevirtual

Implements G4VMaterialExtension.

Definition at line 70 of file G4CrystalExtension.hh.

70{;};

◆ SetCpq()

void G4CrystalExtension::SetCpq ( G4int  p,
G4int  q,
G4double  value 
)

Definition at line 114 of file G4CrystalExtension.cc.

114 {
115 if (p>0 && p<7 && q>0 && q<7) fElReduced[p-1][q-1] = value;
116}

References fElReduced.

◆ SetElReduced()

void G4CrystalExtension::SetElReduced ( const ReducedElasticity mat)

Definition at line 104 of file G4CrystalExtension.cc.

104 {
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}

References fElReduced.

◆ SetMaterial()

void G4CrystalExtension::SetMaterial ( G4Material aMat)
inline

Definition at line 76 of file G4CrystalExtension.hh.

76{fMaterial = aMat;};

References fMaterial.

◆ SetUnitCell()

void G4CrystalExtension::SetUnitCell ( G4CrystalUnitCell aUC)
inline

Definition at line 87 of file G4CrystalExtension.hh.

87{theUnitCell=aUC;}

References theUnitCell.

Field Documentation

◆ fElasticity

Elasticity G4CrystalExtension::fElasticity
protected

Definition at line 100 of file G4CrystalExtension.hh.

Referenced by GetCijkl(), and GetElasticity().

◆ fElReduced

ReducedElasticity G4CrystalExtension::fElReduced
protected

Definition at line 101 of file G4CrystalExtension.hh.

Referenced by GetCpq(), GetElReduced(), SetCpq(), and SetElReduced().

◆ fHash

const std::size_t G4VMaterialExtension::fHash
protectedinherited

Definition at line 79 of file G4VMaterialExtension.hh.

Referenced by G4VMaterialExtension::GetHash().

◆ fMaterial

G4Material* G4CrystalExtension::fMaterial
private

◆ fName

const G4String& G4VMaterialExtension::fName
protectedinherited

Definition at line 76 of file G4VMaterialExtension.hh.

Referenced by G4VMaterialExtension::GetName().

◆ theAtomicBondVector

std::vector<G4AtomicBond*> G4CrystalExtension::theAtomicBondVector
private

Definition at line 172 of file G4CrystalExtension.hh.

Referenced by AddAtomicBond(), GetAtomicBond(), and GetAtomicBondVector().

◆ theCrystalAtomBaseMap

std::map<const G4Element*,G4CrystalAtomBase*> G4CrystalExtension::theCrystalAtomBaseMap
private

Definition at line 123 of file G4CrystalExtension.hh.

Referenced by AddAtomBase(), and GetAtomBase().

◆ theUnitCell

G4CrystalUnitCell* G4CrystalExtension::theUnitCell
private

Definition at line 85 of file G4CrystalExtension.hh.

Referenced by GetAtomPos(), GetUnitCell(), and SetUnitCell().


The documentation for this class was generated from the following files: