7#ifndef gbLAB_Lattice_cpp_
8#define gbLAB_Lattice_cpp_
10#include <Eigen/Eigenvalues>
12#include "../../include/Lattices/LatticeModule.h"
13#include "../../include/Math/GramMatrix.h"
46 ,reciprocalBasis(latticeBasis.inverse().transpose())
56 RLLL rlll(latticeBasis,0.75);
67 const double crossNorm(sqrt(G.determinant()));
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");
83 RLLL rlll(latticeBasis,0.75);
94 const double crossNorm(sqrt(abs(G.determinant())));
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");
110 const double& magnitudeTol,
111 const double& directionTol)
const
117 if ((rld.
cartesian() - d).squaredNorm() > magnitudeTol)
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");
133 const double& magnitudeTol,
134 const double& directionTol)
const
140 if ((rrld.
cartesian() - d).squaredNorm() > magnitudeTol)
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");
168 const bool& useRLLL)
const
170 assert(
this == &l.
lattice &&
"Vectors belong to different Lattices.");
173 std::vector<LatticeDirection<dim>> out;
175 for(
const auto& column : matrix.colwise())
178 if (columnIndex != 0)
184 if(!useRLLL)
return out;
186 Eigen::MatrixXd planeParallelLatticeBasisCartesian(dim,dim-1);
188 for(
auto it=std::next(out.begin()); it!=out.end(); ++it)
190 planeParallelLatticeBasisCartesian.col(index)= (*it).cartesian();
195 const auto& rlll=
RLLL(planeParallelLatticeBasisCartesian, 0.75);
196 planeParallelLatticeBasisCartesian= rlll.reducedBasis();
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)
216 Eigen::MatrixXd pseudoInverse(dim-1,dim);
217 pseudoInverse= (planeParallelLatticeBasisCartesian.transpose() * planeParallelLatticeBasisCartesian).inverse() *
218 (planeParallelLatticeBasisCartesian.transpose());
221 for(
int i=0;i<dim-1;i++)
223 temp= temp + round(out[0].cartesian().dot(pseudoInverse.row(i)))*out[i+1].latticeVector();
234 const bool& useRLLL)
const
236 assert(
this == &l.
lattice &&
"Vectors belong to different Lattices.");
239 std::vector<ReciprocalLatticeDirection<dim>> out;
241 for(
const auto& column : matrix.colwise())
244 if (columnIndex != 0) {
253 if (!useRLLL)
return out;
255 Eigen::MatrixXd directionOrthogonalReciprocalLatticeBasisCartesian(dim,dim-1);
257 for(
auto it=std::next(out.begin()); it!=out.end(); ++it)
259 directionOrthogonalReciprocalLatticeBasisCartesian.col(index)= (*it).cartesian();
263 if constexpr(dim>2) {
264 const auto &rlll =
RLLL(directionOrthogonalReciprocalLatticeBasisCartesian, 0.75);
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);
275 Eigen::MatrixXd pseudoInverse(dim-1,dim);
276 pseudoInverse= (directionOrthogonalReciprocalLatticeBasisCartesian.transpose() * directionOrthogonalReciprocalLatticeBasisCartesian).inverse() *
277 (directionOrthogonalReciprocalLatticeBasisCartesian.transpose());
280 for(
int i=0;i<dim-1;i++)
282 temp= temp + round(out[0].cartesian().dot(pseudoInverse.row(i)))*out[i+1].reciprocalLatticeVector();
304 throw(std::runtime_error(
"The input reciprocal lattice vectors does not belong to the current lattice."));
308 template<
int dim>
template<
int dm>
309 typename std::enable_if<dm==3,std::vector<typename Lattice<dim>::MatrixDimD>>::type
312 std::vector<MatrixDimD> output;
313 std::map<IntScalarType,MatrixDimD> temp;
314 auto basis= planeParallelLatticeBasis(rd);
318 auto b1= basis[1].cartesian();
319 auto b2= basis[2].cartesian();
324 std::vector<std::pair<int,int>> coPrimePairs=
farey(N,
false);
325 for(
const auto& pair : coPrimePairs)
328 auto q2= pair.first*b1+pair.second*b2;
330 if (abs(q2.norm())<epsilon )
continue;
331 double ratio= q2.norm()/q1.norm();
333 double error= ratio-
static_cast<double>(bra.
num)/bra.
den;
334 if (abs(error) > epsilon)
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);
343 rotation= Eigen::AngleAxis<double>(theta,rd.
cartesian().normalized());
345 temp.insert(std::pair<IntScalarType,MatrixDimD>(key,rotation));
348 std::transform(temp.begin(), temp.end(),
349 std::back_inserter(output),
350 [](
const std::pair<IntScalarType,MatrixDimD>& p) {
357 template<
int dim>
template<
int dm>
358 typename std::enable_if<dm==2,std::vector<typename Lattice<dim>::MatrixDimD>>::type
361 std::vector<MatrixDimD> output(generateCoincidentLattices(*
this,maxStrain,maxDen,N));
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,
373 int numberOfConfigurations= 0;
374 const int maxConfigurations= 80;
375 std::vector<MatrixDimD> output;
376 std::map<IntScalarType,MatrixDimD> temp;
379 std::vector<std::pair<int,int>> coPrimePairs=
farey(N,
false);
383 for(
const auto& pair1 : coPrimePairs)
386 pIndices << pair1.first, pair1.second;
388 double ratio= q2.
cartesian().norm() / latticeBasis.col(0).norm();
391 if (abs(maxStrain) < epsilon)
394 double error = ratio -
static_cast<double>(alpha.
num) / alpha.
den;
395 if (abs(error) > epsilon)
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);
405 temp.insert(std::pair<IntScalarType,MatrixDimD>(key,rotation.toRotationMatrix()));
410 RationalApproximations <IntScalarType> alphaSequence(ratio, maxDen, ratio*maxStrain);
415 mn.col(0) = q2 * alpha.d;
416 md.col(0).setConstant(alpha.n);
419 (q2ByAlpha.
cartesian().squaredNorm() - latticeBasis.col(0).squaredNorm()) /
420 latticeBasis.col(0).squaredNorm();
421 if (abs(s1) > maxStrain)
continue;
423 for (
const auto &pair2: coPrimePairs) {
425 qIndices << pair2.first, pair2.second;
427 double ratio2= r2.
cartesian().norm() / latticeBasis.col(1).norm();
431 double s2 = (r2ByBeta.
cartesian().squaredNorm() - latticeBasis.col(1).squaredNorm()) /
432 latticeBasis.col(1).squaredNorm();
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;
440 r2ByBeta.
cartesian() * reciprocalBasis.col(1).transpose();
441 if (F.determinant() < 0)
continue;
443 output.push_back( F );
444 numberOfConfigurations++;
445 std::cout << numberOfConfigurations << std::endl;
446 if (numberOfConfigurations == maxConfigurations)
454 if (abs(maxStrain) < epsilon)
455 std::transform(temp.begin(), temp.end(),
456 std::back_inserter(output),
457 [](
const std::pair<IntScalarType,MatrixDimD>& p) {
463 template<
int dim>
template<
int dm>
464 typename std::enable_if<dm==3,std::vector<LatticeVector<dim>>>::type
469 assert(
this == &boxVector.lattice &&
"Box vectors belong to different lattice.");
474 for (
int i=0; i<dim; ++i) {
475 C.col(i) = boxVectors[i].cartesian();
477 assert(C.determinant() != 0);
480 auto planeParallelBasis= planeParallelLatticeBasis(boxVectors[1].cross(boxVectors[2]),
true);
483 for(
int i=0; i<dim; ++i) {
484 A.col(i) = planeParallelBasis[i].cartesian();
485 mat.col(i)=planeParallelBasis[i].latticeVector();
505 int scale2= round(areaRatio/scale1);
506 int scale0= round(abs(C.determinant()/latticeBasis.determinant()/areaRatio));
508 std::vector<LatticeVector<dim>> output;
513 r0_temp= (boxVectors[1].cross(boxVectors[2])).reciprocalLatticeVector();
514 if (r0_temp.dot(boxVectors[0]) < 0)
516 r1_temp= (boxVectors[2].cross(boxVectors[0])).reciprocalLatticeVector();
517 if (r1_temp.dot(boxVectors[1]) < 0)
519 r2_temp= (boxVectors[0].cross(boxVectors[1])).reciprocalLatticeVector();
520 if (r2_temp.
dot(boxVectors[2]) < 0)
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);
527 for(
int i=0; i<scale0; ++i)
529 for(
int j=0; j<scale1; ++j)
531 for(
int k=0; k<scale2; ++k) {
533 vectorIn_l << i, 0, 0;
541 c0= c0/boxVectors[0].dot(r0);
544 c1= c1/boxVectors[1].dot(r1);
547 c2= c2/boxVectors[2].dot(r2);
548 vector= vector+c0*boxVectors[0]+c1*boxVectors[1]+c2*boxVectors[2];
549 output.push_back(vector);
554 if(!filename.empty()) {
557 if (!file) std::cerr <<
"Unable to open file";
558 file << output.size() << std::endl;
559 file <<
"Lattice=\"";
561 for(
const auto& vector:boxVectors) {
562 file << vector.
cartesian().transpose() <<
" ";
564 file <<
"\" Properties=atom_types:I:1:pos:R:3" << std::endl;
566 for (
const auto &vector: output)
567 file << 1 <<
" " << vector.cartesian().transpose() << std::endl;
575 template<
int dim>
template<
int dm>
576 typename std::enable_if<dm==2,std::vector<LatticeVector<dim>>>::type
581 assert(
this == &boxVector.lattice &&
"Box vectors belong to different lattice.");
586 for (
int i=0; i<dim; ++i) {
587 C.col(i) = boxVectors[i].cartesian();
589 assert(C.determinant() != 0);
590 Lattice<dim> boxLattice(C);
594 LatticeDirection<dim> v0(
598 ReciprocalLatticeVector<dim> temp(*
this);
599 temp << v0.latticeVector()(1),-v0.latticeVector()(0);
600 LatticeDirection<dim> v1(planeParallelLatticeBasis(temp,
true)[0]);
611 int areaRatio= round(abs(boxLattice.latticeBasis.determinant()/
612 (*this).latticeBasis.determinant()));
615 int scale2= round(areaRatio/scale1);
617 std::vector<LatticeVector<dim>> output;
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)
624 r0_temp << boxVectors[1](1),-boxVectors[1](0);
625 if (r0_temp.dot(boxVectors[0]) < 0)
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);
631 for(
int j=0; j<scale1; ++j)
633 for(
int k=0; k<scale2; ++k) {
634 LatticeVector<dim> vector(j*v0.latticeVector()+k*v1.latticeVector());
637 c1= c1/boxVectors[1].dot(r1);
640 c2= c2/boxVectors[0].dot(r0);
641 vector= vector+c1*boxVectors[1]+c2*boxVectors[0];
642 output.push_back(vector);
649 if (!file) std::cerr <<
"Unable to open file";
650 file << output.size() << std::endl;
651 file <<
"Lattice=\"";
653 for(
const auto& vector:boxVectors) {
654 file << vector.cartesian().transpose() <<
" 0 ";
657 file <<
"\" Properties=atom_types:I:1:pos:R:3" << std::endl;
659 for (
const auto &vector: output)
660 file << 1 <<
" " << vector.cartesian().transpose() <<
" 0 " << std::endl;
668 template class Lattice<1>;
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;
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;
685 template class Lattice<4>;
686 template class Lattice<5>;
long long int LongIntType
LatticeDirection< dim > latticeDirection(const VectorDimD &d, const double &tol=FLT_EPSILON) const
Returns the lattice direction along a vector.
typename LatticeCore< dim >::VectorDimD VectorDimD
typename LatticeCore< dim >::MatrixDimI MatrixDimI
std::enable_if< dm==3, std::vector< MatrixDimD > >::type generateCoincidentLattices(const ReciprocalLatticeDirection< dim > &rd, const double &maxDen=100, const int &N=100) const
ReciprocalLatticeDirection< dim > reciprocalLatticeDirection(const VectorDimD &d, const double &tol=FLT_EPSILON) const
Returns the reciprocal lattice direction along a vector.
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 ,...
double interPlanarSpacing(const ReciprocalLatticeDirection< dim > &r) const
Computes the interplanar spacing.
LatticeVector< dim > latticeVector(const VectorDimD &p) const
Returns a lattice vector (in the current lattice) with Cartesian coordinates p.
typename LatticeCore< dim >::MatrixDimD MatrixDimD
RationalLatticeDirection< dim > rationalLatticeDirection(const VectorDimD &d, const typename BestRationalApproximation::LongIntType &maxDen, const double &magnitudeTol, const double &directionTol=FLT_EPSILON) const
Lattice(const MatrixDimD &A, const MatrixDimD &Q=MatrixDimD::Identity())
ReciprocalLatticeVector< dim > reciprocalLatticeVector(const VectorDimD &p) const
Returns a reciprocal lattice vector (in the dual of the current lattice) with Cartesian coordinates p...
const MatrixDimD reciprocalBasis
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 ,...
std::enable_if< dm==3, std::vector< LatticeVector< dim > > >::type box(const std::vector< LatticeVector< dim > > &boxVectors, const std::string &filename="") const
typename LatticeCore< dim >::IntScalarType IntScalarType
typename LatticeCore< dim >::VectorDimI VectorDimI
RationalReciprocalLatticeDirection< dim > rationalReciprocalLatticeDirection(const VectorDimD &d, const typename BestRationalApproximation::LongIntType &maxDen, const double &magnitudeTol, const double &directionTol=FLT_EPSILON) const
const MatrixDimD latticeBasis
const Lattice< dim > & lattice
IntScalarType dot(const ReciprocalLatticeVector< dim > &other) const
VectorDimD cartesian() const
const MatrixType & reducedBasis() const
const Eigen::Matrix< long long int, Eigen::Dynamic, Eigen::Dynamic > & unimodularMatrix() const
std::vector< Rational< IntScalarType > > approximations
const Lattice< dim > & lattice
IntScalarType dot(const LatticeVector< dim > &other) const
VectorDimD cartesian() const
std::vector< std::pair< int, int > > farey(int limit, const bool firstQuadrant)
static Eigen::Matrix< IntScalarType, Eigen::Dynamic, Eigen::Dynamic > ccum(const Eigen::MatrixBase< T > &qin)
static IntScalarType positive_modulo(IntScalarType i, IntScalarType n)
static IntScalarType extended_gcd(IntScalarType a, IntScalarType b, IntScalarType &x, IntScalarType &y)
static IntScalarType gcd(const IntScalarType &a, const IntScalarType &b)
static Eigen::Vector< IntScalarType, Eigen::Dynamic > solveBezout(const Eigen::MatrixBase< T > &a)
const LatticeVector< dim > & latticeVector() const
VectorDimD cartesian() const
VectorDimD cartesian() const
const ReciprocalLatticeVector< dim > & reciprocalLatticeVector() const
Returns a constant reference to the base class (ReciprocalLatticeVector)