oILAB
Loading...
Searching...
No Matches
Lattice.cpp
Go to the documentation of this file.
1/* This file is part of gbLAB.
2 *
3 * gbLAB is distributed without any warranty under the MIT License.
4 */
5
6
7#ifndef gbLAB_Lattice_cpp_
8#define gbLAB_Lattice_cpp_
9
10#include <Eigen/Eigenvalues>
11
12#include "../../include/Lattices/LatticeModule.h"
13#include "../../include/Math/GramMatrix.h"
14#include <iomanip>
15
16namespace oILAB {
17
18// /**********************************************************************/
19// template <int dim>
20// typename Lattice<dim>::MatrixDimD Lattice<dim>::getLatticeBasis(const MatrixDimD &A, const MatrixDimD &Q)
21// {
22//
23// // Check that Q is orthogonal
24// const MatrixDimD QQT(Q * Q.transpose());
25// const double QQTerror((QQT - MatrixDimD::Identity()).norm());
26// if (QQTerror > 10.0 * DBL_EPSILON * dim * dim)
27// {
28// throw std::runtime_error("The rotation matrix Q is not orthogonal: norm(Q*Q^T-I)="+std::to_string(QQTerror)+"\n");
29// }
30//
31// const double Qdet(Q.determinant());
32// if (std::fabs(Q.determinant() - 1.0)>FLT_EPSILON)
33// {
34// throw std::runtime_error("The rotation matrix is not proper: det(Q)="+std::to_string(Qdet)+"\n");
35// }
36//
37// return Q * A;
38// }
39
40
41 /**********************************************************************/
42 /**********************************************************************/
43 template <int dim>
45 /* init */ latticeBasis(Fin*A)
46 /* init */,reciprocalBasis(latticeBasis.inverse().transpose())
47 /* init */,F(Fin)
48 {
49
50 }
51
52 /**********************************************************************/
53 template <int dim>
55 {
56 RLLL rlll(latticeBasis,0.75);
57 auto structureMatrix= rlll.reducedBasis();
58 auto U= rlll.unimodularMatrix();
59 Lattice<dim> reducedLattice(structureMatrix);
60 VectorDimD nd;
61 nd= reducedLattice.reciprocalBasis.transpose()*d;
62 LatticeVector<dim> tempReduced(LatticeCore<dim>::rationalApproximation(nd),reducedLattice);
63 VectorDimI tempRecovered= U*tempReduced;
64 LatticeVector<dim> temp(tempRecovered, *this);
65
66 const GramMatrix<double,2> G(std::array<VectorDimD,2>{temp.cartesian().normalized(),d.normalized()});
67 const double crossNorm(sqrt(G.determinant()));
68 if(crossNorm>tol)
69 {
70 std::cout<<"input direction="<<d.normalized().transpose()<<std::endl;
71 std::cout<<"lattice direction="<<temp.cartesian().normalized().transpose()<<std::endl;
72 std::cout<<"cross product norm="<<std::setprecision(15)<<std::scientific<<crossNorm<<std::endl;
73 std::cout<<"tolerance="<<std::setprecision(15)<<std::scientific<<tol<<std::endl;
74 throw std::runtime_error("LATTICE DIRECTION NOT FOUND\n");
75 }
76 return LatticeDirection<dim>(temp);
77 }
78
79 /**********************************************************************/
80 template <int dim>
82 {
83 RLLL rlll(latticeBasis,0.75);
84 auto structureMatrix= rlll.reducedBasis();
86 Lattice<dim> reducedLattice(structureMatrix);
87 VectorDimD nd;
88 nd= reducedLattice.latticeBasis.transpose()*d;
90 VectorDimI tempRecovered(MatrixDimIExt<IntScalarType,dim>::adjoint(U.matrix()).transpose()* tempReduced);
91 ReciprocalLatticeVector<dim> temp(tempRecovered, *this);
92
93 const GramMatrix<double,2> G(std::array<VectorDimD,2>{temp.cartesian().normalized(),d.normalized()});
94 const double crossNorm(sqrt(abs(G.determinant())));
95 if(crossNorm>tol)
96 {
97 std::cout<<"input direction="<<std::setprecision(15)<<std::scientific<<d.normalized().transpose()<<std::endl;
98 std::cout<<"reciprocal lattice direction="<<std::setprecision(15)<<std::scientific<<temp.cartesian().normalized().transpose()<<std::endl;
99 std::cout<<"cross product norm="<<std::setprecision(15)<<std::scientific<<crossNorm<<std::endl;
100 std::cout<<"tolerance="<<std::setprecision(15)<<std::scientific<<tol<<std::endl;
101 throw std::runtime_error("RECIPROCAL LATTICE DIRECTION NOT FOUND\n");
102 }
104 }
105
106 /**********************************************************************/
107 template <int dim>
109 const typename BestRationalApproximation::LongIntType& maxDen,
110 const double& magnitudeTol,
111 const double& directionTol) const
112 {
113 const LatticeDirection<dim> ld(latticeDirection(d,directionTol));
114 const BestRationalApproximation bra(d.norm() / ld.cartesian().norm(), maxDen);
115 const Rational rat(bra.num, bra.den);
116 const RationalLatticeDirection<dim> rld(rat, ld);
117 if ((rld.cartesian() - d).squaredNorm() > magnitudeTol)
118 {
119 std::cout << "input vector=" << d.transpose() << std::endl;
120 std::cout << "lattice direction=" << ld.cartesian().transpose() << std::endl;
121 std::cout << "rational=" << rat << std::endl;
122 std::cout << "d.norm()/ld.cartesian().norm()=" << d.norm() / ld.latticeVector().norm() << std::endl;
123 throw std::runtime_error("Rational Lattice DirectionType NOT FOUND\n");
124 }
125 return rld;
126 }
127
128 /**********************************************************************/
129
130 template <int dim>
132 const typename BestRationalApproximation::LongIntType &maxDen,
133 const double& magnitudeTol,
134 const double& directionTol) const
135 {
136 const ReciprocalLatticeDirection<dim> rld(reciprocalLatticeDirection(d,directionTol));
137 const BestRationalApproximation bra(d.norm() / rld.cartesian().norm(), maxDen);
138 const Rational rat(bra.num, bra.den);
139 const RationalReciprocalLatticeDirection<dim> rrld(rat, rld);
140 if ((rrld.cartesian() - d).squaredNorm() > magnitudeTol)
141 {
142 std::cout << "input reciprocal vector=" << d.transpose() << std::endl;
143 std::cout << "reciprocal lattice direction=" << rld.cartesian().transpose() << std::endl;
144 std::cout << "rational=" << rat << std::endl;
145 std::cout << "d.norm()/rld.cartesian().norm()=" << d.norm() / rld.reciprocalLatticeVector().norm() << std::endl;
146 throw std::runtime_error("Rational Reciprocal Lattice DirectionType NOT FOUND\n");
147 }
148 return rrld;
149 }
150
151 /**********************************************************************/
152 template <int dim>
154 {
155 return LatticeVector<dim>(p, *this);
156 }
157
158 /**********************************************************************/
159 template <int dim>
164
165 /**********************************************************************/
166 template<int dim>
168 const bool& useRLLL) const
169 {
170 assert(this == &l.lattice && "Vectors belong to different Lattices.");
172 auto matrix= IntegerMath<IntScalarType>::ccum(outOfPlaneVector);
173 std::vector<LatticeDirection<dim>> out;
174 int columnIndex= -1;
175 for(const auto& column : matrix.colwise())
176 {
177 columnIndex++;
178 if (columnIndex != 0)
179 out.push_back(LatticeDirection<dim>(column-column.dot(l.reciprocalLatticeVector())*matrix.col(0),*this));
180 else
181 out.push_back(LatticeDirection<dim>(column,*this));
182 }
183
184 if(!useRLLL) return out;
185
186 Eigen::MatrixXd planeParallelLatticeBasisCartesian(dim,dim-1);
187 int index= 0;
188 for(auto it=std::next(out.begin()); it!=out.end(); ++it)
189 {
190 planeParallelLatticeBasisCartesian.col(index)= (*it).cartesian();
191 index++;
192 }
193
194 if constexpr(dim>2){
195 const auto& rlll= RLLL(planeParallelLatticeBasisCartesian, 0.75);
196 planeParallelLatticeBasisCartesian= rlll.reducedBasis();
197 MatrixDimI outIntegerCoords;
198 for (int i=0; i<dim; ++i)
199 outIntegerCoords.col(i)= out[i].latticeVector();
200 Eigen::Matrix<IntScalarType,dim,dim-1> planeParallelLatticeBasisIntegerCoords= outIntegerCoords.block(0,1,dim,dim-1)*rlll.unimodularMatrix();
201 for (int i=1; i<dim; ++i)
202 out[i] = LatticeDirection<dim>(planeParallelLatticeBasisIntegerCoords.col(i-1),*this);
203 }
204 /*
205 if(dim>2) {
206 planeParallelLatticeBasisCartesian = RLLL(planeParallelLatticeBasisCartesian, 0.75).reducedBasis();
207 index = 1;
208 for (const auto &column: planeParallelLatticeBasisCartesian.colwise()) {
209 out[index] = this->latticeDirection(column);
210 index++;
211 }
212 }
213 */
214
215 // (A^T A)^{-1} A^T
216 Eigen::MatrixXd pseudoInverse(dim-1,dim);
217 pseudoInverse= (planeParallelLatticeBasisCartesian.transpose() * planeParallelLatticeBasisCartesian).inverse() *
218 (planeParallelLatticeBasisCartesian.transpose());
219
220 LatticeVector<dim> temp(*this);
221 for(int i=0;i<dim-1;i++)
222 {
223 temp= temp + round(out[0].cartesian().dot(pseudoInverse.row(i)))*out[i+1].latticeVector();
224 }
225 out[0]= LatticeDirection<dim>(out[0].latticeVector()-temp);
226
227
228 return out;
229 }
230
231 /**********************************************************************/
232 template<int dim>
233 std::vector<ReciprocalLatticeDirection<dim>> Lattice<dim>::directionOrthogonalReciprocalLatticeBasis(const LatticeDirection<dim>& l,
234 const bool& useRLLL) const
235 {
236 assert(this == &l.lattice && "Vectors belong to different Lattices.");
237 auto nonOrthogonalReciprocalVector= IntegerMath<IntScalarType>::solveBezout(l.latticeVector());
238 auto matrix= IntegerMath<IntScalarType>::ccum(nonOrthogonalReciprocalVector);
239 std::vector<ReciprocalLatticeDirection<dim>> out;
240 int columnIndex= -1;
241 for(const auto& column : matrix.colwise())
242 {
243 columnIndex++;
244 if (columnIndex != 0) {
245 VectorDimI temp= column - column.dot(l.latticeVector()) * matrix.col(0);
247 }
248 else {
249 VectorDimI temp = column;
251 }
252 }
253 if (!useRLLL) return out;
254
255 Eigen::MatrixXd directionOrthogonalReciprocalLatticeBasisCartesian(dim,dim-1);
256 int index= 0;
257 for(auto it=std::next(out.begin()); it!=out.end(); ++it)
258 {
259 directionOrthogonalReciprocalLatticeBasisCartesian.col(index)= (*it).cartesian();
260 index++;
261 }
262
263 if constexpr(dim>2) {
264 const auto &rlll = RLLL(directionOrthogonalReciprocalLatticeBasisCartesian, 0.75);
265 MatrixDimI outIntegerCoords;
266 for (int i = 0; i < dim; ++i)
267 outIntegerCoords.col(i) = out[i].reciprocalLatticeVector();
268 Eigen::Matrix<IntScalarType, dim, dim - 1> directionOrthogonalReciprocalLatticeBasisIntegerCoords =
269 outIntegerCoords.block(0, 1, dim, dim - 1) * rlll.unimodularMatrix();
270 for (int i = 1; i < dim; ++i)
272 directionOrthogonalReciprocalLatticeBasisIntegerCoords.col(i - 1), *this);
273 }
274
275 Eigen::MatrixXd pseudoInverse(dim-1,dim);
276 pseudoInverse= (directionOrthogonalReciprocalLatticeBasisCartesian.transpose() * directionOrthogonalReciprocalLatticeBasisCartesian).inverse() *
277 (directionOrthogonalReciprocalLatticeBasisCartesian.transpose());
278
280 for(int i=0;i<dim-1;i++)
281 {
282 temp= temp + round(out[0].cartesian().dot(pseudoInverse.row(i)))*out[i+1].reciprocalLatticeVector();
283 }
284 out[0]= ReciprocalLatticeDirection<dim>(out[0].reciprocalLatticeVector()-temp);
285 return out;
286
287 /*
288 directionOrthogonalReciprocalLatticeBasisCartesian= RLLL(directionOrthogonalReciprocalLatticeBasisCartesian,0.75).reducedBasis();
289 index= 1;
290 for(const auto& column : directionOrthogonalReciprocalLatticeBasisCartesian.colwise())
291 {
292 out[index]= this->reciprocalLatticeDirection(column);
293 index++;
294 }
295 return out;
296 */
297 }
298
299 /**********************************************************************/
300 template<int dim>
302 {
303 if(&(r.lattice) != this)
304 throw(std::runtime_error("The input reciprocal lattice vectors does not belong to the current lattice."));
305 return 1.0/r.cartesian().norm();
306 }
307
308 template<int dim> template<int dm>
309 typename std::enable_if<dm==3,std::vector<typename Lattice<dim>::MatrixDimD>>::type
310 Lattice<dim>::generateCoincidentLattices(const ReciprocalLatticeDirection<dim>& rd, const double& maxDen, const int& N) const
311 {
312 std::vector<MatrixDimD> output;
313 std::map<IntScalarType,MatrixDimD> temp;
314 auto basis= planeParallelLatticeBasis(rd);
315 double epsilon=1e-8;
316 IntScalarType keyScale= 1e6;
317
318 auto b1= basis[1].cartesian();
319 auto b2= basis[2].cartesian();
320
321 // Following the notation in Algorithm 3 of
322 // "Interface dislocations and grain boundary disconnections using Smith normal bicrystallography"
323 auto q1= b1;
324 std::vector<std::pair<int,int>> coPrimePairs= farey(N,false);
325 for(const auto& pair : coPrimePairs)
326 {
327 VectorDimI pIndices;
328 auto q2= pair.first*b1+pair.second*b2;
329
330 if (abs(q2.norm())<epsilon ) continue;
331 double ratio= q2.norm()/q1.norm();
332 BestRationalApproximation bra(ratio,maxDen);
333 double error= ratio-static_cast<double>(bra.num)/bra.den;
334 if (abs(error) > epsilon)
335 continue;
336 else
337 {
338 double cosTheta= b1.dot(q2)/(q1.norm()*q2.norm());
339 if (cosTheta-1>-epsilon) cosTheta= 1.0;
340 if (cosTheta+1<epsilon) cosTheta= -1.0;
341 double theta= acos(cosTheta);
342 MatrixDimD rotation;
343 rotation= Eigen::AngleAxis<double>(theta,rd.cartesian().normalized());
344 IntScalarType key= theta*keyScale;
345 temp.insert(std::pair<IntScalarType,MatrixDimD>(key,rotation));
346 }
347 }
348 std::transform(temp.begin(), temp.end(),
349 std::back_inserter(output),
350 [](const std::pair<IntScalarType,MatrixDimD>& p) {
351 return p.second;
352 });
353 return output;
354 }
355
356
357 template<int dim> template<int dm>
358 typename std::enable_if<dm==2,std::vector<typename Lattice<dim>::MatrixDimD>>::type
359 Lattice<dim>::generateCoincidentLattices(const double& maxStrain, const int& maxDen, const int& N) const
360 {
361 std::vector<MatrixDimD> output(generateCoincidentLattices(*this,maxStrain,maxDen,N));
362 return output;
363 }
364
365
366 template<int dim> template<int dm>
367 typename std::enable_if<dm==2,std::vector<typename Lattice<dim>::MatrixDimD>>::type
369 const double& maxStrain,
370 const int& maxDen,
371 const int& N) const
372 {
373 int numberOfConfigurations= 0;
374 const int maxConfigurations= 80;
375 std::vector<MatrixDimD> output;
376 std::map<IntScalarType,MatrixDimD> temp;
377 MatrixDimI mn(MatrixDimI::Zero());
378 MatrixDimI md(MatrixDimI::Ones());
379 std::vector<std::pair<int,int>> coPrimePairs= farey(N,false);
380 double epsilon=1e-8;
381 IntScalarType keyScale= 1e6;
382
383 for(const auto& pair1 : coPrimePairs)
384 {
385 VectorDimI pIndices;
386 pIndices << pair1.first, pair1.second;
387 LatticeVector<dm> q2(pIndices,undeformedLattice);
388 double ratio= q2.cartesian().norm() / latticeBasis.col(0).norm();
389
390
391 if (abs(maxStrain) < epsilon)
392 {
393 BestRationalApproximation alpha(ratio, maxDen);
394 double error = ratio - static_cast<double>(alpha.num) / alpha.den;
395 if (abs(error) > epsilon)
396 continue;
397 else
398 {
399 double cosTheta = latticeBasis.col(0).dot(q2.cartesian()) / (latticeBasis.col(0).norm() * q2.cartesian().norm());
400 if (cosTheta - 1 > -epsilon) cosTheta = 1.0;
401 if (cosTheta + 1 < epsilon) cosTheta = -1.0;
402 double theta = acos(cosTheta);
403 Eigen::Rotation2D<double> rotation(theta);
404 IntScalarType key= theta*keyScale;
405 temp.insert(std::pair<IntScalarType,MatrixDimD>(key,rotation.toRotationMatrix()));
406 }
407 }
408 else
409 {
410 RationalApproximations <IntScalarType> alphaSequence(ratio, maxDen, ratio*maxStrain);
411 for (const auto& alpha: alphaSequence.approximations)
412 {
413 RationalLatticeDirection<dm> q2ByAlpha(Rational<IntScalarType>(alpha.d, alpha.n), q2);
414
415 mn.col(0) = q2 * alpha.d;
416 md.col(0).setConstant(alpha.n);
417
418 double s1 =
419 (q2ByAlpha.cartesian().squaredNorm() - latticeBasis.col(0).squaredNorm()) /
420 latticeBasis.col(0).squaredNorm();
421 if (abs(s1) > maxStrain) continue;
422
423 for (const auto &pair2: coPrimePairs) {
424 VectorDimI qIndices;
425 qIndices << pair2.first, pair2.second;
426 LatticeVector<dm> r2(qIndices, undeformedLattice);
427 double ratio2= r2.cartesian().norm() / latticeBasis.col(1).norm();
428 RationalApproximations<IntScalarType> betaSequence(ratio2, maxDen,maxStrain*ratio2);
429 for(const auto& beta : betaSequence.approximations) {
430 RationalLatticeDirection<dm> r2ByBeta(Rational<IntScalarType>(beta.d, beta.n), r2);
431 double s2 = (r2ByBeta.cartesian().squaredNorm() - latticeBasis.col(1).squaredNorm()) /
432 latticeBasis.col(1).squaredNorm();
433 double s3 = (q2ByAlpha.cartesian().dot(r2ByBeta.cartesian()) -
434 latticeBasis.col(0).dot(latticeBasis.col(1))) /
435 (latticeBasis.col(0).norm() * latticeBasis.col(1).norm());
436 if (abs(s2) > maxStrain || abs(s3) > maxStrain) continue;
437
438 // calculate the deformation gradient
439 MatrixDimD F = q2ByAlpha.cartesian() * reciprocalBasis.col(0).transpose() +
440 r2ByBeta.cartesian() * reciprocalBasis.col(1).transpose();
441 if (F.determinant() < 0) continue;
442
443 output.push_back( F );
444 numberOfConfigurations++;
445 std::cout << numberOfConfigurations << std::endl;
446 if (numberOfConfigurations == maxConfigurations)
447 return output;
448 }
449 }
450 }
451 }
452 }
453 // sort the angles if rotations are being asked
454 if (abs(maxStrain) < epsilon)
455 std::transform(temp.begin(), temp.end(),
456 std::back_inserter(output),
457 [](const std::pair<IntScalarType,MatrixDimD>& p) {
458 return p.second;
459 });
460 return output;
461 }
462
463 template<int dim> template<int dm>
464 typename std::enable_if<dm==3,std::vector<LatticeVector<dim>>>::type
465 Lattice<dim>::box(const std::vector<LatticeVector<dim>>& boxVectors, const std::string& filename) const
466 {
467 for(const LatticeVector<dim>& boxVector : boxVectors)
468 {
469 assert(this == &boxVector.lattice && "Box vectors belong to different lattice.");
470 }
471
472 // Form the box lattice
473 MatrixDimD C;
474 for (int i=0; i<dim; ++i) {
475 C.col(i) = boxVectors[i].cartesian();
476 }
477 assert(C.determinant() != 0);
478 Lattice<dim> boxLattice(C);
479
480 auto planeParallelBasis= planeParallelLatticeBasis(boxVectors[1].cross(boxVectors[2]),true);
481 MatrixDimD A;
482 MatrixDimI mat;
483 for(int i=0; i<dim; ++i) {
484 A.col(i) = planeParallelBasis[i].cartesian();
485 mat.col(i)=planeParallelBasis[i].latticeVector();
486 }
487
488 // l - lattice with plane parallel basis, where the second and third basis vectors
489 // lie on the plane spanned by boxVectors[1] and boxVectors[2]
490 Lattice<dim> l(A);
491
492 // Form plane parallel vectors v1 and v2
493 // v1 is along boxVector[1]
495 assert(v1.latticeVector()(0)==0);
496 IntScalarType x,y;
498 LatticeVector<dim> temp(l);
499 temp << 0,y,x;
500 LatticeDirection<dim> v2(temp);
501
502 int areaRatio= round((boxLattice.latticeBasis.col(1).cross(boxLattice.latticeBasis.col(2))).norm()/
503 (l.latticeBasis.col(1).cross(l.latticeBasis.col(2))).norm());
504 int scale1= abs(IntegerMath<IntScalarType>::gcd(boxVectors[1]));
505 int scale2= round(areaRatio/scale1);
506 int scale0= round(abs(C.determinant()/latticeBasis.determinant()/areaRatio));
507
508 std::vector<LatticeVector<dim>> output;
509
510 // r0, r1, and r2 are reciprocal vectors perpendicular to areas spanned by boxVectors[0],
511 // boxVectors[1], and boxVectors[2]
512 ReciprocalLatticeVector<dim> r0_temp(*this), r1_temp(*this), r2_temp(*this);
513 r0_temp= (boxVectors[1].cross(boxVectors[2])).reciprocalLatticeVector();
514 if (r0_temp.dot(boxVectors[0]) < 0)
515 r0_temp=-1*r0_temp;
516 r1_temp= (boxVectors[2].cross(boxVectors[0])).reciprocalLatticeVector();
517 if (r1_temp.dot(boxVectors[1]) < 0)
518 r1_temp=-1*r1_temp;
519 r2_temp= (boxVectors[0].cross(boxVectors[1])).reciprocalLatticeVector();
520 if (r2_temp.dot(boxVectors[2]) < 0)
521 r2_temp=-1*r2_temp;
522 assert(r0_temp.dot(boxVectors[1])==0 && r0_temp.dot(boxVectors[2])==0);
523 assert(r1_temp.dot(boxVectors[0])==0 && r1_temp.dot(boxVectors[2])==0);
524 assert(r2_temp.dot(boxVectors[0])==0 && r2_temp.dot(boxVectors[1])==0);
525 ReciprocalLatticeDirection<dim> r0(r0_temp), r1(r1_temp), r2(r2_temp);
526
527 for(int i=0; i<scale0; ++i)
528 {
529 for(int j=0; j<scale1; ++j)
530 {
531 for(int k=0; k<scale2; ++k) {
532 LatticeVector<dim> vectorIn_l(l);
533 vectorIn_l << i, 0, 0;
534 vectorIn_l= vectorIn_l + j*v1.latticeVector()+k*v2.latticeVector();
535 //LatticeVector<dim> vector(latticeVector(vectorIn_l.cartesian()));
536 LatticeVector<dim> vector((mat*vectorIn_l).eval(),*this);
537 //LatticeDirection<dim> v1(MatrixDimIExt<IntScalarType,dim>::adjoint(mat)*boxVectors[1],l);
538
539 int c0= IntegerMath<IntScalarType>::positive_modulo(vector.dot(r0),boxVectors[0].dot(r0)) -
540 vector.dot(r0);
541 c0= c0/boxVectors[0].dot(r0);
542 int c1= IntegerMath<IntScalarType>::positive_modulo(vector.dot(r1),boxVectors[1].dot(r1)) -
543 vector.dot(r1);
544 c1= c1/boxVectors[1].dot(r1);
545 int c2= IntegerMath<IntScalarType>::positive_modulo(vector.dot(r2),boxVectors[2].dot(r2)) -
546 vector.dot(r2);
547 c2= c2/boxVectors[2].dot(r2);
548 vector= vector+c0*boxVectors[0]+c1*boxVectors[1]+c2*boxVectors[2];
549 output.push_back(vector);
550 }
551 }
552 }
553
554 if(!filename.empty()) {
555 std::ofstream file;
556 file.open(filename);
557 if (!file) std::cerr << "Unable to open file";
558 file << output.size() << std::endl;
559 file << "Lattice=\"";
560
561 for(const auto& vector:boxVectors) {
562 file << vector.cartesian().transpose() << " ";
563 }
564 file << "\" Properties=atom_types:I:1:pos:R:3" << std::endl;
565
566 for (const auto &vector: output)
567 file << 1 << " " << vector.cartesian().transpose() << std::endl;
568
569 file.close();
570 }
571 return output;
572 }
573
574
575 template<int dim> template<int dm>
576 typename std::enable_if<dm==2,std::vector<LatticeVector<dim>>>::type
577 Lattice<dim>::box(const std::vector<LatticeVector<dim>>& boxVectors, const std::string& filename) const
578 {
579 for(const LatticeVector<dim>& boxVector : boxVectors)
580 {
581 assert(this == &boxVector.lattice && "Box vectors belong to different lattice.");
582 }
583
584 // Form the box lattice
585 MatrixDimD C;
586 for (int i=0; i<dim; ++i) {
587 C.col(i) = boxVectors[i].cartesian();
588 }
589 assert(C.determinant() != 0);
590 Lattice<dim> boxLattice(C);
591
592 // Form plane parallel vectors v0 and v1
593 //LatticeDirection<dim> v0(this->latticeDirection(boxVectors[0].cartesian()));
594 LatticeDirection<dim> v0(
595 (VectorDimI) boxVectors[0]/IntegerMath<IntScalarType>::gcd(boxVectors[0]),
596 *this);
597
598 ReciprocalLatticeVector<dim> temp(*this);
599 temp << v0.latticeVector()(1),-v0.latticeVector()(0);
600 LatticeDirection<dim> v1(planeParallelLatticeBasis(temp,true)[0]);
601
602 /*
603 IntScalarType x,y;
604 IntegerMath<IntScalarType>::extended_gcd(v0.latticeVector()(0),-v0.latticeVector()(1),x,y);
605 LatticeVector<dim> temp(*this);
606 temp << y,x;
607
608 LatticeDirection<dim> v1(temp);
609 */
610
611 int areaRatio= round(abs(boxLattice.latticeBasis.determinant()/
612 (*this).latticeBasis.determinant()));
613
614 int scale1= abs(IntegerMath<IntScalarType>::gcd(boxVectors[0]));
615 int scale2= round(areaRatio/scale1);
616
617 std::vector<LatticeVector<dim>> output;
618
619 // r0 and r1 are reciprocal vectors perpendicular to boxVectors[1] and boxVectors[0], respectively.
620 ReciprocalLatticeVector<dim> r1_temp(*this), r0_temp(*this);
621 r1_temp << boxVectors[0](1),-boxVectors[0](0);
622 if (r1_temp.dot(boxVectors[1]) < 0)
623 r1_temp=-1*r1_temp;
624 r0_temp << boxVectors[1](1),-boxVectors[1](0);
625 if (r0_temp.dot(boxVectors[0]) < 0)
626 r0_temp=-1*r0_temp;
627 assert(r1_temp.dot(boxVectors[0])==0);
628 assert(r0_temp.dot(boxVectors[1])==0);
629 ReciprocalLatticeDirection<dim> r1(r1_temp), r0(r0_temp);
630
631 for(int j=0; j<scale1; ++j)
632 {
633 for(int k=0; k<scale2; ++k) {
634 LatticeVector<dim> vector(j*v0.latticeVector()+k*v1.latticeVector());
635 int c1= IntegerMath<IntScalarType>::positive_modulo(vector.dot(r1),boxVectors[1].dot(r1)) -
636 vector.dot(r1);
637 c1= c1/boxVectors[1].dot(r1);
638 int c2= IntegerMath<IntScalarType>::positive_modulo(vector.dot(r0),boxVectors[0].dot(r0)) -
639 vector.dot(r0);
640 c2= c2/boxVectors[0].dot(r0);
641 vector= vector+c1*boxVectors[1]+c2*boxVectors[0];
642 output.push_back(vector);
643 }
644 }
645
646 if(!filename.empty()) {
647 std::ofstream file;
648 file.open(filename);
649 if (!file) std::cerr << "Unable to open file";
650 file << output.size() << std::endl;
651 file << "Lattice=\"";
652
653 for(const auto& vector:boxVectors) {
654 file << vector.cartesian().transpose() << " 0 ";
655 }
656 file << " 0 0 1 ";
657 file << "\" Properties=atom_types:I:1:pos:R:3" << std::endl;
658
659 for (const auto &vector: output)
660 file << 1 << " " << vector.cartesian().transpose() << " 0 " << std::endl;
661
662 file.close();
663 }
664 return output;
665 }
666
667
668 template class Lattice<1>;
669
670 template class Lattice<2>;
671 template std::vector<typename Lattice<2>::MatrixDimD> Lattice<2>::generateCoincidentLattices<2>(
672 const double &maxStrain, const int &maxDen, const int &N) const;
673 template std::vector<typename Lattice<2>::MatrixDimD> Lattice<2>::generateCoincidentLattices<2>(
674 const Lattice<2> &undeformedLattice, const double &maxStrain, const int &maxDen, const int &N) const;
675 template std::vector<LatticeVector<2>> Lattice<2>::box<2>(const std::vector<LatticeVector<2>> &boxVectors,
676 const std::string &filename) const;
677
678
679 template class Lattice<3>;
680 template std::vector<typename Lattice<3>::MatrixDimD> Lattice<3>::generateCoincidentLattices<3>(
681 const ReciprocalLatticeDirection<3> &rd, const double &maxDen, const int& N) const;
682 template std::vector<LatticeVector<3>> Lattice<3>::box<3>(const std::vector<LatticeVector<3>> &boxVectors,
683 const std::string &filename) const;
684
685 template class Lattice<4>;
686 template class Lattice<5>;
687 } // namespace oILAB
688#endif
Lattice class.
Definition Lattice.h:31
LatticeDirection< dim > latticeDirection(const VectorDimD &d, const double &tol=FLT_EPSILON) const
Returns the lattice direction along a vector.
Definition Lattice.cpp:54
typename LatticeCore< dim >::VectorDimD VectorDimD
Definition Lattice.h:34
typename LatticeCore< dim >::MatrixDimI MatrixDimI
Definition Lattice.h:37
std::enable_if< dm==3, std::vector< MatrixDimD > >::type generateCoincidentLattices(const ReciprocalLatticeDirection< dim > &rd, const double &maxDen=100, const int &N=100) const
Definition Lattice.cpp:310
ReciprocalLatticeDirection< dim > reciprocalLatticeDirection(const VectorDimD &d, const double &tol=FLT_EPSILON) const
Returns the reciprocal lattice direction along a vector.
Definition Lattice.cpp:81
std::vector< LatticeDirection< dim > > planeParallelLatticeBasis(const ReciprocalLatticeDirection< dim > &l, const bool &useRLLL=false) const
Given a reciprocal lattice direction , this function returns a plane-parallel lattice basis ,...
Definition Lattice.cpp:167
double interPlanarSpacing(const ReciprocalLatticeDirection< dim > &r) const
Computes the interplanar spacing.
Definition Lattice.cpp:301
LatticeVector< dim > latticeVector(const VectorDimD &p) const
Returns a lattice vector (in the current lattice) with Cartesian coordinates p.
Definition Lattice.cpp:153
typename LatticeCore< dim >::MatrixDimD MatrixDimD
Definition Lattice.h:35
RationalLatticeDirection< dim > rationalLatticeDirection(const VectorDimD &d, const typename BestRationalApproximation::LongIntType &maxDen, const double &magnitudeTol, const double &directionTol=FLT_EPSILON) const
Definition Lattice.cpp:108
Lattice(const MatrixDimD &A, const MatrixDimD &Q=MatrixDimD::Identity())
Definition Lattice.cpp:44
ReciprocalLatticeVector< dim > reciprocalLatticeVector(const VectorDimD &p) const
Returns a reciprocal lattice vector (in the dual of the current lattice) with Cartesian coordinates p...
Definition Lattice.cpp:160
const MatrixDimD reciprocalBasis
Definition Lattice.h:42
std::vector< ReciprocalLatticeDirection< dim > > directionOrthogonalReciprocalLatticeBasis(const LatticeDirection< dim > &l, const bool &useRLLL=false) const
Given a lattice direction , this function returns a direction-orthogonal reciprocal lattice basis ,...
Definition Lattice.cpp:233
std::enable_if< dm==3, std::vector< LatticeVector< dim > > >::type box(const std::vector< LatticeVector< dim > > &boxVectors, const std::string &filename="") const
Definition Lattice.cpp:465
typename LatticeCore< dim >::IntScalarType IntScalarType
Definition Lattice.h:38
typename LatticeCore< dim >::VectorDimI VectorDimI
Definition Lattice.h:36
RationalReciprocalLatticeDirection< dim > rationalReciprocalLatticeDirection(const VectorDimD &d, const typename BestRationalApproximation::LongIntType &maxDen, const double &magnitudeTol, const double &directionTol=FLT_EPSILON) const
Definition Lattice.cpp:131
const MatrixDimD latticeBasis
Definition Lattice.h:41
LatticeVector class.
const Lattice< dim > & lattice
IntScalarType dot(const ReciprocalLatticeVector< dim > &other) const
VectorDimD cartesian() const
const MatrixType & reducedBasis() const
Definition RLLL.cpp:185
const Eigen::Matrix< long long int, Eigen::Dynamic, Eigen::Dynamic > & unimodularMatrix() const
Definition RLLL.cpp:191
std::vector< Rational< IntScalarType > > approximations
IntScalarType dot(const LatticeVector< dim > &other) const
std::vector< std::pair< int, int > > farey(int limit, const bool firstQuadrant)
Definition Farey.cpp:5
static Eigen::Matrix< IntScalarType, Eigen::Dynamic, Eigen::Dynamic > ccum(const Eigen::MatrixBase< T > &qin)
static IntScalarType positive_modulo(IntScalarType i, IntScalarType n)
Definition IntegerMath.h:19
static IntScalarType extended_gcd(IntScalarType a, IntScalarType b, IntScalarType &x, IntScalarType &y)
static IntScalarType gcd(const IntScalarType &a, const IntScalarType &b)
Definition IntegerMath.h:26
static Eigen::Vector< IntScalarType, Eigen::Dynamic > solveBezout(const Eigen::MatrixBase< T > &a)
LatticeDirection class.
const LatticeVector< dim > & latticeVector() const
const ReciprocalLatticeVector< dim > & reciprocalLatticeVector() const
Returns a constant reference to the base class (ReciprocalLatticeVector)