oILAB
Loading...
Searching...
No Matches
GbContinuumImplementation.h
Go to the documentation of this file.
1//
2// Created by Nikhil Chandra Admal on 5/27/24.
3//
4
5#ifndef OILAB_GBCONTINUUMIMPLEMENTATION_H
6#define OILAB_GBCONTINUUMIMPLEMENTATION_H
7
8#include<Diff.h>
9#include <numbers>
10
11namespace gbLAB {
12
13 template<int dim>
14 GbContinuum<dim>::GbContinuum(const Eigen::Matrix<double,dim,dim-1>& domain,
15 const std::map<OrderedTuplet<dim+1>,VectorDimD>& xuPairs,
16 const std::array<Eigen::Index,dim-1>& n,
17 const std::map<OrderedTuplet<dim+1>,VectorDimD>& atoms,
18 const bool& verbosity):
19 gbDomain(domain),
20 xuPairs(xuPairs),
21 n(n),
22 bbhat(calculateb(domain,xuPairs,n,atoms)),
23 b(bbhat.first),
24 bhat(bbhat.second)
25 {
26 uAverage.setZero();
27 for (const auto &[x, u]: xuPairs) {
28 uAverage = uAverage + u;
29 }
30 if (xuPairs.size() != 0)
31 uAverage = uAverage / xuPairs.size();
32
33 if (verbosity) {
34 std::cout << "------------------------------------------------------------------------------" << std::endl;
35 std::cout << "Constraints: " << std::endl;
36 }
37 for (const auto &[x, u]: xuPairs) {
38 if (verbosity)
39 std::cout << "x = " << atoms.at(x).transpose() << "; displacement = " << u.transpose() << std::endl;
40 if ((u - displacement(x) - uAverage).norm() > FLT_EPSILON)
41 throw std::runtime_error("GBContinuum construction failed - unable to impose constraints.");
42 }
43 if (verbosity)
44 std::cout << std::endl;
45
46 }
47
48 template<int dim>
49 std::vector<PeriodicFunction<double,dim-1>>
51 {
52 auto bLocal= get_b();
53 auto unitcellLocal= bLocal[0].unitCell;
54
55 // dimension-dependent part
56 Diff<dim-1> dx({1,0},unitcellLocal,n);
57 Diff<dim-1> dy({0,1},unitcellLocal,n);
58
59 std::vector<PeriodicFunction<double,dim-1>> alpha;
60 for(int i=0; i<dim; ++i)
61 {
62 for(int j=0; j<dim-1; ++j)
63 {
64 PeriodicFunction<double,dim-1> alphaij(n, unitcellLocal);;
65 for(int k=0; k<dim-1; ++k)
66 {
67 if (k==j) continue;
68 PeriodicFunction<double,dim-1> dlbi_dxk(n,unitcellLocal);;
69 if (k==0)
70 dx.perform_op(bLocal[i].values.data(), dlbi_dxk.values.data());
71 else if (k==1)
72 dy.perform_op(bLocal[i].values.data(), dlbi_dxk.values.data());
73 if (j==0 && k==1)
74 alphaij.values= alphaij.values + dlbi_dxk.values;
75 else if (j==1 && k==0)
76 alphaij.values= alphaij.values - dlbi_dxk.values;
77 }
78 alpha.push_back(alphaij);
79 }
80 }
81 return alpha;
82 }
83
84 template<int dim>
85 std::vector<PeriodicFunction<double,dim-1>>
87 {
88 Eigen::Matrix<double,dim,dim-1> orthogonalVectors;
89 Eigen::Matrix<double,dim,dim-1> unitCell= b[0].unitCell;
90
91 // Form a rotation matrix to transform to a 2D coordinate system containing the grain boundary
92 for(int i=0; i<dim-1; ++i)
93 {
94 orthogonalVectors.col(i)= unitCell.col(i);
95 for(int j=0; j<i; ++j)
96 orthogonalVectors.col(i)-= unitCell.col(i).dot(orthogonalVectors.col(j))* orthogonalVectors.col(j);
97 orthogonalVectors.col(i)= orthogonalVectors.col(i).normalized();
98 }
99 Rotation<dim> rotation(orthogonalVectors);
100 assert((rotation*unitCell).row(2).norm()<FLT_EPSILON);
101
102 // transform the slip vectors b to the local coordinate system
103 std::vector<PeriodicFunction<double, dim - 1>> bLocal;
104 Eigen::Matrix<double,dim-1,dim-1> unitcellLocal= (rotation*unitCell).block(0,0,dim-1,dim-1);
105 for(int i= 0; i<dim; ++i) {
106 PeriodicFunction<double, dim - 1> temp(n,unitcellLocal);
107 for (int j = 0; j < dim; ++j)
108 temp.values = temp.values + rotation(i, j) * b[j].values;
109 bLocal.push_back(temp);
110 }
111 return bLocal;
112 }
113
114 template<int dim>
115 PeriodicFunction<double,dim-1>
116 GbContinuum<dim>::get_pi(const Eigen::Matrix<double,dim,dim-1>& domain,
117 const std::array<Eigen::Index,dim-1>& n,
118 const VectorDimD& point)
119 {
120 Eigen::Matrix<double,dim,dim-1> basisVectors(domain.transpose().completeOrthogonalDecomposition().pseudoInverse());
121 VectorDimD normal(domain.col(0).cross(domain.col(1)));
122 normal.normalize();
123
124 //DisplacementKernel f(normal,10);
125 // change 1e6 entry
126 DisplacementKernel f(normal);
128 return PeriodicFunction<double,dim-1>(n,domain,piFunction);
129 }
130
131
132 template<int dim>
134 GbContinuum<dim>::get_pihat(const Eigen::Matrix<double,dim,dim-1>& domain,
135 const std::array<Eigen::Index,dim-1>& n,
136 const VectorDimD& point)
137 {
138 Eigen::Matrix<double,dim,dim-1> basisVectors(domain.transpose().completeOrthogonalDecomposition().pseudoInverse());
139 VectorDimD normal(domain.col(0).cross(domain.col(1)));
140 normal.normalize();
141
142 std::vector<LatticeFunction<std::complex<double>,dim-1>> pihat;
143 //DisplacementKernel f(normal,10);
144
145 ShiftedDisplacementKernelFT pihatFunction(point,normal);
146 return LatticeFunction<std::complex<double>,dim-1> (n,basisVectors,pihatFunction);
147 }
148
149 template<int dim>
151 GbContinuum<dim>::calculateb(const Eigen::Matrix<double,dim,dim-1>& domain,
152 const std::map<OrderedTuplet<dim+1>,VectorDimD>& xuPairs,
153 const std::array<Eigen::Index,dim-1>& n,
154 const std::map<OrderedTuplet<dim+1>,VectorDimD>& atoms)
155 {
156 if (HhatInvComponents.size() ==0 ) HhatInvComponents= getHhatInvComponents(domain,n);
157 Eigen::Matrix<double,dim,dim-1> basisVectors(domain.transpose().completeOrthogonalDecomposition().pseudoInverse());
158 // lb0 - read as "local b0"
159 std::vector<PeriodicFunction<double,dim-1>> lb0, lb;
160 std::vector<LatticeFunction<std::complex<double>,dim-1>> lb0hat;
161 for(int i=0; i<dim; ++i)
162 {
163 lb0.push_back(PeriodicFunction<double,dim-1>(n,domain));
164 lb.push_back(PeriodicFunction<double,dim-1>(n,domain));
165 lb0hat.push_back(LatticeFunction<std::complex<double>,dim-1>(n,basisVectors));
166 }
167 std::vector<LatticeFunction<std::complex<double>,dim-1>> lbhat(lb0hat);
168
169
170 if(piPeriodicFunctions.empty()) {
171 for (const auto& [key, value]: atoms) {
172 // the cross product of the domain vectors has to be parallel to nA
173 Eigen::Matrix<double,dim,dim-1> basisVectors(domain.transpose().completeOrthogonalDecomposition().pseudoInverse());
174 VectorDimD normal(domain.col(0).cross(domain.col(1)));
175 normal.normalize();
176
177 VectorDimD perturbedValue= value;
178 if(abs(value.dot(normal)) < FLT_EPSILON) // belongs to the GB
179 {
180 if (key(dim)==1)
181 perturbedValue= value - (value.dot(normal) + FLT_EPSILON) * normal;
182 else if (key(dim)==2)
183 perturbedValue= value - (value.dot(normal) - FLT_EPSILON) * normal;
184 }
185
186 else if (key(dim)==-1 || key(dim)==-2) // belongs to (lattice 1 & grain 2) || (lattice 2 & grain 1)
187 {
188 perturbedValue= value - 2 * value.dot(normal) * normal;
189 }
190 piPeriodicFunctions.insert({key, get_pi(domain, n, perturbedValue)});
191 pihatLatticeFunctions.insert({key, get_pihat(domain, n, perturbedValue)});
192
193 /*
194 if(abs(value.dot(normal)) < FLT_EPSILON && key(dim)==1) // belongs to lattice 1 and on the GB
195 {
196 VectorDimD perturbedValue= value - (value.dot(normal) + FLT_EPSILON) * normal;
197 piPeriodicFunctions.insert({key, get_pi(domain, n, perturbedValue)});
198 pihatLatticeFunctions.insert({key, get_pihat(domain, n, perturbedValue)});
199 }
200 else if(abs(value.dot(normal)) < FLT_EPSILON && key(dim)==2) // belongs to lattice 2 and on the GB
201 {
202 VectorDimD perturbedValue= value - (value.dot(normal) - FLT_EPSILON) * normal;
203 piPeriodicFunctions.insert({key, get_pi(domain, n, perturbedValue)});
204 pihatLatticeFunctions.insert({key, get_pihat(domain, n, perturbedValue)});
205 }
206 else {
207 piPeriodicFunctions.insert({key, get_pi(domain, n, value)});
208 pihatLatticeFunctions.insert({key, get_pihat(domain, n, value)});
209 }
210 */
211 }
212 }
213
214 // matrices P and u
215 const int numberOfCslPoints= xuPairs.size();
216 Eigen::Matrix<std::complex<double>,Eigen::Dynamic,dim> uMatrix(numberOfCslPoints,dim);
217 Eigen::Matrix<std::complex<double>,Eigen::Dynamic,Eigen::Dynamic> P(numberOfCslPoints,numberOfCslPoints);
218
219
220 VectorDimD uAverage;
221 uAverage.setZero();
222 for(const auto& [xi,ui] : xuPairs)
223 {
224 uAverage= uAverage + ui;
225 }
226 if (xuPairs.size() != 0)
227 uAverage= uAverage/xuPairs.size();
228
229 int row= -1;
230 for(const auto& [xi,ui] : xuPairs)
231 {
232 row++;
233 //uMatrix.row(row)= ui;
234 uMatrix.row(row)= ui-uAverage;
235 int col= -1;
236 for(const auto& [xj,uj] : xuPairs)
237 {
238 col++;
239 P(row, col) = pihatLatticeFunctions.at(xj).dot(pihatLatticeFunctions.at(xi));
240 }
241
242 }
243
244 // Solve P \alpha = u
245 Eigen::Matrix<std::complex<double>,Eigen::Dynamic,dim> alpha(numberOfCslPoints,dim);
246 Eigen::LLT<Eigen::Matrix<std::complex<double>,Eigen::Dynamic,Eigen::Dynamic>> llt;
247 llt.compute(P);
248 alpha= llt.solve(uMatrix);
249
250 // calculate lb0 and lb0hat
251
252 for(int i=0; i<dim; ++i)
253 {
254 lb0hat[i].values.setZero();
255 int j= -1;
256 for(const auto& [xj,uj] : xuPairs) {
257 j++;
258 lb0hat[i].values = lb0hat[i].values + pihatLatticeFunctions.at(xj).values * alpha(j,i);
259 }
260 }
261 for(int i=0; i<dim; ++i) {
262 lb0[i].values.setZero();
263 int j = -1;
264 for (const auto &[xj, uj]: xuPairs) {
265 j++;
266 lb0[i].values = lb0[i].values + piPeriodicFunctions.at(xj).values * alpha(j,i).real();
267 }
268 }
269
270 // calculate lagrange multipliers, lbhat, and lb
271 Eigen::VectorXcd uFlattened;
272 uFlattened= uMatrix.reshaped();
273 Eigen::VectorXcd lm(dim*numberOfCslPoints);
274 Eigen::MatrixXcd M(dim*numberOfCslPoints,dim*numberOfCslPoints);
275 M.setZero();
276
277 for (int i=0; i<dim; ++i)
278 {
279 int j= -1;
280 for (const auto& [xj,uj] : xuPairs)
281 {
282 ++j;
283 int row= i*numberOfCslPoints + j;
284 for (int l=0; l<dim; ++l)
285 {
286 int k= -1;
287 for(const auto& [xk,uk] : xuPairs)
288 {
289 k++;
290 int col= l*numberOfCslPoints + k;
291
292 std::complex<double> temp;
293
294 if (i==0 && l==0) temp = (HhatInvComponents[0]*pihatLatticeFunctions.at(xk)).dot(pihatLatticeFunctions.at(xj));
295 if (i==1 && l==1) temp = (HhatInvComponents[1]*pihatLatticeFunctions.at(xk)).dot(pihatLatticeFunctions.at(xj));
296 if (i==2 && l==2) temp = (HhatInvComponents[2]*pihatLatticeFunctions.at(xk)).dot(pihatLatticeFunctions.at(xj));
297 if ((i==1 && l==2) || (i==2 && l==1)) temp = (HhatInvComponents[3]*pihatLatticeFunctions.at(xk)).dot(pihatLatticeFunctions.at(xj));
298 if ((i==0 && l==2) || (i==2 && l==0)) temp = (HhatInvComponents[4]*pihatLatticeFunctions.at(xk)).dot(pihatLatticeFunctions.at(xj));
299 if ((i==0 && l==1) || (i==1 && l==0)) temp = (HhatInvComponents[5]*pihatLatticeFunctions.at(xk)).dot(pihatLatticeFunctions.at(xj));
300 M(row, col) = temp;
301 }
302 }
303 }
304 }
305
306 // Solve M lm = u2
307 Eigen::LLT<Eigen::Matrix<std::complex<double>,Eigen::Dynamic,Eigen::Dynamic>> llt2;
308 llt2.compute(M);
309 lm= llt2.solve(uFlattened);
310
311 lbhat[0].values.setZero();
312 lbhat[1].values.setZero();
313 lbhat[2].values.setZero();
314
315 auto lmMatrix= lm.reshaped(numberOfCslPoints,dim);
316 // temp = \sum_k p_k * lm_k
317 std::vector<LatticeFunction<std::complex<double>,dim-1>> temp;
318
319 for (int i=0; i<dim; ++i) {
320 temp.push_back(LatticeFunction<std::complex<double>, dim - 1>(n, basisVectors));
321 int k= -1;
322 for (const auto& [xk,uk] : xuPairs) {
323 k++;
324 temp[i].values = temp[i].values + pihatLatticeFunctions.at(xk).values * lmMatrix(k, i);
325 }
326 }
327
328 lbhat[0].values= HhatInvComponents[0].values * temp[0].values +
329 HhatInvComponents[5].values * temp[1].values +
330 HhatInvComponents[4].values * temp[2].values;
331
332 lbhat[1].values= HhatInvComponents[5].values * temp[0].values +
333 HhatInvComponents[1].values * temp[1].values +
334 HhatInvComponents[3].values * temp[2].values;
335
336 lbhat[2].values= HhatInvComponents[4].values * temp[0].values +
337 HhatInvComponents[3].values * temp[1].values +
338 HhatInvComponents[2].values * temp[2].values;
339
340 for (int i=0; i<dim; ++i)
341 lb[i].values = lbhat[i].ifft().values.real();
342
343 //return std::make_pair(lb,lbhat);
344 return std::make_pair(lb,lbhat);
345 }
346
347 template<int dim>
349 {
350 VectorDimD u;
351
352 // u = f \star b
353 for(int i=0; i<dim; ++i)
354 u(i)= (bhat[i].dot(pihatLatticeFunctions.at(t))).real();
355
356 return u;
357
358 }
359
360
361 template<int dim>
363 {
364 VectorDimD u;
365
366 Eigen::Matrix<double,dim,dim-1> basisVectors(gbDomain.transpose().completeOrthogonalDecomposition().pseudoInverse());
367 VectorDimD normal(gbDomain.col(0).cross(gbDomain.col(1)));
368 normal.normalize();
369
370
371 ShiftedDisplacementKernelFT pihatFunction(t,normal);
372 auto lf= LatticeFunction<std::complex<double>,dim-1> (n,basisVectors,pihatFunction);
373
374 // u = f \star b
375 for(int i=0; i<dim; ++i) {
376 u(i) = bhat[i].dot(lf).real();
377 //assert(abs(bhat[i].dot(lf).imag()) < FLT_EPSILON);
378 }
379
380
381
382 return u;
383
384 }
385
386 template<int dim>
387 std::vector<LatticeFunction<std::complex<double>,dim-1>>
388 GbContinuum<dim>::getHhatInvComponents(const Eigen::Matrix<double, dim,dim-1>& domain,
389 const std::array<Eigen::Index,dim-1>& n)
390 {
391 Eigen::Matrix<double,Eigen::Dynamic,dim-1> basisVectors(domain.transpose().completeOrthogonalDecomposition().pseudoInverse());
392 std::vector<LatticeFunction<std::complex<double>,dim-1>> output;
393
394 output.push_back(LatticeFunction<std::complex<double>,dim-1>(n,basisVectors,HhatInvFunction(0,0,domain)));
395 output.push_back(LatticeFunction<std::complex<double>,dim-1>(n,basisVectors,HhatInvFunction(1,1,domain)));
396 output.push_back(LatticeFunction<std::complex<double>,dim-1>(n,basisVectors,HhatInvFunction(2,2,domain)));
397 output.push_back(LatticeFunction<std::complex<double>,dim-1>(n,basisVectors,HhatInvFunction(1,2,domain)));
398 output.push_back(LatticeFunction<std::complex<double>,dim-1>(n,basisVectors,HhatInvFunction(0,2,domain)));
399 output.push_back(LatticeFunction<std::complex<double>,dim-1>(n,basisVectors,HhatInvFunction(0,1,domain)));
400 return output;
401 }
402
403
404 /*----------------------------*/
405 DisplacementKernel::DisplacementKernel(const Eigen::Vector<double, Eigen::Dynamic>& _normal):
406 Function<DisplacementKernel, double>(),
407 normal(_normal)
408 {}
409 double DisplacementKernel::operator()(const Eigen::Vector<double, Eigen::Dynamic> &x) const
410 {
411 return -x.dot(normal) / (4 * std::numbers::pi * (std::pow(x.norm(), 3)));
412 }
413 /*----------------------------*/
415 const Eigen::VectorXd& normal):
416 x(x),normal(normal.normalized())
417 { }
418
419 std::complex<double> ShiftedDisplacementKernelFT::operator()(const Eigen::VectorXd& xi) const
420 {
421 double xNormalComponent= x.dot(normal);
422 auto xProjected= (x-xNormalComponent*normal).eval();
423 double xiBar= xi.norm();
424 std::complex<double> output= -0.5*std::exp(-2*std::numbers::pi* std::complex<double>(0,1) *
425 xProjected.dot(xi)
426 ) * std::exp(-2*std::numbers::pi*xiBar*abs(xNormalComponent));
427 if (abs(xNormalComponent) < DBL_EPSILON)
428 return output;
429 else
430 return output*xNormalComponent/abs(xNormalComponent);
431 }
432 /*----------------------------*/
433 HhatInvFunction::HhatInvFunction(const int& t, const int& i, const Eigen::MatrixXd& domain) :
434 t(t),i(i), e1(domain.col(0).normalized()), e3(domain.col(1).normalized())
435 { }
436 std::complex<double> HhatInvFunction::operator()(const Eigen::VectorXd& xi) const
437 {
438 if (xi.isZero())
439 return std::complex<double>(0,0);
440 std::complex<double> output(0,0);
441 double xi1= xi.dot(e1);
442 double xi3= xi.dot(e3);
443 double xi2= (xi-xi1*e1-xi3*e3).norm();
444
445 output= tensorHhat(t,i,(Eigen::VectorXd(3) << xi1,xi2,xi3).finished());
446 output= -std::pow(2*std::numbers::pi,2) * output;
447 return output;
448 }
449 /*----------------------------*/
450
451}
452
453#endif
double operator()(const Eigen::Vector< double, Eigen::Dynamic > &x) const
DisplacementKernel(const Eigen::Vector< double, Eigen::Dynamic > &_normal)
Eigen::Vector< double, Eigen::Dynamic > normal
Definition GbContinuum.h:20
VectorDimD displacement(const OrderedTuplet< dim+1 > &x) const
VectorDimD uAverage
Definition GbContinuum.h:86
typename std::pair< std::vector< PeriodicFunction< double, dim-1 > >, std::vector< LatticeFunction< std::complex< double >, dim-1 > > > FunctionFFTPair
Definition GbContinuum.h:48
static PeriodicFunction< double, dim-1 > get_pi(const Eigen::Matrix< double, dim, dim-1 > &domain, const std::array< Eigen::Index, dim-1 > &n, const VectorDimD &point)
GbContinuum(const Eigen::Matrix< double, dim, dim-1 > &domain, const std::map< OrderedTuplet< dim+1 >, VectorDimD > &xuPairs, const std::array< Eigen::Index, dim-1 > &n, const std::map< OrderedTuplet< dim+1 >, VectorDimD > &atoms, const bool &verbosity=false)
std::map< OrderedTuplet< dim+1 >, VectorDimD > atoms
Definition GbContinuum.h:85
static GbLatticeFunctions getHhatInvComponents(const Eigen::Matrix< double, dim, dim-1 > &domain, const std::array< Eigen::Index, dim-1 > &n)
std::vector< PeriodicFunction< double, dim-1 > > get_b() const
typename LatticeCore< dim >::VectorDimD VectorDimD
Definition GbContinuum.h:46
static LatticeFunction< std::complex< double >, dim-1 > get_pihat(const Eigen::Matrix< double, dim, dim-1 > &domain, const std::array< Eigen::Index, dim-1 > &n, const VectorDimD &point)
const std::map< OrderedTuplet< dim+1 >, VectorDimD > xuPairs
Definition GbContinuum.h:80
static FunctionFFTPair calculateb(const Eigen::Matrix< double, dim, dim-1 > &domain, const std::map< OrderedTuplet< dim+1 >, VectorDimD > &xuPairs, const std::array< Eigen::Index, dim-1 > &n, const std::map< OrderedTuplet< dim+1 >, VectorDimD > &points)
std::vector< PeriodicFunction< double, dim-1 > > get_alpha() const
static std::complex< double > tensorHhat(const int &t, const int &i, const Eigen::VectorXd &xi)
const Eigen::VectorXd e3
Definition GbContinuum.h:37
HhatInvFunction(const int &t, const int &i, const Eigen::MatrixXd &domain)
std::complex< double > operator()(const Eigen::VectorXd &xi) const
const Eigen::VectorXd e1
Definition GbContinuum.h:37
std::complex< double > dot(const LatticeFunction< std::complex< double >, dim > &other) const
const Eigen::Matrix< double, Eigen::Dynamic, dim > unitCell
Eigen::Vector< double, Eigen::Dynamic > normal
Definition GbContinuum.h:28
ShiftedDisplacementKernelFT(const Eigen::VectorXd &x, const Eigen::VectorXd &normal)
std::complex< double > operator()(const Eigen::VectorXd &xi) const
Eigen::Vector< double, Eigen::Dynamic > x
Definition GbContinuum.h:28