DAMASK_EICMD/python/damask/_orientation.py

236 lines
10 KiB
Python
Raw Normal View History

import numpy as np
2019-05-30 23:32:55 +05:30
2020-03-19 19:49:11 +05:30
from . import Lattice
from . import Rotation
2019-02-21 17:06:27 +05:30
class Orientation: # ToDo: make subclass of lattice and Rotation
2020-03-22 21:33:28 +05:30
"""
Crystallographic orientation.
A crystallographic orientation contains a rotation and a lattice.
"""
__slots__ = ['rotation','lattice']
def __repr__(self):
"""Report lattice type and orientation."""
return self.lattice.__repr__()+'\n'+self.rotation.__repr__()
def __init__(self, rotation, lattice):
"""
New orientation from rotation and lattice.
Parameters
----------
rotation : Rotation
Rotation specifying the lattice orientation.
lattice : Lattice
Lattice type of the crystal.
"""
if isinstance(lattice, Lattice):
self.lattice = lattice
else:
self.lattice = Lattice(lattice) # assume string
if isinstance(rotation, Rotation):
self.rotation = rotation
else:
2020-05-16 14:47:12 +05:30
self.rotation = Rotation.from_quaternion(rotation) # assume quaternion
2020-03-22 21:33:28 +05:30
2020-03-22 21:33:28 +05:30
def disorientation(self,
other,
SST = True,
symmetries = False):
"""
Disorientation between myself and given other orientation.
Rotation axis falls into SST if SST == True.
2020-04-12 19:04:29 +05:30
Currently requires same symmetry for both orientations.
Look into A. Heinz and P. Neumann 1991 for cases with differing sym.
2020-03-22 21:33:28 +05:30
"""
if self.lattice.symmetry != other.lattice.symmetry:
raise NotImplementedError('disorientation between different symmetry classes not supported yet.')
mySymEqs = self.equivalentOrientations() if SST else self.equivalentOrientations([0]) # take all or only first sym operation
otherSymEqs = other.equivalentOrientations()
for i,sA in enumerate(mySymEqs):
aInv = sA.rotation.inversed()
for j,sB in enumerate(otherSymEqs):
b = sB.rotation
r = b*aInv
for k in range(2):
r.inverse()
2020-05-21 03:20:08 +05:30
breaker = self.lattice.symmetry.inFZ(r.as_Rodrigues(vector=True)) \
and (not SST or other.lattice.symmetry.inDisorientationSST(r.as_Rodrigues(vector=True)))
2020-03-22 21:33:28 +05:30
if breaker: break
if breaker: break
if breaker: break
return (Orientation(r,self.lattice), i,j, k == 1) if symmetries else r # disorientation ...
# ... own sym, other sym,
# self-->other: True, self<--other: False
2019-10-23 03:01:27 +05:30
def inFZ_vec(self):
2020-06-28 22:33:06 +05:30
"""Check if orientations fall into Fundamental Zone."""
if not self.rotation.shape:
return self.lattice.symmetry.inFZ(self.rotation.as_Rodrigues(vector=True))
else:
return [self.lattice.symmetry.inFZ(\
2020-06-28 22:33:06 +05:30
self.rotation.as_Rodrigues(vector=True)[l]) for l in range(self.rotation.shape[0])]
def inFZ(self):
return self.lattice.symmetry.inFZ(self.rotation.as_Rodrigues(vector=True))
@property
2020-06-28 22:33:06 +05:30
def equivalent_vec(self):
"""
Return orientations which are symmetrically equivalent.
One dimension (length according to symmetrically equivalent orientations)
is added to the left of the rotation array.
"""
2020-06-28 22:33:06 +05:30
s = self.lattice.symmetry.symmetry_operations #24 lines (sym) x 4 columns (quat)
s = s.reshape(s.shape[:1]+(1,)*len(self.rotation.shape)+(4,)) #reshape zo (24,1,4)
s = Rotation(np.broadcast_to(s,s.shape[:1]+self.rotation.quaternion.shape))
2020-06-28 22:33:06 +05:30
r = np.broadcast_to(self.rotation.quaternion,s.shape[:1]+self.rotation.quaternion.shape) #(24,NumRots,4)
r = Rotation(r) #(24, NumRot)
return self.__class__(s@r,self.lattice)
2020-03-22 21:33:28 +05:30
def equivalentOrientations(self,members=[]):
"""List of orientations which are symmetrically equivalent."""
try:
iter(members) # asking for (even empty) list of members?
except TypeError:
return self.__class__(self.lattice.symmetry.symmetryOperations(members)*self.rotation,self.lattice) # no, return rotation object
else:
return [self.__class__(q*self.rotation,self.lattice) \
for q in self.lattice.symmetry.symmetryOperations(members)] # yes, return list of rotations
def relatedOrientations_vec(self,model):
"""List of orientations related by the given orientation relationship."""
2020-06-28 22:33:06 +05:30
h = self.lattice.relationOperations(model)
rot= h['rotations']
op=np.array([o.as_quaternion() for o in rot])
s = op.reshape(op.shape[:1]+(1,)*len(self.rotation.shape)+(4,))
s = Rotation(np.broadcast_to(s,s.shape[:1]+self.rotation.quaternion.shape))
r = np.broadcast_to(self.rotation.quaternion,s.shape[:1]+self.rotation.quaternion.shape)
r = Rotation(r)
return self.__class__(s@r,h['lattice'])
2020-03-22 21:33:28 +05:30
def relatedOrientations(self,model):
"""List of orientations related by the given orientation relationship."""
r = self.lattice.relationOperations(model)
return [self.__class__(o*self.rotation,r['lattice']) for o in r['rotations']]
2019-10-23 03:01:27 +05:30
2020-06-28 22:33:06 +05:30
@property
def reduced_vec(self):
"""Transform orientation to fall into fundamental zone according to symmetry."""
2020-06-29 21:55:45 +05:30
equi= self.equivalent_vec.rotation #equivalent orientations
r= 1 if not self.rotation.shape else equi.shape[1] #number of rotations
num_equi=equi.shape[0] #number of equivalente orientations
quat= np.reshape( equi.as_quaternion(), (r*num_equi,4) ,order='F') #equivalents are listed in intiuitive order
boolean=Orientation(quat, self.lattice).inFZ_vec() #check which ones are in FZ
if sum(boolean) == r:
return self.__class__(quat[boolean],self.lattice)
else:
print('More than 1 equivalent orientation has been found for an orientation')
index=np.empty(r, dtype=int)
for l,h in enumerate(range(0,r*num_equi, num_equi)):
index[l]=np.where(boolean[h:h+num_equi])[0][0] + (l*num_equi) #get first index that is true then go check to next orientation
return self.__class__(quat[index],self.lattice)
2020-06-28 22:33:06 +05:30
2020-03-22 21:33:28 +05:30
def reduced(self):
"""Transform orientation to fall into fundamental zone according to symmetry."""
for me in self.equivalentOrientations():
2020-05-21 03:20:08 +05:30
if self.lattice.symmetry.inFZ(me.rotation.as_Rodrigues(vector=True)): break
2019-02-21 17:06:27 +05:30
2020-03-22 21:33:28 +05:30
return self.__class__(me.rotation,self.lattice)
2020-03-22 21:33:28 +05:30
def inversePole(self,
axis,
proper = False,
SST = True):
"""Axis rotated according to orientation (using crystal symmetry to ensure location falls into SST)."""
if SST: # pole requested to be within SST
for i,o in enumerate(self.equivalentOrientations()): # test all symmetric equivalent quaternions
pole = o.rotation@axis # align crystal direction to axis
2020-03-22 21:33:28 +05:30
if self.lattice.symmetry.inSST(pole,proper): break # found SST version
else:
pole = self.rotation@axis # align crystal direction to axis
2020-03-22 21:33:28 +05:30
return (pole,i if SST else 0)
2020-03-22 21:33:28 +05:30
def IPFcolor(self,axis):
"""TSL color of inverse pole figure for given axis."""
color = np.zeros(3,'d')
2020-03-22 21:33:28 +05:30
for o in self.equivalentOrientations():
pole = o.rotation@axis # align crystal direction to axis
2020-03-22 21:33:28 +05:30
inSST,color = self.lattice.symmetry.inSST(pole,color=True)
if inSST: break
2020-03-22 21:33:28 +05:30
return color
2020-06-29 21:55:45 +05:30
def IPF_color(self,axis):
2020-06-28 22:35:10 +05:30
"""TSL color of inverse pole figure for given axis. Not for hex or triclinic lattices."""
2020-06-28 22:33:06 +05:30
eq = self.equivalent_vec
pole = eq.rotation @ np.broadcast_to(axis/np.linalg.norm(axis),eq.rotation.shape+(3,))
in_SST, color = self.lattice.symmetry.in_SST(pole,color=True)
# ignore duplicates (occur for highly symmetric orientations)
found = np.zeros_like(in_SST[1],dtype=bool)
c = np.empty(color.shape[1:])
for s in range(in_SST.shape[0]):
c = np.where(np.expand_dims(np.logical_and(in_SST[s],~found),-1),color[s],c)
found = np.logical_or(in_SST[s],found)
return c
2020-03-22 21:33:28 +05:30
@staticmethod
def fromAverage(orientations,
weights = []):
"""Create orientation from average of list of orientations."""
2020-06-23 02:44:58 +05:30
# further read: Orientation distribution analysis in deformed grains, https://doi.org/10.1107/S0021889801003077
2020-03-22 21:33:28 +05:30
if not all(isinstance(item, Orientation) for item in orientations):
raise TypeError("Only instances of Orientation can be averaged.")
2020-03-22 21:33:28 +05:30
closest = []
ref = orientations[0]
for o in orientations:
closest.append(o.equivalentOrientations(
ref.disorientation(o,
SST = False, # select (o[ther]'s) sym orientation
symmetries = True)[2]).rotation) # with lowest misorientation
2020-03-22 21:33:28 +05:30
return Orientation(Rotation.fromAverage(closest,weights),ref.lattice)
2019-04-17 12:59:49 +05:30
2020-03-22 21:33:28 +05:30
def average(self,other):
"""Calculate the average rotation."""
return Orientation.fromAverage([self,other])