Source code for schrodinger.application.matsci.elasticity.strain

"""
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.

Copyright Schrodinger, LLC. All rights reserved.
"""

# Based on pymatgen/analysis/elasticity/strain.py (LEGAL-413)
# last updated from upstream: 8767b89  on Aug 7, 2017

import collections.abc
import itertools

import numpy as np
import scipy

from schrodinger.application.matsci.elasticity.tensors import SquareTensor
from schrodinger.application.matsci.elasticity.tensors import symmetry_reduce
from schrodinger.application.matsci.nano import xtal


[docs]class Deformation(SquareTensor): """ Subclass of SquareTensor that describes the deformation gradient tensor """ 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(Deformation, cls).__new__(cls, deformation_gradient) return obj.view(cls)
[docs] def is_independent(self, tol=1e-8): """ checks to determine whether the deformation is independent """ return len(self.get_perturbed_indices(tol)) == 1
[docs] def get_perturbed_indices(self, tol=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)
[docs] 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() vecs = np.array(xtal.get_vectors_from_chorus(def_struct)) new_vecs = np.array([np.dot(self, vec) for vec in vecs]) xtal.set_pbc_properties(def_struct, new_vecs.flat) frac = xtal.trans_cart_to_frac_from_vecs(def_struct.getXYZ(), *vecs) new_cart = xtal.trans_frac_to_cart_from_vecs(frac, *new_vecs) def_struct.setXYZ(new_cart) return def_struct
[docs] @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)
[docs]class DeformedStructureSet(collections.abc.Sequence): """ class that generates a set of independently deformed structures that can be used to calculate linear stress-strain response """ NORM_INDICES = [(0, 0), (1, 1), (2, 2)] NORM_STRAINS = [-0.01, -0.005, 0.005, 0.01] SHEAR_INDICES = [(0, 1), (0, 2), (1, 2)] SHEAR_STRAINS = [-0.06, -0.03, 0.03, 0.06]
[docs] def __init__(self, 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. :param structure: structure to undergo deformation :type structure: Structure :param: norm_strains: strain values to apply to each normal mode. :type norm_strains: list[float] :param shear_strains: strain values to apply to each shear mode. :type shear_strains: list[float] :param symmetry: whether or not to use symmetry reduction. :type symmetry: bool """ norm_strains = norm_strains or self.NORM_STRAINS shear_strains = shear_strains or self.SHEAR_STRAINS self.undeformed_structure = structure self.strains = [] self.deformations = [] self.def_structs = [] self.sym_dict = {} # Generate deformations for ind in self.NORM_INDICES: for amount in norm_strains: strain = Strain.from_index_amount(ind, amount) self.strains.append(strain) self.deformations.append(strain.deformation_matrix) for ind in self.SHEAR_INDICES: for amount in shear_strains: strain = Strain.from_index_amount(ind, amount) self.strains.append(strain) self.deformations.append(strain.deformation_matrix) # Perform symmetry reduction if specified if symmetry: self.sym_dict = symmetry_reduce(self.deformations, structure) self.deformations = list(self.sym_dict.keys()) self.deformed_structures = [ defo.apply_to_structure(structure) for defo in self.deformations ]
def __iter__(self): return iter(self.deformed_structures)
[docs] def __len__(self): return len(self.deformed_structures)
def __getitem__(self, ind): return self.deformed_structures[ind]
[docs]class Strain(SquareTensor): """ Subclass of SquareTensor that describes the Green-Lagrange strain tensor. """ def __new__(cls, strain_matrix, dfm=None, dfm_shape="upper"): """ 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(Strain, cls).__new__(cls, strain_matrix, vscale=vscale) if dfm is None: obj._dfm = convert_strain_to_deformation(obj, dfm_shape) else: dfm = Deformation(dfm) gls_test = 0.5 * (np.dot(dfm.trans, dfm) - np.eye(3)) if (gls_test - obj > 1e-10).any(): raise ValueError("Strain and deformation gradients " "do not match!") obj._dfm = Deformation(dfm) 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._dfm = getattr(obj, "_dfm", None) self._vscale = getattr(obj, "_vscale", None)
[docs] @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)), dfm)
[docs] @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. :param idx: index to be perturbed, can be voigt or full-tensor notation :type idx: tuple or integer :param amount: amount to perturb selected index :type amount: float """ if np.array(idx).ndim == 0: v = np.zeros(6) v[idx] = amount return cls.from_voigt(v) elif np.array(idx).ndim == 1: v = np.zeros((3, 3)) for i in itertools.permutations(idx): v[i] = amount return cls(v) else: raise ValueError("Index must either be 2-tuple or integer " "corresponding to full-tensor or voigt index")
@property def deformation_matrix(self): """ returns the deformation matrix """ return self._dfm @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)
[docs]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: :param strain: (3x3 array-like) strain matrix :param shape: method for determining deformation, supports: - "upper" produces an upper triangular defo - "lower" produces a lower triangular defo - "symmetric" produces a symmetric defo :type shape: str """ 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)