Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CrystalUnitCell.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// G4CrystalUnitCell
29//
30// Class Description:
31//
32
33#ifndef G4CrystalUnitCell_H
34#define G4CrystalUnitCell_H 1
35
36#include "G4ThreeVector.hh"
37#include "globals.hh"
38
41
42#include <vector>
43
45{
46 public:
47 G4CrystalUnitCell(G4double sizeA, G4double sizeB, G4double sizeC, G4double alpha, G4double beta,
48 G4double gamma, G4int spacegroup);
49
50 virtual ~G4CrystalUnitCell() = default;
51
52 inline G4int GetSpaceGroup() const { return theSpaceGroup; };
53 inline void SetSpaceGroup(G4int aInt) { theSpaceGroup = aInt; };
54
57
58 const G4ThreeVector& GetBasis(G4int idx) const;
59 const G4ThreeVector& GetUnitBasis(G4int idx) const;
60 inline G4ThreeVector GetSize() const { return theSize; }
61 inline G4ThreeVector GetAngle() const { return theAngle; }
62
63 // return theUnitBase[2] vector
65
66 const G4ThreeVector& GetRecBasis(G4int idx) const;
67 const G4ThreeVector& GetRecUnitBasis(G4int idx) const;
68 inline G4ThreeVector GetRecSize() const { return theRecSize; }
69 inline G4ThreeVector GetRecAngle() const { return theRecAngle; }
70
71 // Methods to populate atom position in the lattice from the basis
72 // and the unit basis
73 G4bool FillAtomicUnitPos(G4ThreeVector& pos, std::vector<G4ThreeVector>& vecout);
74 G4bool FillAtomicPos(G4ThreeVector& pos, std::vector<G4ThreeVector>& vecout);
75
76 // Methods to populate elasticity and reduced elasticity tensors
77 G4bool FillElReduced(G4double Cij[6][6]);
78
79 // The volumes of the cell
80 G4double ComputeCellVolume(); // compute and store the volume
81
82 inline G4double GetVolume() const { return theVolume; } // get the stored volume
83 inline G4double GetRecVolume() const { return theRecVolume; } // get the stored volume
84
85 // Squared Reciprocal and direct interplanar spacing
87 G4int l); // squared interplanar spacing
88
90 G4int l); // squared reciprocal interplanar spacing
91
92 G4double GetIntCosAng(G4int h1, G4int k1, G4int l1, G4int h2, G4int k2,
93 G4int l2); // cosine of the angle between two planes
94
95 protected:
96 G4ThreeVector nullVec; // Size and angles of the crystalline unit cell
97
98 G4ThreeVector theSize; // cell sizes
99 G4ThreeVector theAngle; // cell angles
100 G4ThreeVector theUnitBasis[3]; // Basis unit vectors in direct orientation
101 G4ThreeVector theBasis[3]; // Basis vectors in direct orientation
102 //
103 // Reciprocal size and angles of the crystalline unit cell
104 //
105 G4ThreeVector theRecSize; // reciprocal cell sizes
106 G4ThreeVector theRecAngle; // reciprocal cell angles
107 G4ThreeVector theRecUnitBasis[3]; // Basis unit vectors in reciprocal orientation
108 G4ThreeVector theRecBasis[3]; // Basis vectors in reciprocal orientation
109
110 private:
113 G4bool FillAmorphous(G4double Cij[6][6]) const;
114 G4bool FillCubic(G4double Cij[6][6]) const;
115 G4bool FillTetragonal(G4double Cij[6][6]) const;
116 G4bool FillOrthorhombic(G4double Cij[6][6]) const;
117 G4bool FillRhombohedral(G4double Cij[6][6]) const;
118 G4bool FillMonoclinic(G4double Cij[6][6]) const;
119 G4bool FillTriclinic(G4double Cij[6][6]) const;
120 G4bool FillHexagonal(G4double Cij[6][6]) const;
121
122 G4bool ReflectElReduced(G4double Cij[6][6]) const;
123
124 private:
125 G4int theSpaceGroup; //
126
127 G4double cosa, cosb, cosg;
128 G4double sina, sinb, sing;
129 G4double cosar, cosbr, cosgr;
130
131 G4double theVolume; // the cell volume
132 G4double theRecVolume; // the cell volume
133};
134
135#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4int GetSpaceGroup() const
G4ThreeVector theBasis[3]
const G4ThreeVector & GetRecUnitBasis(G4int idx) const
G4ThreeVector GetUnitBasisTrigonal()
const G4ThreeVector & GetRecBasis(G4int idx) const
G4double GetIntSp2(G4int h, G4int k, G4int l)
G4bool FillElReduced(G4double Cij[6][6])
virtual ~G4CrystalUnitCell()=default
G4ThreeVector GetRecAngle() const
G4CrystalUnitCell(G4double sizeA, G4double sizeB, G4double sizeC, G4double alpha, G4double beta, G4double gamma, G4int spacegroup)
G4bool FillAtomicPos(G4ThreeVector &pos, std::vector< G4ThreeVector > &vecout)
const G4ThreeVector & GetUnitBasis(G4int idx) const
G4ThreeVector theRecUnitBasis[3]
G4double GetRecIntSp2(G4int h, G4int k, G4int l)
G4double GetVolume() const
G4ThreeVector theRecBasis[3]
G4bool FillAtomicUnitPos(G4ThreeVector &pos, std::vector< G4ThreeVector > &vecout)
G4ThreeVector GetRecSize() const
G4double GetRecVolume() const
G4ThreeVector GetSize() const
G4ThreeVector theUnitBasis[3]
theBravaisLatticeType GetBravaisLattice()
G4ThreeVector theRecAngle
void SetSpaceGroup(G4int aInt)
theLatticeSystemType GetLatticeSystem()
G4ThreeVector GetAngle() const
const G4ThreeVector & GetBasis(G4int idx) const
G4double GetIntCosAng(G4int h1, G4int k1, G4int l1, G4int h2, G4int k2, G4int l2)