oILAB
Loading...
Searching...
No Matches
Gb.cpp
Go to the documentation of this file.
1//
2// Created by Nikhil Chandra Admal on 11/5/22.
3//
4#ifndef OILAB_GB_CPP
5#define OILAB_GB_CPP
6
7#include "../../include/Lattices/Gb.h"
8
9namespace oILAB {
10template <int dim>
13 : /* init */ bc(bc)
14 /* init */,
15 nA(&n.lattice == &(bc.A) ? n
16 : bc.getReciprocalLatticeDirectionInA(
17 -1 * n.reciprocalLatticeVector()))
18 /* init */,
19 nB(&n.lattice == &(bc.B) ? n
20 : bc.getReciprocalLatticeDirectionInB(
21 -1 * n.reciprocalLatticeVector()))
22 /* init */,
23 basisT(getBasisT(bc, n))
24 /* init */,
25 T(Lattice<dim>(bc.dscl.latticeBasis *
26 getBasisT(bc, n).template cast<double>())) {
27 if (&n.lattice != &(bc.A) && &n.lattice != &(bc.B))
28 throw std::runtime_error(
29 "The normal does not belong to the reciprocal lattices of A and B \n");
30 }
31 catch(std::runtime_error& e)
32 {
33 std::cout << e.what() << std::endl;
34 throw(std::runtime_error("GB construction failed. "));
35 }
36
37 template<int dim>
39 {
40 ReciprocalLatticeDirection<dim> dir= bc.getReciprocalLatticeDirectionInC(nA.reciprocalLatticeVector());
41 double cslPlaneSpacing= dir.planeSpacing();
42 double step= bc.shiftTensorA(d).cartesian().dot(nA.cartesian().normalized());
43 return std::remainder(step,cslPlaneSpacing);
44 }
45
46 template<int dim>
48 {
49 ReciprocalLatticeDirection<dim> dir= bc.getReciprocalLatticeDirectionInC(nA.reciprocalLatticeVector());
50 double cslPlaneSpacing= dir.planeSpacing();
51 double step= (bc.shiftTensorA(d)-bc.shiftTensorB(d)).cartesian().dot(nA.cartesian().normalized())/2.0;
52 return std::remainder(step,cslPlaneSpacing);
53 }
54
55 template<int dim>
57 {
58 ReciprocalLatticeDirection<dim> dir= bc.getReciprocalLatticeDirectionInC(nA.reciprocalLatticeVector());
59 double cslPlaneSpacing= dir.planeSpacing();
60 double step= bc.shiftTensorB(d).cartesian().dot(nB.cartesian().normalized());
61 return std::remainder(step,cslPlaneSpacing);
62 }
63
64 template<int dim> template<int dm>
65 typename std::enable_if<dm==2 || dm==3,std::vector<LatticeVector<dim>>>::type
66 Gb<dim>::box(std::vector<LatticeVector<dim>>& boxVectors,
67 const double& orthogonality,
68 const int& dsclFactor,
69 std::string filename,
70 bool orient) const
71 {
72 for (auto iter= std::next(boxVectors.begin()); iter < boxVectors.end(); iter++)
73 assert((*iter).dot(bc.getReciprocalLatticeDirectionInC(nA.reciprocalLatticeVector())) == 0 &&
74 "Box vectors not parallel to the grain boundary.");
75
76 auto config= bc.box(boxVectors,orthogonality,dsclFactor);
77 std::vector<LatticeVector<dim>> configuration;
78 for (const LatticeVector<dim>& latticeVector : config)
79 {
80 if (&(latticeVector.lattice) == &bc.A && latticeVector.dot(nA)<=0)
81 configuration.push_back(latticeVector);
82 if (&(latticeVector.lattice) == &bc.B && latticeVector.dot(nB)<=0)
83 configuration.push_back(latticeVector);
84 if (&(latticeVector.lattice) == &bc.csl || &(latticeVector.lattice) == &bc.dscl)
85 configuration.push_back(latticeVector);
86 }
87
88 // form the rotation matrix used to orient the system
89 MatrixDimD rotation= Eigen::Matrix<double,dim,dim>::Identity();;
90 Eigen::Matrix<double,dim,dim-1> orthogonalVectors;
91 if (orient) {
92 if (dim==3) {
93 orthogonalVectors.col(0) = boxVectors[1].cartesian().normalized();
94 orthogonalVectors.col(1) = boxVectors[2].cartesian().normalized();
95 }
96 else if (dim==2)
97 orthogonalVectors.col(0)= boxVectors[1].cartesian().normalized();
98
99 rotation=Rotation<dim>(orthogonalVectors);
100 }
101 //assert((rotation*rotation.transpose()).template isApprox(Eigen::Matrix<double,dim,dim>::Identity())
102 assert((rotation*rotation.transpose()).isApprox(Eigen::Matrix<double,dim,dim>::Identity())
103 && "Cannot orient the grain boundary. The GB plane box vectors are not orthogonal.");
104
105 if(!filename.empty()) {
106 std::ofstream file;
107 file.open(filename);
108 if (!file) std::cerr << "Unable to open file";
109 file << configuration.size() << std::endl;
110 file << "Lattice=\"";
111
112 LatticeVector<dim> origin(-1*boxVectors[0]);
113 if (dim == 2) {
114 file << (rotation * 2 * boxVectors[0].cartesian()).transpose() << " 0 ";
115 file << (rotation * boxVectors[1].cartesian()).transpose() << " 0 ";
116 file << " 0 0 1 ";
117 file << "\" Properties=atom_types:I:1:pos:R:3:radius:R:1 PBC=\"F T T\" origin=\"";
118 file << (rotation * origin.cartesian()).transpose() << " 0.0\"" << std::endl;
119 for (const auto &vector: configuration)
120 if (&(vector.lattice) == &bc.A)
121 file << 1 << " " << (rotation * vector.cartesian()).transpose() << " " << 0.0 << " "
122 << 0.05 << std::endl;
123 else if (&(vector.lattice) == &bc.B)
124 file << 2 << " " << (rotation * vector.cartesian()).transpose() << " " << 0.0 << " "
125 << 0.05 << std::endl;
126 else if (&(vector.lattice) == &bc.csl)
127 file << 3 << " " << (rotation * vector.cartesian()).transpose() << " " << 0.0 << " " << 0.2
128 << std::endl;
129 else
130 file << 4 << " " << (rotation * vector.cartesian()).transpose() << " " << 0.0 << " "
131 << 0.01 << std::endl;
132 } else if (dim == 3) {
133 file << (rotation * 2 * boxVectors[0].cartesian()).transpose() << " ";
134 file << (rotation * boxVectors[1].cartesian()).transpose() << " ";
135 file << (rotation * boxVectors[2].cartesian()).transpose() << " ";
136 file << "\" Properties=atom_types:I:1:pos:R:3:radius:R:1 PBC=\"F T T\" origin=\"";
137 file << (rotation * origin.cartesian()).transpose() << "\"" << std::endl;
138
139 for (const auto &vector: configuration)
140 if (&(vector.lattice) == &bc.A)
141 file << 1 << " " << (rotation * vector.cartesian()).transpose() << " " << 0.05
142 << std::endl;
143 else if (&(vector.lattice) == &bc.B)
144 file << 2 << " " << (rotation * vector.cartesian()).transpose() << " " << 0.05
145 << std::endl;
146 else if (&(vector.lattice) == &bc.csl)
147 file << 3 << " " << (rotation * vector.cartesian()).transpose() << " " << 0.2 << std::endl;
148 else
149 file << 4 << " " << (rotation * vector.cartesian()).transpose() << " " << 0.01
150 << std::endl;
151 }
152 file.close();
153 }
154 return configuration;
155 }
156
157 template<int dim> template<int dm>
158 typename std::enable_if<dm==2,LatticeVector<dim>>::type
160 {
161 LatticeVector<dm> axisAxnA(nA.reciprocalLatticeVector().cross().latticeVector());
162 return (LatticeVector<dm>(bc.getLatticeDirectionInC(axisAxnA).latticeVector()));
163 }
164
165 template<int dim> template<int dm>
166 typename std::enable_if<dm==3,LatticeVector<dim>>::type
168 {
169 assert(abs(nA.cartesian().dot(axis.cartesian())) < FLT_EPSILON);
170 ReciprocalLatticeDirection<dm> axisA(bc.A.reciprocalLatticeVector(axis.cartesian()));
171 LatticeVector<dm> axisAxnA(axisA.reciprocalLatticeVector().cross(nA.reciprocalLatticeVector()).latticeVector());
172 return (LatticeVector<dm>(bc.getLatticeDirectionInC(axisAxnA).latticeVector()));
173 }
174
175
176 template<int dim>
178 {
179 MatrixDimI output;
180 MatrixDimI Q = 2*bc.LambdaA-Eigen::Matrix<IntScalarType,dim,dim>::Identity();
181 auto nC= bc.getReciprocalLatticeDirectionInC(this->nB.reciprocalLatticeVector());
182
185
186 IntScalarType detMN= (bc.M*bc.N).template cast<double>().determinant();
187 Eigen::DiagonalMatrix<IntScalarType,dim> nCDiagonal(nC.reciprocalLatticeVector());
188 RationalMatrix<dim> temp(adjM*adjN*Q*nCDiagonal, Eigen::Matrix<IntScalarType,dim,dim>::Constant(2*detMN));
189
190 IntScalarType alpha= IntegerMath<IntScalarType>::gcd(temp.integerMatrix.diagonal().eval());
191 IntScalarType beta= temp.mu;
193
194 assert(IntegerMath<IntScalarType>::gcd(alpha,beta) == 1);
195
196 //auto basis= bc.dscl.planeParallelLatticeBasis(m,true);
197 auto basis= bc.dscl.planeParallelLatticeBasis(m,false);
198 for(int i=0; i<dim; ++i) {
199 output.col(i) = basis[i].latticeVector();
200 if (i==0) output.col(i)= beta*output.col(i);
201 }
202 return output;
203 }
204
205 template<int dim>
207 {
209 IntScalarType det= basisT.template cast<double>().determinant();
210 // basisTA
211 // T --------> dscl
212 if (&(v.lattice) == &(bc.csl) && IntegerMath<IntScalarType>::gcd(v)%2 == 0)
213 {
214 VectorDimI integerCoordinates= adj * bc.getLatticeVectorInD(v);
215 assert(IntegerMath<IntScalarType>::gcd(integerCoordinates) % det == 0);
216 integerCoordinates= integerCoordinates/det;
217
218 auto temp= LatticeVector<dim>(integerCoordinates,T);
219 if (temp.cartesian().dot(v.cartesian()) < 0) temp= -1*temp;
220 return temp;
221 }
222 else
223 throw(std::runtime_error("The input lattice vector should belong to twice the CSL lattice"));
224
225 }
226
227 template<int dim>
229 {
230 // basisTA
231 // T --------> dscl
232 if (&(v.lattice) == &(bc.dscl))
233 return ReciprocalLatticeVector<dim>((basisT.transpose()*v).eval(),T);
234 else
235 throw(std::runtime_error("The input reciprocal lattice vector should belong to the reciprocal lattice of the DSCL"));
236
237 }
238
239 template<int dim>
244
245 template class Gb<2>;
247 template std::vector<LatticeVector<2>> Gb<2>::box<2>(std::vector<LatticeVector<2>>& boxVectors,
248 const double& orthogonality,
249 const int& dsclFactor,
250 std::string filename,
251 bool orient) const;
252
253 template class Gb<3>;
255 template std::vector<LatticeVector<3>> Gb<3>::box<3>(std::vector<LatticeVector<3>>& boxVectors,
256 const double& orthogonality,
257 const int& dsclFactor,
258 std::string filename,
259 bool orient) const;
260
261 template class Gb<4>;
262 template class Gb<5>;
263 } // namespace oILAB
264
265#endif //OILAB_GB_CPP
const MatrixDimI N
Integer matrix that connects the bases and of lattices and the CSL , respectively: .
Definition BiCrystal.h:63
const MatrixDimI LambdaA
Shift tensor describes the shift in the CSL when lattice is shifted by a DSCL vector.
Definition BiCrystal.h:119
ReciprocalLatticeDirection< dim > getReciprocalLatticeDirectionInC(const ReciprocalLatticeVector< dim > &v) const
const MatrixDimI M
Integer matrix that connects the bases and of lattices and the CSL , respectively: .
Definition BiCrystal.h:58
const Lattice< dim > dscl
DCSL lattice .
Definition BiCrystal.h:106
Definition Gb.h:13
std::enable_if< dm==2, LatticeVector< dim > >::type getPeriodVector(const ReciprocalLatticeVector< dim > &axis) const
Definition Gb.cpp:159
ReciprocalLatticeDirection< dim > getReciprocalLatticeDirectionInT(const ReciprocalLatticeVector< dim > &v) const
Definition Gb.cpp:240
LatticeVector< dim > getLatticeVectorInT(const LatticeVector< dim > &v) const
Definition Gb.cpp:206
double stepHeightB(const LatticeVector< dim > &d) const
Computes the step height of a disconnection formed by displacing lattice by a Burgers vector .
Definition Gb.cpp:56
double stepHeight(const LatticeVector< dim > &d) const
Computes the step height of a disconnection formed by displacing lattice by and by
Definition Gb.cpp:47
Gb(const BiCrystal< dim > &bc, const ReciprocalLatticeDirection< dim > &n)
Constructs a grain boundary of a given orientation in a bicrystal.
Definition Gb.cpp:11
MatrixDimI getBasisT(const BiCrystal< dim > &bc, const ReciprocalLatticeDirection< dim > &n)
Definition Gb.cpp:177
typename LatticeCore< dim >::MatrixDimI MatrixDimI
Definition Gb.h:17
typename LatticeCore< dim >::IntScalarType IntScalarType
Definition Gb.h:18
typename LatticeCore< dim >::VectorDimI VectorDimI
Definition Gb.h:14
const BiCrystal< dim > & bc
Definition Gb.h:28
ReciprocalLatticeVector< dim > getReciprocalLatticeVectorInT(const ReciprocalLatticeVector< dim > &v) const
Definition Gb.cpp:228
std::enable_if< dm==2||dm==3, std::vector< LatticeVector< dim > > >::type box(std::vector< LatticeVector< dim > > &boxVectors, const double &orthogonality, const int &dsclFactor, std::string filename="", bool orient=false) const
Definition Gb.cpp:66
double stepHeightA(const LatticeVector< dim > &d) const
Computes the step height of a disconnection formed by displacing lattice by a Burgers vector .
Definition Gb.cpp:38
typename LatticeCore< dim >::MatrixDimD MatrixDimD
Definition Gb.h:16
Lattice class.
Definition Lattice.h:31
LatticeVector class.
const Lattice< dim > & lattice
VectorDimD cartesian() const
static Eigen::Matrix< IntScalarType, dim, dim > adjoint(const Eigen::Matrix< IntScalarType, dim, dim > &input)
const MatrixDimI & integerMatrix
const IntScalarType & mu
static IntScalarType gcd(const IntScalarType &a, const IntScalarType &b)
Definition IntegerMath.h:26
const ReciprocalLatticeVector< dim > & reciprocalLatticeVector() const
Returns a constant reference to the base class (ReciprocalLatticeVector)
double planeSpacing() const
Returns the spacing between two consecutive lattice planes.