oILAB
Loading...
Searching...
No Matches
GbContinuum.h
Go to the documentation of this file.
1//
2// Created by Nikhil Chandra Admal on 5/16/24.
3//
4
5#ifndef OILAB_GBPLASTICITY_H
6#define OILAB_GBPLASTICITY_H
7
8#include "LatticeCore.h"
9#include "Eigen/Dense"
10#include "../Math/PeriodicFunction.h"
11#include "LatticeFunction.h"
12#include "../Math/Function.h"
13#include "GbMaterialTensors.h"
14#include "OrderedTuplet.h"
15
16namespace oILAB {
17
18class DisplacementKernel : public Function<DisplacementKernel, double> {
19private:
20 Eigen::Vector<double, Eigen::Dynamic> normal;
21
22public:
23 explicit DisplacementKernel(
24 const Eigen::Vector<double, Eigen::Dynamic> &_normal);
25 double operator()(const Eigen::Vector<double, Eigen::Dynamic> &x) const;
26 };
27
28 class ShiftedDisplacementKernelFT : public Function<ShiftedDisplacementKernelFT, std::complex<double>>{
29 private:
30 Eigen::Vector<double, Eigen::Dynamic> x, normal;
31 public:
32 explicit ShiftedDisplacementKernelFT(const Eigen::VectorXd& x, const Eigen::VectorXd& normal);
33 std::complex<double> operator()(const Eigen::VectorXd& xi) const;
34 };
35
36 class HhatInvFunction: public Function<HhatInvFunction, std::complex<double>>, GbMaterialTensors
37 {
38 private:
39 const Eigen::VectorXd e1,e3;
40 public:
41 const int t,i;
42 explicit HhatInvFunction(const int& t, const int& i, const Eigen::MatrixXd& domain);
43 std::complex<double> operator()(const Eigen::VectorXd& xi) const;
44 };
45
46 template<int dim>
49 using FunctionFFTPair= typename std::pair<std::vector<PeriodicFunction<double,dim-1>>,
50 std::vector<LatticeFunction<std::complex<double>,dim-1>>>;
51 using GbLatticeFunctions= typename std::vector<LatticeFunction<std::complex<double>,dim-1>> ;
52
53 private:
54
56 //static FunctionFFTPair pipihat;
57 static thread_local std::map<OrderedTuplet<dim+1>,PeriodicFunction<double, dim - 1>> piPeriodicFunctions;
58 static thread_local std::map<OrderedTuplet<dim+1>,LatticeFunction<std::complex<double>, dim - 1>> pihatLatticeFunctions;
59 // GBMesostateEnsemble should generate the bicrystal (member variable <OrderedTuplet,VectorDimD>) and pass it as a reference to each mesostate
60 // pipihat should be map from OrderedTuplet to FunctionFFTPair. should be computed once in calculateb
61 // at the same time, compute pipihat once
62 // change xuPairs type to <Tiplet,VectorDimD>
63
65 static FunctionFFTPair calculateb(const Eigen::Matrix<double, dim,dim-1>& domain,
67 const std::array<Eigen::Index,dim-1>& n,
68 const std::map<OrderedTuplet<dim+1>,VectorDimD>& points);
69 static GbLatticeFunctions getHhatInvComponents(const Eigen::Matrix<double, dim,dim-1>& domain,
70 const std::array<Eigen::Index,dim-1>& n);
71 static PeriodicFunction<double,dim-1>get_pi(const Eigen::Matrix<double,dim,dim-1>& domain,
72 const std::array<Eigen::Index,dim-1>& n,
73 const VectorDimD& point);
74 static LatticeFunction<std::complex<double>,dim-1>get_pihat(const Eigen::Matrix<double,dim,dim-1>& domain,
75 const std::array<Eigen::Index,dim-1>& n,
76 const VectorDimD& point);
77
78 public:
79
80
81 const Eigen::Matrix<double,dim,dim-1> gbDomain;
82 const std::map<OrderedTuplet<dim+1>,VectorDimD> xuPairs;
83 std::array<Eigen::Index,dim-1> n;
84 std::vector<PeriodicFunction<double,dim-1>> b;
85
86 std::vector<LatticeFunction<std::complex<double>,dim-1>> bhat;
87 std::map<OrderedTuplet<dim+1>,VectorDimD> atoms;
89
90 double energy;
91
92 GbContinuum(const Eigen::Matrix<double, dim,dim-1>& domain,
93 const std::map<OrderedTuplet<dim+1>, VectorDimD>& xuPairs,
94 const std::array<Eigen::Index,dim-1>& n,
95 const std::map<OrderedTuplet<dim+1>,VectorDimD>& atoms,
96 const bool& verbosity=false);
97
99
100 VectorDimD displacement(const VectorDimD& x) const;
101
102
107 std::vector<PeriodicFunction<double,dim-1>> get_b() const;
108
114 std::vector<PeriodicFunction<double,dim-1>> get_alpha() const;
115
116 static void reset(){
117 std::map<OrderedTuplet<dim+1>,PeriodicFunction<double, dim - 1>>().swap(piPeriodicFunctions);
118 std::map<OrderedTuplet<dim+1>,LatticeFunction<std::complex<double>, dim - 1>>().swap(pihatLatticeFunctions);
120
121 //HhatInvComponents.clear();
122 //piPeriodicFunctions.clear();
123 //pihatLatticeFunctions.clear();
124 }
125 };
126
127
128 template<int dim>
129 thread_local std::vector<LatticeFunction<std::complex<double>,dim-1>> GbContinuum<dim>::HhatInvComponents;
130
131 template<int dim>
132 thread_local std::map<OrderedTuplet<dim+1>,PeriodicFunction<double, dim - 1>> GbContinuum<dim>::piPeriodicFunctions;
133
134 template<int dim>
135 thread_local std::map<OrderedTuplet<dim+1>,LatticeFunction<std::complex<double>, dim - 1>> GbContinuum<dim>::pihatLatticeFunctions;
136
137 } // namespace oILAB
138
140#endif //OILAB_GBPLASTICITY_H
double operator()(const Eigen::Vector< double, Eigen::Dynamic > &x) const
Eigen::Vector< double, Eigen::Dynamic > normal
Definition GbContinuum.h:20
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)
VectorDimD uAverage
Definition GbContinuum.h:88
std::vector< LatticeFunction< std::complex< double >, dim-1 > > bhat
Definition GbContinuum.h:86
const Eigen::Matrix< double, dim, dim-1 > gbDomain
Definition GbContinuum.h:81
static thread_local GbLatticeFunctions HhatInvComponents
Definition GbContinuum.h:55
typename std::vector< LatticeFunction< std::complex< double >, dim-1 > > GbLatticeFunctions
Definition GbContinuum.h:51
std::vector< PeriodicFunction< double, dim-1 > > b
Definition GbContinuum.h:84
static thread_local std::map< OrderedTuplet< dim+1 >, PeriodicFunction< double, dim - 1 > > piPeriodicFunctions
Definition GbContinuum.h:57
typename LatticeCore< dim >::VectorDimD VectorDimD
Definition GbContinuum.h:48
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)
typename std::pair< std::vector< PeriodicFunction< double, dim-1 > >, std::vector< LatticeFunction< std::complex< double >, dim-1 > > > FunctionFFTPair
Definition GbContinuum.h:50
static thread_local std::map< OrderedTuplet< dim+1 >, LatticeFunction< std::complex< double >, dim - 1 > > pihatLatticeFunctions
Definition GbContinuum.h:58
std::map< OrderedTuplet< dim+1 >, VectorDimD > atoms
Definition GbContinuum.h:87
FunctionFFTPair bbhat
Definition GbContinuum.h:64
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
std::array< Eigen::Index, dim-1 > n
Definition GbContinuum.h:83
std::vector< PeriodicFunction< double, dim-1 > > get_alpha() const
VectorDimD displacement(const OrderedTuplet< dim+1 > &x) const
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)
const std::map< OrderedTuplet< dim+1 >, VectorDimD > xuPairs
Definition GbContinuum.h:82
static void reset()
const Eigen::VectorXd e1
Definition GbContinuum.h:39
const Eigen::VectorXd e3
Definition GbContinuum.h:39
std::complex< double > operator()(const Eigen::VectorXd &xi) const
Eigen::Vector< double, Eigen::Dynamic > normal
Definition GbContinuum.h:30
Eigen::Vector< double, Eigen::Dynamic > x
Definition GbContinuum.h:30
std::complex< double > operator()(const Eigen::VectorXd &xi) const
Eigen::Matrix< double, dim, 1 > VectorDimD
Definition LatticeCore.h:20