oILAB
Loading...
Searching...
No Matches
BiCrystal.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_BiCrystal_cpp_
8#define gbLAB_BiCrystal_cpp_
9
10#include <LatticeModule.h>
11#include <numbers>
12
13namespace gbLAB
14{
15
16 template <int dim>
19 {
21 for(int i=0;i<dim;++i)
22 {
23 const auto& dii(sd.matrixD()(i,i));
24 M(i,i)=dii/IntegerMath<IntScalarType>::gcd(rm.mu,dii);
25 }
26 return M;
27 }
28
29 template <int dim>
32 {
34 for(int i=0;i<dim;++i)
35 {
36 const auto& dii=sd.matrixD()(i,i);
37 N(i,i)=rm.mu/IntegerMath<IntScalarType>::gcd(rm.mu,dii);
38 }
39 return N;
40 }
41
42 template <int dim>
44 const typename BiCrystal<dim>::MatrixDimI& N)
45 {
47 for(int col=0; col<dim; ++col)
48 {
50 for (int i = 0; i < dim; ++i)
51 {
52 IntScalarType a = N(i, i);
53 IntScalarType b = -M(i, i);
54 IntScalarType c = -(i==col);
56 }
57 LambdaA.col(col)= M*y;
58 }
59 return LambdaA;
60 }
61
62 template <int dim>
64 const typename BiCrystal<dim>::MatrixDimI& N)
65 {
66 typename BiCrystal<dim>::MatrixDimI LambdaB;
67 for(int col=0; col<dim; ++col)
68 {
70 for (int i = 0; i < dim; ++i)
71 {
72 IntScalarType a = N(i, i);
73 IntScalarType b = -M(i, i);
74 IntScalarType c = (i==col);
76 }
77 LambdaB.col(col)= N*x;
78 }
79 return LambdaB;
80 }
81
82 template <int dim>
84 const Lattice<dim>& B,
86 const typename BiCrystal<dim>::MatrixDimI& M,
87 const typename BiCrystal<dim>::MatrixDimI& N,
88 const bool& useRLLL)
89 {
90 // The transition matrix is T=P/sigma, where P=rm.integerMatrix is
91 // an integer matrix and sigma=rm.sigma is an integer
92 // The integer matrix P can be decomposed as P=X*D*Y using the Smith decomposition.
93 // X and Y are unimodular, and D is diagonal with D(k,k) dividing D(k+1,k+1)
94 // The decomposition also computes the matices U and V such that D=U*P*V
95 // From T=inv(A)*B=P/sigma=X*D*Y/sigma=X*D*inv(V)/sigma, we have
96 // B1*(sigma*I)=A1*D
97 // where
98 // B1=B*V
99 // A1=A*X
100 // Since V and X are unimodular matrices, B1 and A1 are new bases
101 // of the lattices B and A, respectively. Moreover, since
102 // (sigma*I) and D are diagonal, the columns of B1 and A1 are
103 // proportional, with rational proportionality factors different for each column.
104 // For example, columns "i" read
105 // b1_i*sigma=a1_i*D(i,i)
106 // Therefore, the i-th primitive vectors of the CSL is
107 // c_i=b1_i*sigma/gcd(sigma,D(i,i))=a1_i*D(i,i)/gcd(sigma,D(i,i))
108 // or, in matrix form
109 // C=B1*N=A1*M, that is
110 // C=B*V*N=A*X*M
111 // where M=diag(D(i,i)/gcd(sigma,D(i,i))) and
112 // N=diag(sigma/gcd(sigma,D(i,i))) and
113
114 const auto C1(A.latticeBasis*(sd.matrixX()*M).template cast<double>());
115 const auto C2(B.latticeBasis*(sd.matrixV()*N).template cast<double>());
116 if ((C1-C2).norm()/C1.norm()>FLT_EPSILON || (C1-C2).norm()/C2.norm()>FLT_EPSILON)
117 {
118 throw std::runtime_error("CSL calculation failed.\n");
119 }
120
121 if(useRLLL)
122 {
123 return RLLL(0.5*(C1+C2),0.75).reducedBasis();
124 }
125 else
126 {
127 return 0.5*(C1+C2);
128 }
129 }
130
131 template <int dim>
133 const Lattice<dim>& B,
134 const SmithDecomposition<dim>& sd,
135 const typename BiCrystal<dim>::MatrixDimI& M,
136 const typename BiCrystal<dim>::MatrixDimI& N,
137 const bool& useRLLL)
138 {
139
140 const auto D1(A.latticeBasis*sd.matrixX().template cast<double>()*N.template cast<double>().inverse());
141 const auto D2(B.latticeBasis*sd.matrixV().template cast<double>()*M.template cast<double>().inverse());
142 if ((D1-D2).norm()/D1.norm()>FLT_EPSILON || (D1-D2).norm()/D2.norm()>FLT_EPSILON)
143 {
144 throw std::runtime_error("DSCL calculation failed.\n");
145 }
146 if(useRLLL)
147 {
148 return RLLL(0.5*(D1+D2),0.75).reducedBasis();
149 }
150 else
151 {
152 return 0.5*(D1+D2);
153 }
154 }
155
156
157 template <int dim>
159 const Lattice<dim>& B_in,
160 const bool& useRLLL) try :
161 /* init */ RationalMatrix<dim>(A_in.reciprocalBasis.transpose()*B_in.latticeBasis)
162 /* init */,SmithDecomposition<dim>(this->integerMatrix)
163 /* init */,A(A_in)
164 /* init */,B(B_in)
165 /* init */,M(getM(*this,*this))
166 /* init */,N(getN(*this,*this))
167 /* init */,sigmaA(round(M.template cast<double>().determinant()))
168 /* init */,sigmaB(round(N.template cast<double>().determinant()))
169 /* init */,sigma(std::abs(sigmaA)==std::abs(sigmaB)? std::abs(sigmaA) : 0)
170 /* init */, csl(getCSLBasis (A,B,*this,M,N,useRLLL),MatrixDimD::Identity())
171 /* init */,dscl(getDSCLBasis(A,B,*this,M,N,useRLLL),MatrixDimD::Identity())
172 /* init */,Ap(A.latticeBasis*this->matrixX().template cast<double>())
173 /* init */,Bp(B.latticeBasis*this->matrixV().template cast<double>())
174 /* init */,LambdaA(getLambdaA(M,N))
175 /* init */,LambdaB(getLambdaB(M,N))
176 {
177
178 if(true)
179 {//verify that CSL can be obtained as multiple of A and B
180
181 const MatrixDimD tempA(A.reciprocalBasis.transpose()*csl.latticeBasis);
182 if ((tempA-tempA.array().round().matrix()).norm()/tempA.norm()>FLT_EPSILON)
183 {
184 //std::cout << (tempA-tempA.array().round().matrix()).norm()/tempA.norm() << std::endl;
185 throw std::runtime_error("CSL is not a multiple of lattice A.\n");
186 }
187
188 const MatrixDimD tempB(B.reciprocalBasis.transpose()*csl.latticeBasis);
189 if ((tempB-tempB.array().round().matrix()).norm()/tempB.norm()>FLT_EPSILON)
190 {
191 //std::cout << (tempB-tempB.array().round().matrix()).norm()/tempB.norm() << std::endl;
192 throw std::runtime_error("CSL is not a multiple of lattice B\n");
193 }
194 }
195
196 if(true)
197 {//verify that A and B are multiples of DSCL
198
199 const MatrixDimD tempA(dscl.reciprocalBasis.transpose()*A.latticeBasis);
200 if ((tempA-tempA.array().round().matrix()).norm()/tempA.norm()>FLT_EPSILON)
201 {
202 //std::cout << (tempA-tempA.array().round().matrix()).norm()/tempA.norm() << std::endl;
203 throw std::runtime_error("Lattice A is not a multiple of the DSCL\n");
204 }
205
206 const MatrixDimD tempB(dscl.reciprocalBasis.transpose()*B.latticeBasis);
207 if ((tempB-tempB.array().round().matrix()).norm()/tempB.norm()>FLT_EPSILON)
208 {
209 //std::cout << (tempB-tempB.array().round().matrix()).norm()/tempB.norm() << std::endl;
210 throw std::runtime_error("Lattice B is not a multiple of the DSCL\n");
211 }
212 }
213
214 if(true)
215 {//verify LambdaA + LambdaB = I
216 if (!(LambdaA+LambdaB).isIdentity())
217 {
218 throw std::runtime_error("LambdaA + LambdaB != I\n");
219 }
220 }
221
222 }
223
224 catch(std::runtime_error& e)
225 {
226 std::cout << e.what() << std::endl;
227 throw(std::runtime_error("Bicrystal construction failed. "));
228 }
229
230 template<int dim>
232 {
233 VectorDimI integerCoordinates;
234
235 if(&(v.lattice) == &(this->A))
236 return v;
237 else if(&(v.lattice) == &(this->csl))
238 // U*M*v
239 integerCoordinates= this->matrixX() * M * v;
240 else
241 throw(std::runtime_error("The input lattice vector should belong "
242 "to lattice A or the CSL"));
243 auto temp= LatticeVector<dim>(integerCoordinates,A);
244 if (temp.cartesian().dot(v.cartesian()) < 0) temp= -1*temp;
245 return temp;
246 }
247
248 template<int dim>
250 {
251 VectorDimI integerCoordinates;
252
253 if(&(v.lattice) == &(this->B))
254 return v;
255 else if(&(v.lattice) == &(this->csl))
256 // V*N*v
257 integerCoordinates= this->matrixV() * N * v;
258 else
259 throw(std::runtime_error("The input lattice vector should belong "
260 "to lattice B or the CSL"));
261 auto temp= LatticeVector<dim>(integerCoordinates,B);
262 if (temp.cartesian().dot(v.cartesian()) < 0) temp= -1*temp;
263 return temp;
264 }
265
266 template<int dim>
268 {
269 VectorDimI integerCoordinates;
270
273
274 if(&(v.lattice) == &(this->A))
275 // N*inv(U)*v
276 integerCoordinates = N * adjX * v;
277 else if(&(v.lattice) == &(this->B))
278 // M*inv(V)*v
279 integerCoordinates = M * adjV * v;
280 else if(&(v.lattice) == &(this->csl))
281 // N*M*v
282 integerCoordinates = N * M * v;
283 else if(&(v.lattice) == &(this->dscl))
284 return LatticeVector<dim>(v);
285 else
286 throw(std::runtime_error("The input lattice vector should belong to one of the four lattices of the bicrystal"));
287
288 auto temp= LatticeVector<dim>(integerCoordinates,dscl);
289 if (temp.cartesian().dot(v.cartesian()) < 0) temp= -1*temp;
290
291 //return LatticeVector<dim>(integerCoordinates,dscl);
292 return temp;
293 }
294
295
296 template<int dim>
298 {
299 VectorDimI integerCoordinates;
300
306
307 if(&(v.lattice) == &(this->A))
308 // inv(M)*inv(U)*v
309 integerCoordinates = adjM * adjX * v;
310 else if(&(v.lattice) == &(this->B))
311 // inv(N)*inv(V)*v
312 integerCoordinates = adjN * adjV * v;
313 else if(&(v.lattice) == &(this->csl))
314 return LatticeDirection<dim>(v);
315 else if(&(v.lattice) == &(this->dscl))
316 integerCoordinates = adjMN* v;
317 else
318 throw(std::runtime_error("The input reciprocal lattice vector should belong to one of the four reciprocal lattices of the bicrystal"));
319
320 auto temp= LatticeVector<dim>(integerCoordinates,csl);
321 if (temp.cartesian().dot(v.cartesian()) < 0) temp= -1*temp;
322 return LatticeDirection<dim>(temp);
323 }
324 template<int dim>
326 {
327 return LatticeDirection<dim>(getLatticeVectorInD(v));
328 }
329
330 template<int dim>
332 {
333 VectorDimI integerCoordinates;
334
337
338 if(&(rv.lattice) == &(this->A))
340 else if(&(rv.lattice) == &(this->B))
341 // U^-T*inverse(M)*N*V^T
342 integerCoordinates= adjX.transpose() * adjM * N * (this->matrixV()).transpose() * rv;
343 else if(&(rv.lattice) == &(this->csl))
344 // U^-T*inverse(M)
345 integerCoordinates= adjX.transpose() * adjM * rv;
346 else if(&(rv.lattice) == &(this->dscl))
347 // U^-T*N*rv
348 integerCoordinates = adjX.transpose() * N * rv;
349 else
350 throw(std::runtime_error("The input reciprocal lattice vector should belong to one of the four reciprocal lattices of the bicrystal"));
351
352 auto temp= ReciprocalLatticeVector<dim>(integerCoordinates,A);
353 if (temp.cartesian().dot(rv.cartesian()) < 0) temp= -1*temp;
355 }
356 template<int dim>
358 {
359 VectorDimI integerCoordinates;
360
365
366 if(&(rv.lattice) == &(this->A))
367 // V^-T*inverse(N)*M*U^T
368 integerCoordinates= adjV.transpose() * adjN * M * (this->matrixX()).transpose() * rv;
369 else if(&(rv.lattice) == &(this->B))
371 else if(&(rv.lattice) == &(this->csl))
372 // V^-T*inverse(N)*rv
373 integerCoordinates= adjV.transpose() * adjN * rv;
374 else if(&(rv.lattice) == &(this->dscl))
375 // V^-T*M*rv
376 integerCoordinates= adjV.transpose() * M * rv;
377 else
378 throw(std::runtime_error("The input reciprocal lattice vector should belong to one of the four reciprocal lattices of the bicrystal"));
379
380 auto temp= ReciprocalLatticeVector<dim>(integerCoordinates,B);
381 if (temp.cartesian().dot(rv.cartesian()) < 0) temp= -1*temp;
383 }
384 template<int dim>
386 {
387 VectorDimI integerCoordinates;
388
389 if(&(rv.lattice) == &(this->A))
390 // M*U^T*rv
391 integerCoordinates= M*this->matrixX().transpose()*rv;
392 else if(&(rv.lattice) == &(this->B))
393 // N*V^T*rv
394 integerCoordinates= N*this->matrixV().transpose()*rv;
395 else if(&(rv.lattice) == &(this->csl))
397 else if(&(rv.lattice) == &(this->dscl))
398 // M*N*rv
399 integerCoordinates= M*N*rv;
400 else
401 throw(std::runtime_error("The input reciprocal lattice vector should belong to one of the four reciprocal lattices of the bicrystal"));
402
403 auto temp= ReciprocalLatticeVector<dim>(integerCoordinates,csl);
404 if (temp.cartesian().dot(rv.cartesian()) < 0) temp= -1*temp;
406 }
407 template<int dim>
409 {
410 VectorDimI integerCoordinates;
411
414
415 if(&(rv.lattice) == &(this->A))
416 integerCoordinates= adjN * (this->matrixX().transpose()) * rv;
417 else if(&(rv.lattice) == &(this->B))
418 integerCoordinates= adjM * (this->matrixV().transpose()) * rv;
419 else if(&(rv.lattice) == &(this->csl))
420 integerCoordinates= adjN * adjM * rv;
421 else
422 throw(std::runtime_error("The input reciprocal lattice vector should belong to one of the four reciprocal lattices of the bicrystal"));
423
424 auto temp= ReciprocalLatticeVector<dim>(integerCoordinates,dscl);
425 if (temp.cartesian().dot(rv.cartesian()) < 0) temp= -1*temp;
427 }
428
429 template<int dim>
431 {
432 if(&d.lattice != &this->dscl)
433 throw(std::runtime_error("Input vector is not a DSCL vectors"));
434 return LatticeVector<dim>((LambdaA*d).eval(),d.lattice);
435 }
436
437 template<int dim>
439 {
440 if(&d.lattice != &this->dscl)
441 throw(std::runtime_error("Input vector is not a DSCL vectors"));
442 return LatticeVector<dim>((LambdaB*d).eval(),d.lattice);
443 }
444
445 template<int dim> template<int dm>
446 typename std::enable_if<dm==2 || dm==3,std::map<typename BiCrystal<dim>::IntScalarType,Gb<dm>>>::type
448 {
449 if (&d.lattice != &A && &d.lattice != &B)
450 throw std::runtime_error("The tilt axis does not belong to lattices A and B \n");
451 std::vector<Gb<dm>> gbVec;
452 double epsilon=1e-8;
453 int count= -1;
454 IntScalarType keyScale= 1e6;
455 auto basis= d.lattice.directionOrthogonalReciprocalLatticeBasis(d,true);
456 if (dm==3)
457 {
458 for (int i = -div; i <= div; ++i)
459 {
460 for (int j = -div; j <= div; ++j)
461 {
462 if (i==0 && j==0) continue;
463 count++;
464 ReciprocalLatticeVector<dm> rv = i * basis[1].reciprocalLatticeVector() + j * basis[2].reciprocalLatticeVector();
465 try
466 {
467 Gb<dm> gb(*this, rv);
468 gbVec.push_back(gb);
469 }
470 catch(std::runtime_error& e)
471 {
472 std::cout << e.what() << std::endl;
473 std::cout << "Unable to form GB with normal = " << rv << std::endl;
474 std::cout << "moving on to next inclination" << std::endl;
475 }
476 }
477 }
478 }
479 else if(dm==2)
480 {
481 auto rv= basis[0].reciprocalLatticeVector();
482 gbVec.push_back(Gb<dm>(*this, rv));
483 }
484 std::map<IntScalarType,Gb<dm>> gbSet;
485 for(const Gb<dim>& gb:gbVec)
486 {
487 double cosAngle;
488 cosAngle= gb.nA.cartesian().normalized().dot(gbVec[0].nA.cartesian().normalized());
489 if (cosAngle-1>-epsilon) cosAngle= 1.0;
490 if (cosAngle+1<epsilon) cosAngle= -1.0;
491
492 double angle= acos(cosAngle);
493 IntScalarType key= angle*keyScale;
494 gbSet.insert(std::pair<IntScalarType,Gb<dm>>(key,gb));
495 }
496 return gbSet;
497 }
498
499 template<int dim> template<int dm>
500 typename std::enable_if<dm==2 || dm==3,std::vector<LatticeVector<dim>>>::type
501 BiCrystal<dim>::box(std::vector<LatticeVector<dim>>& boxVectors,
502 const double& orthogonality,
503 const int& dsclFactor,
504 std::string filename,
505 bool orient) const
506 {
507 assert(orthogonality>=0.0 && orthogonality<=1.0 &&
508 "The \"orthogonality\" parameter should be between 0.0 and 1.0");
509 assert(dsclFactor>=1 &&
510 "The \"dsclFactor\" should be greater than 1.");
511 assert(boxVectors.size()==dim);
512 for(const auto& boxVector : boxVectors)
513 {
514 assert(&csl == &boxVector.lattice &&
515 "Box vectors do not belong to the CSL.");
516 }
517
518 // Form the box lattice
519 MatrixDimD C;
520 for (int i=0; i<dim; ++i) {
521 C.col(i) = boxVectors[i].cartesian();
522 }
523 assert(abs(C.determinant()) > FLT_EPSILON && "Box volume is equal to zero.");
524
525
526 // Adjust boxVector[0] such that it is as orthogonal as possible to boxVector[1]
527 auto boxVectorTemp= boxVectors[0];
529 if (dim==2)
530 nC= boxVectors[1].cross();
531 if (dim==3)
532 nC= boxVectors[1].cross(boxVectors[2]);
533 auto basis= csl.planeParallelLatticeBasis(nC,true);
534
535 int planesToExplore= nC.stacking();
536 MatrixDimI boxLatticeIndices;
537 boxLatticeIndices.col(0)= boxVectors[0];
538 for (int i=1; i<dim; ++i)
539 boxLatticeIndices.col(i)= boxVectors[i]/IntegerMath<long long int>::gcd(boxVectors[i]);
540 double minDotProduct= std::numbers::pi/2;
541 int minStep;
542
543 auto boxVectorUpdated(boxVectors[0]);
544
545 for(int i=0;i<planesToExplore;++i)
546 {
547 int sign= boxVectors[0].cartesian().dot(basis[0].cartesian())>0? 1 : -1;
548 boxVectorTemp= boxVectors[0]+i*sign*basis[0].latticeVector();
549 boxLatticeIndices.col(0)= boxVectorTemp;
550 Lattice<dim> boxLattice(csl.latticeBasis*boxLatticeIndices.template cast<double>());
551 ReciprocalLatticeVector<dim> rC(boxLattice);
552 rC(0)=1;
553 VectorDimI temp=
554 boxLatticeIndices*boxLattice.planeParallelLatticeBasis(ReciprocalLatticeDirection<dim>(rC),true)[0].latticeVector();
555 boxVectorTemp= LatticeVector<dim>(temp,csl);
556 double dotProduct= abs(acos(boxVectorTemp.cartesian().normalized().dot(rC.cartesian().normalized())-FLT_EPSILON));
557 if(dotProduct<minDotProduct) {
558 minDotProduct = dotProduct;
559 boxVectorUpdated= boxVectorTemp;
560 if (dotProduct < (1-orthogonality)*std::numbers::pi/2)
561 break;
562 }
563
564 }
565 boxVectors[0]=boxVectorUpdated;
566 C.col(0)= boxVectors[0].cartesian();
567
568
569 // form the rotation matrix used to orient the system
570 MatrixDimD rotation= Eigen::Matrix<double,dim,dim>::Identity();;
571 Eigen::Matrix<double,dim,dim-1> orthogonalVectors;
572 if (orient) {
573 if (dim==3) {
574 orthogonalVectors.col(0) = C.col(1).normalized();
575 orthogonalVectors.col(1) = C.col(2).normalized();
576 }
577 else if (dim==2)
578 orthogonalVectors.col(0)= C.col(1).normalized();
579
580 rotation=Rotation<dim>(orthogonalVectors);
581 }
582 //assert((rotation*rotation.transpose()).template isApprox(Eigen::Matrix<double,dim,dim>::Identity())
583 //assert((rotation*rotation.transpose()).isApprox(Eigen::Matrix<double,dim,dim>::Identity())
584 assert((rotation*rotation.transpose()).isApprox(Eigen::Matrix<double,dim,dim>::Identity(),1e-6)
585 && "Cannot orient the grain boundary. Box vectors are not orthogonal.");
586
587
588 std::vector<LatticeVector<dim>> configurationA, configurationB, configurationC, configurationD;
589 std::vector<LatticeVector<dim>> configuration;
590
591 std::vector<LatticeVector<dim>> boxVectorsInA, boxVectorsInB, boxVectorsInD;
592 // calculate boxVectors in A, B, and D
593 for(const auto& boxVector : boxVectors) {
594 boxVectorsInA.push_back(getLatticeVectorInA(boxVector));
595 boxVectorsInB.push_back(getLatticeVectorInB(boxVector));
596 boxVectorsInD.push_back(getLatticeVectorInD(boxVector));
597 }
598
599 // prepare boxVectors for D
600 auto dsclVector=getLatticeDirectionInD(boxVectors[0]).latticeVector();
601 auto nD= getReciprocalLatticeDirectionInD(nC.reciprocalLatticeVector());
602 if(abs((dsclFactor*dsclVector).dot(nD)) < abs(boxVectorsInD[0].dot(nD)))
603 boxVectorsInD[0]= dsclFactor*dsclVector;
604
605 std::vector<LatticeVector<dim>> boxVectorsForA(boxVectorsInA),
606 boxVectorsForB(boxVectorsInB),
607 boxVectorsForC(boxVectors),
608 boxVectorsForD(boxVectorsInD);
609 boxVectorsForA[0]=2*boxVectorsInA[0];
610 boxVectorsForB[0]=2*boxVectorsInB[0];
611 boxVectorsForC[0]=2*boxVectors[0];
612 boxVectorsForD[0]=2*boxVectorsInD[0];
613
614 configurationA= A.box(boxVectorsForA);
615 configurationB= B.box(boxVectorsForB);
616 configurationC= csl.box(boxVectorsForC);
617 configurationD= dscl.box(boxVectorsForD);
618
619 LatticeVector<dim> origin(-1*boxVectors[0]);
620 for(auto& vector : configurationA)
621 vector= vector + LatticeVector<dim>(-1*boxVectorsInA[0]);
622 for(auto& vector : configurationB)
623 vector= vector + LatticeVector<dim>(-1*boxVectorsInB[0]);
624 for(auto& vector : configurationC)
625 vector= vector + origin;
626 for(auto& vector : configurationD)
627 vector= vector + LatticeVector<dim>(-1*boxVectorsInD[0]);
628
629 configuration= configurationA;
630 configuration.insert(configuration.end(),configurationB.begin(),configurationB.end());
631 configuration.insert(configuration.end(),configurationC.begin(),configurationC.end());
632 configuration.insert(configuration.end(),configurationD.begin(),configurationD.end());
633
634 if(!filename.empty()) {
635 std::ofstream file;
636 file.open(filename);
637 if (!file) std::cerr << "Unable to open file";
638 file << configuration.size() << std::endl;
639 file << "Lattice=\"";
640
641 if (dim==2) {
642 file << (rotation*boxVectorsForC[0].cartesian()).transpose() << " 0 ";
643 file << (rotation*boxVectorsForC[1].cartesian()).transpose() << " 0 ";
644 file << " 0 0 1 ";
645 file << "\" Properties=atom_types:I:1:pos:R:3:radius:R:1 PBC=\"F T T\" origin=\"";
646 file << (rotation*origin.cartesian()).transpose() << " 0.0\"" << std::endl;
647 for (const auto &vector: configurationA)
648 file << 1 << " " << (rotation*vector.cartesian()).transpose() << " " << 0.0 << " " << 0.05 << std::endl;
649 for (const auto &vector: configurationB)
650 file << 2 << " " << (rotation*vector.cartesian()).transpose() << " " << 0.0 << " " << 0.05 << std::endl;
651 for (const auto &vector: configurationC)
652 file << 3 << " " << (rotation*vector.cartesian()).transpose() << " " << 0.0 << " " << 0.2 << std::endl;
653 for (const auto &vector: configurationD)
654 file << 4 << " " << (rotation*vector.cartesian()).transpose() << " " << 0.0 << " " << 0.01 << std::endl;
655 }
656 else if (dim==3){
657 file << (rotation*boxVectorsForC[0].cartesian()).transpose() << " ";
658 file << (rotation*boxVectorsForC[1].cartesian()).transpose() << " ";
659 file << (rotation*boxVectorsForC[2].cartesian()).transpose() << " ";
660 file << "\" Properties=atom_types:I:1:pos:R:3:radius:R:1 PBC=\"F T T\" origin=\"";
661 file << (rotation*origin.cartesian()).transpose() << "\"" << std::endl;
662
663 for (const auto &vector: configurationA)
664 file << 1 << " " << (rotation*vector.cartesian()).transpose() << " " << 0.05 << std::endl;
665 for (const auto &vector: configurationB)
666 file << 2 << " " << (rotation*vector.cartesian()).transpose() << " " << 0.05 << std::endl;
667 for (const auto &vector: configurationC)
668 file << 3 << " " << (rotation*vector.cartesian()).transpose() << " " << 0.2 << std::endl;
669 for (const auto &vector: configurationD)
670 file << 4 << " " << (rotation*vector.cartesian()).transpose() << " " << 0.01 << std::endl;
671 }
672
673
674 file.close();
675 }
676 return configuration;
677 }
678
679// template class BiCrystal<1>;
680 template class BiCrystal<2>;
681 template std::map<BiCrystal<2>::IntScalarType, Gb<2>>
683 template std::vector<LatticeVector<2>>
684 BiCrystal<2>::box<2>(std::vector<LatticeVector<2>> &boxVectors,
685 const double &orthogonality, const int &dsclFactor,
686 std::string filename, bool orient) const;
687
688 template class BiCrystal<3>;
689 template std::map<BiCrystal<3>::IntScalarType, Gb<3>>
691 template std::vector<LatticeVector<3>>
692 BiCrystal<3>::box<3>(std::vector<LatticeVector<3>> &boxVectors,
693 const double &orthogonality, const int &dsclFactor,
694 std::string filename, bool orient) const;
695
696 template class BiCrystal<4>;
697 template class BiCrystal<5>;
698
699} // end namespace
700#endif
701
const Lattice< dim > & A
Definition BiCrystal.h:56
LatticeVector< dim > shiftTensorB(const LatticeVector< dim > &d) const
LatticeVector< dim > getLatticeVectorInB(const LatticeVector< dim > &v) const
static MatrixDimI getLambdaA(const MatrixDimI &M, const MatrixDimI &N)
Definition BiCrystal.cpp:43
LatticeVector< dim > getLatticeVectorInD(const LatticeVector< dim > &v) const
ReciprocalLatticeDirection< dim > getReciprocalLatticeDirectionInD(const ReciprocalLatticeVector< dim > &v) const
LatticeVector< dim > shiftTensorA(const LatticeVector< dim > &d) const
const Lattice< dim > dscl
DCSL lattice .
Definition BiCrystal.h:90
const MatrixDimI LambdaB
Shift tensor describes the shift in the CSL when lattice is shifted by a DSCL vector.
Definition BiCrystal.h:108
const Lattice< dim > & B
Definition BiCrystal.h:57
ReciprocalLatticeDirection< dim > getReciprocalLatticeDirectionInA(const ReciprocalLatticeVector< dim > &v) const
LatticeCore< dim >::IntScalarType IntScalarType
Definition BiCrystal.h:30
ReciprocalLatticeDirection< dim > getReciprocalLatticeDirectionInC(const ReciprocalLatticeVector< dim > &v) const
LatticeCore< dim >::VectorDimI VectorDimI
Definition BiCrystal.h:33
LatticeVector< dim > getLatticeVectorInA(const LatticeVector< dim > &v) const
ReciprocalLatticeDirection< dim > getReciprocalLatticeDirectionInB(const ReciprocalLatticeVector< dim > &v) const
BiCrystal(const Lattice< dim > &A, const Lattice< dim > &B, const bool &useRLLL=false)
Constructs a bicrystal from two lattices and by computing the parallel bases for lattices ,...
static MatrixDimD getDSCLBasis(const Lattice< dim > &A, const Lattice< dim > &B, const SmithDecomposition< dim > &sd, const MatrixDimI &M, const MatrixDimI &N, const bool &useRLLL)
LatticeDirection< dim > getLatticeDirectionInC(const LatticeVector< dim > &v) const
const Lattice< dim > csl
CSL lattice .
Definition BiCrystal.h:86
static MatrixDimI getN(const RationalMatrix< dim > &rm, const SmithDecomposition< dim > &sd)
Definition BiCrystal.cpp:30
LatticeDirection< dim > getLatticeDirectionInD(const LatticeVector< dim > &v) const
LatticeCore< dim >::MatrixDimI MatrixDimI
Definition BiCrystal.h:34
std::enable_if< dm==2||dm==3, std::map< IntScalarType, Gb< dm > > >::type generateGrainBoundaries(const LatticeDirection< dim > &d, int div=30) const
Given a tilt axis , that belongs to lattices or , this function generate a set of tilt GBs....
static MatrixDimD getCSLBasis(const Lattice< dim > &A, const Lattice< dim > &B, const SmithDecomposition< dim > &sd, const MatrixDimI &M, const MatrixDimI &N, const bool &useRLLL)
Definition BiCrystal.cpp:83
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
LatticeCore< dim >::MatrixDimD MatrixDimD
Definition BiCrystal.h:32
static MatrixDimI getLambdaB(const MatrixDimI &M, const MatrixDimI &N)
Definition BiCrystal.cpp:63
static MatrixDimI getM(const RationalMatrix< dim > &rm, const SmithDecomposition< dim > &sd)
Definition BiCrystal.cpp:17
const MatrixDimI LambdaA
Shift tensor describes the shift in the CSL when lattice is shifted by a DSCL vector.
Definition BiCrystal.h:103
Definition Gb.h:16
Lattice class.
Definition Lattice.h:34
const MatrixDimD latticeBasis
Definition Lattice.h:46
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:169
LatticeVector class.
VectorDimD cartesian() const
const Lattice< dim > & lattice
static Eigen::Matrix< IntScalarType, dim, dim > adjoint(const Eigen::Matrix< IntScalarType, dim, dim > &input)
const MatrixType & reducedBasis() const
Definition RLLL.cpp:192
const IntScalarType & mu
std::enable_if< dm==2, LatticeDirection< dim > >::type cross(const ReciprocalLatticeVector< dim > &other) const
const MatrixNi & matrixD() const
const MatrixNi & matrixV() const
const MatrixNi & matrixX() const
static IntScalarType gcd(const IntScalarType &a, const IntScalarType &b)
Definition IntegerMath.h:33
static void solveDiophantine2vars(IntScalarType a, IntScalarType b, IntScalarType c, IntScalarType &x, IntScalarType &y)
LatticeDirection class.
int stacking() const
Returns the number of planes in the stacking sequence.
const ReciprocalLatticeVector< dim > & reciprocalLatticeVector() const
Returns a constant reference to the base class (ReciprocalLatticeVector)