File indexing completed on 2025-01-18 09:58:05
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033 #ifndef G4CrystalUnitCell_H
0034 #define G4CrystalUnitCell_H 1
0035
0036 #include "G4ThreeVector.hh"
0037 #include "globals.hh"
0038
0039 #include "G4CrystalBravaisLattices.h"
0040 #include "G4CrystalLatticeSystems.h"
0041
0042 #include <vector>
0043
0044 class G4CrystalUnitCell
0045 {
0046 public:
0047 G4CrystalUnitCell(G4double sizeA, G4double sizeB, G4double sizeC, G4double alpha, G4double beta,
0048 G4double gamma, G4int spacegroup);
0049
0050 virtual ~G4CrystalUnitCell() = default;
0051
0052 inline G4int GetSpaceGroup() const { return theSpaceGroup; };
0053 inline void SetSpaceGroup(G4int aInt) { theSpaceGroup = aInt; };
0054
0055 theLatticeSystemType GetLatticeSystem() { return GetLatticeSystem(theSpaceGroup); }
0056 theBravaisLatticeType GetBravaisLattice() { return GetBravaisLattice(theSpaceGroup); }
0057
0058 const G4ThreeVector& GetBasis(G4int idx) const;
0059 const G4ThreeVector& GetUnitBasis(G4int idx) const;
0060 inline G4ThreeVector GetSize() const { return theSize; }
0061 inline G4ThreeVector GetAngle() const { return theAngle; }
0062
0063
0064 G4ThreeVector GetUnitBasisTrigonal();
0065
0066 const G4ThreeVector& GetRecBasis(G4int idx) const;
0067 const G4ThreeVector& GetRecUnitBasis(G4int idx) const;
0068 inline G4ThreeVector GetRecSize() const { return theRecSize; }
0069 inline G4ThreeVector GetRecAngle() const { return theRecAngle; }
0070
0071
0072
0073 G4bool FillAtomicUnitPos(G4ThreeVector& pos, std::vector<G4ThreeVector>& vecout);
0074 G4bool FillAtomicPos(G4ThreeVector& pos, std::vector<G4ThreeVector>& vecout);
0075
0076
0077 G4bool FillElReduced(G4double Cij[6][6]);
0078
0079
0080 G4double ComputeCellVolume();
0081
0082 inline G4double GetVolume() const { return theVolume; }
0083 inline G4double GetRecVolume() const { return theRecVolume; }
0084
0085
0086 G4double GetIntSp2(G4int h, G4int k,
0087 G4int l);
0088
0089 G4double GetRecIntSp2(G4int h, G4int k,
0090 G4int l);
0091
0092 G4double GetIntCosAng(G4int h1, G4int k1, G4int l1, G4int h2, G4int k2,
0093 G4int l2);
0094
0095 protected:
0096 G4ThreeVector nullVec;
0097
0098 G4ThreeVector theSize;
0099 G4ThreeVector theAngle;
0100 G4ThreeVector theUnitBasis[3];
0101 G4ThreeVector theBasis[3];
0102
0103
0104
0105 G4ThreeVector theRecSize;
0106 G4ThreeVector theRecAngle;
0107 G4ThreeVector theRecUnitBasis[3];
0108 G4ThreeVector theRecBasis[3];
0109
0110 private:
0111 theLatticeSystemType GetLatticeSystem(G4int aGroup);
0112 theBravaisLatticeType GetBravaisLattice(G4int aGroup);
0113 G4bool FillAmorphous(G4double Cij[6][6]) const;
0114 G4bool FillCubic(G4double Cij[6][6]) const;
0115 G4bool FillTetragonal(G4double Cij[6][6]) const;
0116 G4bool FillOrthorhombic(G4double Cij[6][6]) const;
0117 G4bool FillRhombohedral(G4double Cij[6][6]) const;
0118 G4bool FillMonoclinic(G4double Cij[6][6]) const;
0119 G4bool FillTriclinic(G4double Cij[6][6]) const;
0120 G4bool FillHexagonal(G4double Cij[6][6]) const;
0121
0122 G4bool ReflectElReduced(G4double Cij[6][6]) const;
0123
0124 private:
0125 G4int theSpaceGroup;
0126
0127 G4double cosa, cosb, cosg;
0128 G4double sina, sinb, sing;
0129 G4double cosar, cosbr, cosgr;
0130
0131 G4double theVolume;
0132 G4double theRecVolume;
0133 };
0134
0135 #endif