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