-
Notifications
You must be signed in to change notification settings - Fork 854
/
strain.py
263 lines (219 loc) · 8.63 KB
/
strain.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
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
# Copyright (c) Pymatgen Development Team.
# Distributed under the terms of the MIT License.
"""
This module provides classes and methods used to describe deformations and
strains, including applying those deformations to structure objects and
generating deformed structure sets for further calculations.
"""
import collections
import itertools
import numpy as np
import scipy
from pymatgen.core.lattice import Lattice
from pymatgen.core.structure import Structure
from pymatgen.core.tensors import SquareTensor, symmetry_reduce
__author__ = "Joseph Montoya"
__copyright__ = "Copyright 2012, The Materials Project"
__credits__ = "Maarten de Jong, Mark Asta, Anubhav Jain"
__version__ = "1.0"
__maintainer__ = "Joseph Montoya"
__email__ = "montoyjh@lbl.gov"
__status__ = "Production"
__date__ = "July 24, 2018"
class Deformation(SquareTensor):
"""
Subclass of SquareTensor that describes the deformation gradient tensor
"""
symbol = "d"
def __new__(cls, deformation_gradient):
"""
Create a Deformation object. Note that the constructor uses __new__
rather than __init__ according to the standard method of subclassing
numpy ndarrays.
Args:
deformation_gradient (3x3 array-like): the 3x3 array-like
representing the deformation gradient
"""
obj = super().__new__(cls, deformation_gradient)
return obj.view(cls)
def is_independent(self, tol: float = 1e-8):
"""
checks to determine whether the deformation is independent
"""
return len(self.get_perturbed_indices(tol)) == 1
def get_perturbed_indices(self, tol: float = 1e-8):
"""
Gets indices of perturbed elements of the deformation gradient,
i. e. those that differ from the identity
"""
indices = list(zip(*np.where(abs(self - np.eye(3)) > tol)))
return indices
@property
def green_lagrange_strain(self):
"""
calculates the euler-lagrange strain from
the deformation gradient
"""
return Strain.from_deformation(self)
def apply_to_structure(self, structure):
"""
Apply the deformation gradient to a structure.
Args:
structure (Structure object): the structure object to
be modified by the deformation
"""
def_struct = structure.copy()
old_latt = def_struct.lattice.matrix
new_latt = np.transpose(np.dot(self, np.transpose(old_latt)))
def_struct.lattice = Lattice(new_latt)
return def_struct
@classmethod
def from_index_amount(cls, matrixpos, amt):
"""
Factory method for constructing a Deformation object
from a matrix position and amount
Args:
matrixpos (tuple): tuple corresponding the matrix position to
have a perturbation added
amt (float): amount to add to the identity matrix at position
matrixpos
"""
f = np.identity(3)
f[matrixpos] += amt
return cls(f)
class DeformedStructureSet(collections.abc.Sequence):
"""
class that generates a set of independently deformed structures that
can be used to calculate linear stress-strain response
"""
def __init__(self, structure: Structure, norm_strains=None, shear_strains=None, symmetry=False):
"""
constructs the deformed geometries of a structure. Generates
m + n deformed structures according to the supplied parameters.
Args:
structure (Structure): structure to undergo deformation
norm_strains (list of floats): strain values to apply
to each normal mode.
shear_strains (list of floats): strain values to apply
to each shear mode.
symmetry (bool): whether or not to use symmetry reduction.
"""
norm_strains = norm_strains or [-0.01, -0.005, 0.005, 0.01]
shear_strains = shear_strains or [-0.06, -0.03, 0.03, 0.06]
self.undeformed_structure = structure
self.deformations: list[Deformation] = []
self.def_structs: list[Structure] = []
# Generate deformations
for ind in [(0, 0), (1, 1), (2, 2)]:
for amount in norm_strains:
strain = Strain.from_index_amount(ind, amount)
self.deformations.append(strain.get_deformation_matrix())
for ind in [(0, 1), (0, 2), (1, 2)]:
for amount in shear_strains:
strain = Strain.from_index_amount(ind, amount)
self.deformations.append(strain.get_deformation_matrix())
# Perform symmetry reduction if specified
if symmetry:
self.sym_dict = symmetry_reduce(self.deformations, structure)
self.deformations = list(self.sym_dict)
self.deformed_structures = [defo.apply_to_structure(structure) for defo in self.deformations]
def __iter__(self):
return iter(self.deformed_structures)
def __len__(self):
return len(self.deformed_structures)
def __getitem__(self, ind):
return self.deformed_structures[ind]
class Strain(SquareTensor):
"""
Subclass of SquareTensor that describes the Green-Lagrange strain tensor.
"""
symbol = "e"
def __new__(cls, strain_matrix):
"""
Create a Strain object. Note that the constructor uses __new__
rather than __init__ according to the standard method of
subclassing numpy ndarrays. Note also that the default constructor
does not include the deformation gradient
Args:
strain_matrix (3x3 array-like): the 3x3 array-like
representing the Green-Lagrange strain
"""
vscale = np.ones((6,))
vscale[3:] *= 2
obj = super().__new__(cls, strain_matrix, vscale=vscale)
if not obj.is_symmetric():
raise ValueError(
"Strain objects must be initialized "
"with a symmetric array or a voigt-notation "
"vector with six entries."
)
return obj.view(cls)
def __array_finalize__(self, obj):
if obj is None:
return
self.rank = getattr(obj, "rank", None)
self._vscale = getattr(obj, "_vscale", None)
@classmethod
def from_deformation(cls, deformation):
"""
Factory method that returns a Strain object from a deformation
gradient
Args:
deformation (3x3 array-like):
"""
dfm = Deformation(deformation)
return cls(0.5 * (np.dot(dfm.trans, dfm) - np.eye(3)))
@classmethod
def from_index_amount(cls, idx, amount):
"""
Like Deformation.from_index_amount, except generates
a strain from the zero 3x3 tensor or voigt vector with
the amount specified in the index location. Ensures
symmetric strain.
Args:
idx (tuple or integer): index to be perturbed, can be voigt or
full-tensor notation
amount (float): amount to perturb selected index
"""
if np.array(idx).ndim == 0:
v = np.zeros(6)
v[idx] = amount
return cls.from_voigt(v)
if np.array(idx).ndim == 1:
v = np.zeros((3, 3))
for i in itertools.permutations(idx):
v[i] = amount
return cls(v)
raise ValueError("Index must either be 2-tuple or integer corresponding to full-tensor or voigt index")
def get_deformation_matrix(self, shape="upper"):
"""
returns the deformation matrix
"""
return convert_strain_to_deformation(self, shape=shape)
@property
def von_mises_strain(self):
"""
Equivalent strain to Von Mises Stress
"""
eps = self - 1 / 3 * np.trace(self) * np.identity(3)
return np.sqrt(np.sum(eps * eps) * 2 / 3)
def convert_strain_to_deformation(strain, shape="upper"):
"""
This function converts a strain to a deformation gradient that will
produce that strain. Supports three methods:
Args:
strain (3x3 array-like): strain matrix
shape: (string): method for determining deformation, supports
"upper" produces an upper triangular defo
"lower" produces a lower triangular defo
"symmetric" produces a symmetric defo
"""
strain = SquareTensor(strain)
ftdotf = 2 * strain + np.eye(3)
if shape == "upper":
result = scipy.linalg.cholesky(ftdotf)
elif shape == "symmetric":
result = scipy.linalg.sqrtm(ftdotf)
else:
raise ValueError('shape must be "upper" or "symmetric"')
return Deformation(result)