-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathelastic2D_tools.py
198 lines (171 loc) · 8.5 KB
/
elastic2D_tools.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
import numpy as np
import matplotlib.pyplot as plt
def gen_strained_lats(lattice, f=0.02, N=7):
"""
Generates strained lattices for later post-processing into the second-order
elastic constants by fitting stress-strain data points to the generalized
Hooke's Law.
Args:
lattice (Numpy array) : lattice corresponding to a 2D material
(which includes a vvacuum layer)
f (float) : fraction of the lattice deformation
N (int) : number of stress-strain data points
Returns:
tuple of (Numpy array, list):
1. Numpy array of strain values
2. list of deformed lattices for 2D elastic constant calculations
"""
ds = np.linspace(-f, f, N)
eps_11_deformed_lattices = []
# apply eps_11 (for C1111, C2211, C1211)
for d in ds:
strain = np.zeros((3, 3))
strain[0, 0] = d
deformed_lattice = ((np.eye(3) + strain) @ lattice.T).T
eps_11_deformed_lattices.append(deformed_lattice)
eps_22_deformed_lattices = []
# apply eps_22 (for C1122, C2222, C1222)
for d in ds:
strain = np.zeros((3, 3))
strain[1, 1] = d
deformed_lattice = ((np.eye(3) + strain) @ lattice.T).T
eps_22_deformed_lattices.append(deformed_lattice)
eps_12_deformed_lattices = []
# apply eps_12 (for C1112, C2212, C1212)
for d in ds:
strain = np.zeros((3, 3))
strain[0, 1] = d / 2
strain[1, 0] = d / 2
deformed_lattice = ((np.eye(3) + strain) @ lattice.T).T
eps_12_deformed_lattices.append(deformed_lattice)
return ds, [eps_11_deformed_lattices, eps_22_deformed_lattices, eps_12_deformed_lattices]
def process_2D_elastic_constants(ds, stresses, c=None, symmetrize=True, plot=False):
"""
Fits the stress-strain data points to a quadratic polynomial and extracts
the linear coefficient as the elastic constants. Note that the units of the
returned elastic constants is the units of the stresses times the units of
the c variable (if given). The distinction between the Voigt and Mandel notation
of the returned elastic constants can be found in
Marcin Maździarz 2019 2D Mater.6 048001
(https://iopscience.iop.org/article/10.1088/2053-1583/ab2ef3).
Args:
ds (Numpy array) : strain array output by gen_strained_lats()
stresses (list) : list of stresses in the same format as the list of
deformed lattices from gen_strained_lats()
c (float) : height of the lattice (c-axis length)
symmetrize (bool): whether the elastic constants matrix is symmetrized
Returns:
tuple of Numpy arrays: (Cs_voigt, Cs_mandel, Cs_res)
Cs_voigt and Cs_mandel are the elastic constants in Voigt and Mandel (tensor)
notation while Cs_res shows the stress-strain fitting error.
"""
stresses_eps11, stresses_eps22, stresses_eps12 = [], [], []
# symmetrize stresses (in case)
for stress in stresses[0]:
stresses_eps11.append((stress + stress.T) / 2)
for stress in stresses[1]:
stresses_eps22.append((stress + stress.T) / 2)
for stress in stresses[2]:
stresses_eps12.append((stress + stress.T) / 2)
# for C1111, C2211, C1211
C1111_fit = np.polyfit(ds, np.array([stress[0, 0] for stress in stresses_eps11]), 2, full=True)
C1111 = C1111_fit[0][-2]
C1111_res = C1111_fit[1][0]
C2211_fit = np.polyfit(ds, np.array([stress[1, 1] for stress in stresses_eps11]), 2, full=True)
C2211 = C2211_fit[0][-2]
C2211_res = C2211_fit[1][0]
C1211_fit = np.polyfit(ds, np.array([stress[0, 1] for stress in stresses_eps11]), 2, full=True)
C1211 = C1211_fit[0][-2]
C1211_res = C1211_fit[1][0]
# for C1122, C2222, C1222
C1122_fit = np.polyfit(ds, np.array([stress[0, 0] for stress in stresses_eps22]), 2, full=True)
C1122 = C1122_fit[0][-2]
C1122_res = C1122_fit[1][0]
C2222_fit = np.polyfit(ds, np.array([stress[1, 1] for stress in stresses_eps22]), 2, full=True)
C2222 = C2222_fit[0][-2]
C2222_res = C2222_fit[1][0]
C1222_fit = np.polyfit(ds, np.array([stress[0, 1] for stress in stresses_eps22]), 2, full=True)
C1222 = C1222_fit[0][-2]
C1222_res = C1222_fit[1][0]
# for C1112, C2212, C1212
C1112_fit = np.polyfit(ds, np.array([stress[0, 0] for stress in stresses_eps12]), 2, full=True)
C1112 = C1112_fit[0][-2]
C1112_res = C1112_fit[1][0]
C2212_fit = np.polyfit(ds, np.array([stress[1, 1] for stress in stresses_eps12]), 2, full=True)
C2212 = C2212_fit[0][-2]
C2212_res = C2212_fit[1][0]
C1212_fit = np.polyfit(ds, np.array([stress[0, 1] for stress in stresses_eps12]), 2, full=True)
C1212 = C1212_fit[0][-2]
C1212_res = C1212_fit[1][0]
if plot:
fig, axs = plt.subplots(3, 3, sharex=True, gridspec_kw={'wspace': 0.5, 'hspace': 0})
axs[0, 0].plot(ds, np.array([stress[0, 0] for stress in stresses_eps11]), 'bx')
axs[0, 0].plot(ds, np.poly1d(C1111_fit[0])(ds), '--r')
axs[0, 1].plot(ds, np.array([stress[0, 0] for stress in stresses_eps22]), 'bx')
axs[0, 1].plot(ds, np.poly1d(C1122_fit[0])(ds), '--r')
axs[0, 2].plot(ds, np.array([stress[0, 0] for stress in stresses_eps12]), 'bx')
axs[0, 2].plot(ds, np.poly1d(C1112_fit[0])(ds), '--r')
axs[1, 0].plot(ds, np.array([stress[1, 1] for stress in stresses_eps11]), 'bx')
axs[1, 0].plot(ds, np.poly1d(C2211_fit[0])(ds), '--r')
axs[1, 1].plot(ds, np.array([stress[1, 1] for stress in stresses_eps22]), 'bx')
axs[1, 1].plot(ds, np.poly1d(C2222_fit[0])(ds), '--r')
axs[1, 2].plot(ds, np.array([stress[1, 1] for stress in stresses_eps12]), 'bx')
axs[1, 2].plot(ds, np.poly1d(C2212_fit[0])(ds), '--r')
axs[2, 0].plot(ds, np.array([stress[0, 1] for stress in stresses_eps11]), 'bx')
axs[2, 0].plot(ds, np.poly1d(C1211_fit[0])(ds), '--r')
axs[2, 1].plot(ds, np.array([stress[0, 1] for stress in stresses_eps22]), 'bx')
axs[2, 1].plot(ds, np.poly1d(C1222_fit[0])(ds), '--r')
axs[2, 2].plot(ds, np.array([stress[0, 1] for stress in stresses_eps12]), 'bx')
axs[2, 2].plot(ds, np.poly1d(C1212_fit[0])(ds), '--r')
plt.show()
# possible future TODO: additional symmetrization based on lattice type
if symmetrize:
C1122_sym = (C1122 + C2211) / 2
C1122 = C2211 = C1122_sym
C1112_sym = (C1112 + C1211) / 2
C1112 = C1211 = C1112_sym
C2212_sym = (C2212 + C1222) / 2
C2212 = C1222 = C2212_sym
Cs_mandel = np.array([[C1111, C1122, np.sqrt(2) * C1112],
[C2211, C2222, np.sqrt(2) * C2212],
[np.sqrt(2) * C1211, np.sqrt(2) * C1222, 2 * C1212]])
Cs_voigt = np.array([[C1111, C1122, C1112],
[C2211, C2222, C2212],
[C1211, C1222, C1212]])
Cs_res = np.array([[C1111_res, C1122_res, C1112_res],
[C2211_res, C2222_res, C2212_res],
[C1211_res, C1222_res, C1212_res]])
if c is not None:
return Cs_voigt * c, Cs_mandel * c, np.sqrt(Cs_res) * c
else:
return Cs_voigt, Cs_mandel, np.sqrt(Cs_res)
def check_elastic_stability(Cs):
"""
Assesses the elastic stability using the criterion that the eigenvalues
of the elastic constants in Mandel notation must be positive.
Ref: Marcin Maździarz 2019 2D Mater.6 048001
(https://iopscience.iop.org/article/10.1088/2053-1583/ab2ef3).
Args:
Cs (Numpy array): 2D elastic constants in the second-rank tensor
or Mandel notation
Returns:
Tuple of (boo, Numpy array): (Stable?, eigenvalues of Cs)
"""
evals = np.linalg.eigvals(Cs)
return np.all(evals > 0), evals
def planar_bulk_modulus(Cs, extensibility=False):
"""
Computes the planar bulk modulus or areal extensibility as defined and derived
in https://doi.org/10.1016/j.carbon.2019.10.041.
Args:
Cs (Numpy array) : 2D elastic constants in Voigt notation
extensibility (bool): Whether to return the areal extensibility instead
Returns:
float: planar bulk modulus (or optionally areal extensibility)
"""
Ss = np.linalg.inv(Cs) # compliance matrix
k = Ss[0, 0] + Ss[1, 1] + 2 * Ss[0, 1]
if extensibility:
return k
else:
return 1 / k