Alamo
CrystalPlastic.H
Go to the documentation of this file.
1//
2// A basic and relatively untested implementation of finite-deformation FCC crystal plasticity.
3// Use with caution.
4//
5// Inherits from :ref:`Model::Solid::Finite::PseudoLinear::Cubic` to provide elastic response.
6//
7// Plastic flow modeled as:
8//
9// .. math::
10//
11// \mathbf{F}^p\mathbf{F}^{p-1} =
12// \sum \dot{\gamma}_n \mathbf{a}_n\otimes\mathbf{N}_n
13//
14// where :math:`\mathbf{F}=\mathbf{F}^e\mathbf{F}^p` and :math:`\gamma_n` are
15// slips on system :math:`n`.
16//
17// Power law rate hardening is used in the integration of slips:
18//
19// .. math::
20//
21// \dot{\gamma}_n = \dot{\gamma}_0 \Big(\frac{\tau}{\tau_{crss}}\Big)^m \operatorname{sign}(\tau)
22//
23//
24
25#ifndef MODEL_SOLID_FINITE_CRYSTALPLASTIC_H_
26#define MODEL_SOLID_FINITE_CRYSTALPLASTIC_H_
27
28
29#include "IO/ParmParse.H"
30#include "Model/Solid/Solid.H"
31#include "PseudoLinear/Cubic.H"
32
33namespace Model
34{
35namespace Solid
36{
37namespace Finite
38{
40{
41public:
43
44 virtual ~CrystalPlastic() {};
45
46 Set::Scalar W(const Set::Matrix & F) const override
47 {
48 return PseudoLinear::Cubic::W(F * Set::reduce(Fp).inverse());
49 }
50 Set::Matrix DW(const Set::Matrix & F) const override
51 {
52 Set::Matrix Fpinv = Set::reduce(Fp).inverse();
53 return PseudoLinear::Cubic::DW(F * Fpinv) * Fpinv.transpose();
54 }
56 {
57 Set::Matrix Fpinv = Set::reduce(Fp).inverse();
58 auto DDW = PseudoLinear::Cubic::DDW(F * Fpinv);
60 for (int i = 0; i < AMREX_SPACEDIM; i++)
61 for (int J = 0; J < AMREX_SPACEDIM; J++)
62 for (int k = 0; k < AMREX_SPACEDIM; k++)
63 for (int L = 0; L < AMREX_SPACEDIM; L++)
64 {
65 ret(i,J,k,L) = 0.0;
66 for (int Q = 0; Q < AMREX_SPACEDIM; Q++)
67 for (int R = 0; R < AMREX_SPACEDIM; R++)
68 ret(i,J,k,L) += DDW(i,Q,k,R)*Fpinv(J,Q)*Fpinv(L,R);
69 }
70 return ret;
71 }
72
73 virtual void Advance(Set::Scalar dt, Set::Matrix /*F*/, Set::Matrix P, Set::Scalar time) override
74 {
75 if (time < tstart) return;
76
77 q.normalize();
78 Set::Matrix3d R = q.toRotationMatrix();
79 std::array<std::pair<Set::Vector3d,Set::Vector3d>,12> ss = slipSystems(R);
80
81 Set::Matrix3d L = Set::Matrix3d::Zero();
82
83 for (int n = 0; n < 12; n++)
84 {
85 Set::Vector3d A = ss[n].first;
86 Set::Vector3d N = ss[n].second;
87 Set::Scalar tau = A.dot(Set::expand(P)*N);
88 Set::Scalar sign = tau>0 ? 1.0 : -1.0;
89
90 Set::Scalar gammadot
91 = gammadot0 * std::pow(fabs(tau/tau_crss[n]),m_rate_inv) * sign;
92
93 L += gammadot * A * N.transpose();
94
95 gamma[n] += gammadot * dt;
96 }
97
98 Fp += L * Fp * dt;
99 }
100
101
102public:
103 //Set::Matrix4<AMREX_SPACEDIM,Set::Sym::MajorMinor> CC;
104 //static constexpr KinematicVariable kinvar = KinematicVariable::F;
105 //Set::Scalar C11=1.68, C12=1.21, C44=0.75;
106 //mutable Set::Quaternion q = Set::Quaternion(1.0, 0.0, 0.0, 0.0);
107
108 // Plastic state variables
109 Set::Matrix3d Fp = Set::Matrix3d::Identity();
110 Eigen::Matrix<Set::Scalar,12,1> gamma = Eigen::Matrix<Set::Scalar,12,1>::Zero();
111
112 // Parameters
114 //Set::Scalar tau_crss = NAN;
115 Eigen::Matrix<Set::Scalar,12,1> tau_crss = Eigen::Matrix<Set::Scalar,12,1>::Zero();
117
118 // Control variables
120
121private:
122 static std::array<std::pair<Set::Vector3d,Set::Vector3d>,12> slipSystems(Set::Matrix3d R)
123 {
124 const std::array<Set::Vector3d,12> n0 =
125 {
126 Set::Vector3d::Constant( 1, 1, 1) / sqrt(3.0),
127 Set::Vector3d::Constant( 1, 1, 1) / sqrt(3.0),
128 Set::Vector3d::Constant( 1, 1, 1) / sqrt(3.0),
129 Set::Vector3d::Constant(-1, 1,-1) / sqrt(3.0),
130 Set::Vector3d::Constant(-1, 1,-1) / sqrt(3.0),
131 Set::Vector3d::Constant(-1, 1,-1) / sqrt(3.0),
132 Set::Vector3d::Constant( 1,-1,-1) / sqrt(3.0),
133 Set::Vector3d::Constant( 1,-1,-1) / sqrt(3.0),
134 Set::Vector3d::Constant( 1,-1,-1) / sqrt(3.0),
135 Set::Vector3d::Constant(-1,-1, 1) / sqrt(3.0),
136 Set::Vector3d::Constant(-1,-1, 1) / sqrt(3.0),
137 Set::Vector3d::Constant(-1,-1, 1) / sqrt(3.0)
138 };
139 const std::array<Set::Vector3d,12> a0 =
140 {
141 Set::Vector3d::Constant( 1, 0,-1) / sqrt(2.0), // 1 1 1
142 Set::Vector3d::Constant( 0,-1, 1) / sqrt(2.0), // 1 1 1
143 Set::Vector3d::Constant( 1,-1, 0) / sqrt(2.0), // 1 1 1
144 Set::Vector3d::Constant( 1, 0,-1) / sqrt(2.0), // -1 1 -1
145 Set::Vector3d::Constant( 1, 1, 0) / sqrt(2.0), // -1 1 -1
146 Set::Vector3d::Constant( 0, 1, 1) / sqrt(2.0), // -1 1 -1
147 Set::Vector3d::Constant( 1, 1, 0) / sqrt(2.0), // 1 -1 -1
148 Set::Vector3d::Constant( 0,-1, 1) / sqrt(2.0), // 1 -1 -1
149 Set::Vector3d::Constant( 1, 0, 1) / sqrt(2.0), // 1 -1 -1
150 Set::Vector3d::Constant( 0, 1, 1) / sqrt(2.0), // -1 -1 1
151 Set::Vector3d::Constant( 1, 0, 1) / sqrt(2.0), // -1 -1 1
152 Set::Vector3d::Constant( 1,-1, 0) / sqrt(2.0) // -1 -1 1
153 };
154
155
156 static std::array<std::pair<Set::Vector3d,Set::Vector3d>,12> ret;
157 for (int n = 0 ; n < 12 ; n ++)
158 {
159 ret[n].first = R*a0[n];
160 ret[n].second = R*n0[n];
161 }
162 return ret;
163 }
164
165public:
167 {
168 CrystalPlastic ret;
169 ret.C11 = 0.0;
170 ret.C12 = 0.0;
171 ret.C44 = 0.0;
172 ret.q = Set::Quaternion(0.0,0.0,0.0,0.0);
173 ret.Fp = Set::Matrix3d::Zero();
174 ret.gamma = Eigen::Matrix<Set::Scalar,12,1>::Zero();
175 ret.m_rate_inv = 0.0;
176 ret.tau_crss = Eigen::Matrix<Set::Scalar,12,1>::Zero();
177 ret.gammadot0 = 0.0;
178 return ret;
179 }
181 {
183 }
185 {
187 ret.C11 = C11;
188 ret.C12 = C12;
189 ret.C44 = C44;
190 ret.q = Set::Quaternion::UnitRandom();
191 return ret;
192 }
193 static void Parse(CrystalPlastic & value, IO::ParmParse & pp)
194 {
196 // TODO Add inputs for the CP values
197
198 std::vector<Set::Scalar> tau_crss;
199 // Critical resolved shear stress :math:`\tau_{crss}`
200 pp_queryarr("tau_crss",tau_crss);
201 if (tau_crss.size() == 1)
202 for (int n = 0; n < 12; n++) value.tau_crss[n] = tau_crss[0];
203 else if (tau_crss.size() == 12)
204 for (int n = 0; n < 12; n++) value.tau_crss[n] = tau_crss[n];
205 else
206 Util::Exception(INFO,"Expected either 1 or 12 values for tau_crss but got ",tau_crss.size());
207
208 // Rate hardening coefficient :math:`\dot{\gamma}_0`
209 pp_query_default("gammadot0",value.gammadot0,1.0);
210
211 // Inverse of the hardening exponent :math:`\frac{1}{m}`
212 pp_query_default("m_rate_inv",value.m_rate_inv,0.5);
213
214 // Time to activate plastic slip
215 pp_query_default("tstart",value.tstart,0.0);
216 }
217
218#define OP_CLASS CrystalPlastic
219#define OP_VARS X(C11) X(C12) X(C44) X(q) X(Fp) X(gamma) X(m_rate_inv) X(tau_crss) X(gammadot0) X(tstart)
221};
223
224
225
226
227
228}
229}
230}
231
232
233template<>
235{
236 return 12;
237}
238
239template<>
241{
242 if (i==0) return name + "gamma1";
243 if (i==1) return name + "gamma2";
244 if (i==2) return name + "gamma3";
245 if (i==3) return name + "gamma4";
246 if (i==4) return name + "gamma5";
247 if (i==5) return name + "gamma6";
248 if (i==6) return name + "gamma7";
249 if (i==7) return name + "gamma8";
250 if (i==8) return name + "gamma9";
251 if (i==9) return name + "gamma10";
252 if (i==10) return name + "gamma11";
253 if (i==11) return name + "gamma12";
254 return name;
255}
256
257template<>
259{
260 for (amrex::MFIter mfi(a_dst, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
261 {
262 const amrex::Box& bx = mfi.growntilebox(amrex::IntVect(a_nghost));
263 if (bx.ok())
264 {
265 amrex::Array4<const Model::Solid::Finite::CrystalPlastic> const & src = ((*this)[a_lev])->array(mfi);
266 amrex::Array4<Set::Scalar> const & dst = a_dst.array(mfi);
267 for (int n = 0; n < 12; n++)
268 {
269 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
270 dst(i,j,k,a_dstcomp + n) = src(i,j,k).gamma(n);
271 });
272 }
273 }
274 }
275}
276
277
278
279#endif
#define pp_queryarr(...)
Definition ParmParse.H:103
#define pp_query_default(...)
Definition ParmParse.H:100
#define INFO
Definition Util.H:20
virtual void Advance(Set::Scalar dt, Set::Matrix, Set::Matrix P, Set::Scalar time) override
Set::Scalar W(const Set::Matrix &F) const override
static std::array< std::pair< Set::Vector3d, Set::Vector3d >, 12 > slipSystems(Set::Matrix3d R)
Eigen::Matrix< Set::Scalar, 12, 1 > tau_crss
static void Parse(CrystalPlastic &value, IO::ParmParse &pp)
Set::Matrix DW(const Set::Matrix &F) const override
Set::Matrix4< AMREX_SPACEDIM, Set::Sym::Major > DDW(const Set::Matrix &F) const override
static CrystalPlastic Random(Set::Scalar C11, Set::Scalar C12, Set::Scalar C44)
Eigen::Matrix< Set::Scalar, 12, 1 > gamma
Set::Scalar W(const Set::Matrix &F) const override
Definition Cubic.H:22
Set::Matrix DW(const Set::Matrix &F) const override
Definition Cubic.H:30
Set::Matrix4< AMREX_SPACEDIM, Set::Sym::Major > DDW(const Set::Matrix &F) const override
Definition Cubic.H:47
static void Parse(Cubic &value, IO::ParmParse &pp)
Definition Cubic.H:115
std::string Name(int) const
Definition Set.H:73
void Copy(int, amrex::MultiFab &, int, int) const
Definition Set.H:68
int NComp() const
Definition Set.H:72
amrex::Real Scalar
Definition Base.H:19
Eigen::Matrix3d Matrix3d
Definition Base.H:24
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition Base.H:20
Eigen::Vector3d Vector3d
Definition Base.H:21
AMREX_FORCE_INLINE Set::Matrix3d expand(const Set::Matrix &in)
Definition Base.H:34
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, AMREX_SPACEDIM > Matrix
Definition Base.H:23
AMREX_FORCE_INLINE Set::Matrix reduce(const Set::Matrix3d &in)
Definition Base.H:28
Set::Scalar Random()
Definition Set.cpp:9
void Exception(std::string file, std::string func, int line, Args const &... args)
Definition Util.H:205