2021-06-03 13:48:27 +05:30
|
|
|
|
import copy
|
2023-02-24 00:19:08 +05:30
|
|
|
|
from typing import Optional, Union, TypeVar
|
2021-02-11 01:38:48 +05:30
|
|
|
|
|
2015-05-08 19:44:44 +05:30
|
|
|
|
import numpy as np
|
2019-05-30 23:32:55 +05:30
|
|
|
|
|
2023-08-17 01:31:20 +05:30
|
|
|
|
from ._typehints import FloatSequence, IntSequence, CrystalFamily, BravaisLattice
|
2020-03-19 19:49:11 +05:30
|
|
|
|
from . import Rotation
|
2021-06-08 01:19:04 +05:30
|
|
|
|
from . import Crystal
|
2020-11-10 01:50:56 +05:30
|
|
|
|
from . import util
|
2020-11-16 03:44:46 +05:30
|
|
|
|
from . import tensor
|
2021-04-29 23:14:21 +05:30
|
|
|
|
|
2022-02-13 05:54:02 +05:30
|
|
|
|
MyType = TypeVar('MyType', bound='Orientation')
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
2021-06-08 01:19:04 +05:30
|
|
|
|
class Orientation(Rotation,Crystal):
|
2020-03-22 21:33:28 +05:30
|
|
|
|
"""
|
2020-11-10 01:50:56 +05:30
|
|
|
|
Representation of crystallographic orientation as combination of rotation and either crystal family or Bravais lattice.
|
|
|
|
|
|
2021-04-29 23:14:21 +05:30
|
|
|
|
The crystal family is one of:
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
- triclinic
|
|
|
|
|
- monoclinic
|
|
|
|
|
- orthorhombic
|
|
|
|
|
- tetragonal
|
|
|
|
|
- hexagonal
|
|
|
|
|
- cubic
|
2020-03-22 21:33:28 +05:30
|
|
|
|
|
2020-11-10 01:50:56 +05:30
|
|
|
|
and enables symmetry-related operations such as
|
|
|
|
|
"equivalent", "reduced", "disorientation", "IPF_color", or "to_SST".
|
|
|
|
|
|
2021-04-29 23:14:21 +05:30
|
|
|
|
The Bravais lattice is given in the Pearson notation:
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
2021-03-27 14:40:35 +05:30
|
|
|
|
- triclinic
|
|
|
|
|
- aP : primitive
|
|
|
|
|
|
2021-07-05 02:55:00 +05:30
|
|
|
|
- monoclinic
|
2021-03-27 14:40:35 +05:30
|
|
|
|
- mP : primitive
|
|
|
|
|
- mS : base-centered
|
|
|
|
|
|
|
|
|
|
- orthorhombic
|
|
|
|
|
- oP : primitive
|
|
|
|
|
- oS : base-centered
|
|
|
|
|
- oI : body-centered
|
|
|
|
|
- oF : face-centered
|
|
|
|
|
|
|
|
|
|
- tetragonal
|
|
|
|
|
- tP : primitive
|
|
|
|
|
- tI : body-centered
|
|
|
|
|
|
|
|
|
|
- hexagonal
|
|
|
|
|
- hP : primitive
|
|
|
|
|
|
|
|
|
|
- cubic
|
|
|
|
|
- cP : primitive
|
|
|
|
|
- cI : body-centered
|
|
|
|
|
- cF : face-centered
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
and inherits the corresponding crystal family.
|
|
|
|
|
Specifying a Bravais lattice, compared to just the crystal family,
|
|
|
|
|
extends the functionality of Orientation objects to include operations such as
|
|
|
|
|
"Schmid", "related", or "to_pole" that require a lattice type and its parameters.
|
|
|
|
|
|
|
|
|
|
Examples
|
|
|
|
|
--------
|
|
|
|
|
An array of 3 x 5 random orientations reduced to the fundamental zone of tetragonal symmetry:
|
|
|
|
|
|
2021-07-25 23:01:48 +05:30
|
|
|
|
>>> import damask
|
2021-08-09 03:26:54 +05:30
|
|
|
|
>>> o=damask.Orientation.from_random(shape=(3,5),family='tetragonal').reduced
|
2020-03-22 21:33:28 +05:30
|
|
|
|
|
|
|
|
|
"""
|
|
|
|
|
|
2023-03-03 00:16:00 +05:30
|
|
|
|
@util.extend_docstring(adopted_parameters=Crystal.__init__)
|
2020-11-10 01:50:56 +05:30
|
|
|
|
def __init__(self,
|
2022-01-28 14:41:34 +05:30
|
|
|
|
rotation: Union[FloatSequence, Rotation] = np.array([1.,0.,0.,0.]),
|
|
|
|
|
*,
|
2022-11-23 02:56:15 +05:30
|
|
|
|
family: Optional[CrystalFamily] = None,
|
2023-08-17 01:31:20 +05:30
|
|
|
|
lattice: Optional[BravaisLattice] = None,
|
2022-11-23 02:56:15 +05:30
|
|
|
|
a: Optional[float] = None, b: Optional[float] = None, c: Optional[float] = None,
|
|
|
|
|
alpha: Optional[float] = None, beta: Optional[float] = None, gamma: Optional[float] = None,
|
2021-12-14 21:35:00 +05:30
|
|
|
|
degrees: bool = False):
|
2020-03-22 21:33:28 +05:30
|
|
|
|
"""
|
2021-03-27 14:40:35 +05:30
|
|
|
|
New orientation.
|
2020-03-22 21:33:28 +05:30
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2023-02-21 20:57:06 +05:30
|
|
|
|
rotation : list, numpy.ndarray, or Rotation, optional
|
2020-11-10 01:50:56 +05:30
|
|
|
|
Unit quaternion in positive real hemisphere.
|
|
|
|
|
Use .from_quaternion to perform a sanity check.
|
|
|
|
|
Defaults to no rotation.
|
2020-03-22 21:33:28 +05:30
|
|
|
|
|
|
|
|
|
"""
|
2021-06-06 23:19:29 +05:30
|
|
|
|
Rotation.__init__(self,rotation)
|
2021-06-08 01:19:04 +05:30
|
|
|
|
Crystal.__init__(self,family=family, lattice=lattice,
|
2021-06-06 23:19:29 +05:30
|
|
|
|
a=a,b=b,c=c, alpha=alpha,beta=beta,gamma=gamma, degrees=degrees)
|
2020-03-22 21:33:28 +05:30
|
|
|
|
|
|
|
|
|
|
2021-12-14 21:35:00 +05:30
|
|
|
|
def __repr__(self) -> str:
|
2022-07-08 21:36:41 +05:30
|
|
|
|
"""
|
|
|
|
|
Return repr(self).
|
|
|
|
|
|
2023-02-21 20:57:06 +05:30
|
|
|
|
Give short, human-readable summary.
|
2022-07-08 21:36:41 +05:30
|
|
|
|
|
|
|
|
|
"""
|
2022-02-23 19:07:23 +05:30
|
|
|
|
return util.srepr([Crystal.__repr__(self),
|
|
|
|
|
Rotation.__repr__(self)])
|
|
|
|
|
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
2022-02-13 05:54:02 +05:30
|
|
|
|
def __copy__(self: MyType,
|
2022-11-23 02:56:15 +05:30
|
|
|
|
rotation: Union[None, FloatSequence, Rotation] = None) -> MyType:
|
2022-07-08 21:36:41 +05:30
|
|
|
|
"""
|
|
|
|
|
Return deepcopy(self).
|
|
|
|
|
|
|
|
|
|
Create deep copy.
|
|
|
|
|
|
|
|
|
|
"""
|
2021-06-03 13:48:27 +05:30
|
|
|
|
dup = copy.deepcopy(self)
|
|
|
|
|
if rotation is not None:
|
2021-11-01 03:07:54 +05:30
|
|
|
|
dup.quaternion = Rotation(rotation).quaternion
|
2021-06-03 13:48:27 +05:30
|
|
|
|
return dup
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
copy = __copy__
|
|
|
|
|
|
2020-06-30 17:25:09 +05:30
|
|
|
|
|
2021-12-14 21:35:00 +05:30
|
|
|
|
|
2022-01-26 19:39:09 +05:30
|
|
|
|
def __eq__(self,
|
|
|
|
|
other: object) -> bool:
|
2020-11-10 01:50:56 +05:30
|
|
|
|
"""
|
2022-07-08 21:36:41 +05:30
|
|
|
|
Return self==other.
|
|
|
|
|
|
2022-07-08 21:37:07 +05:30
|
|
|
|
Test equality of other.
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
|
|
|
|
other : Orientation
|
|
|
|
|
Orientation to check for equality.
|
2020-04-24 10:22:09 +05:30
|
|
|
|
|
2020-03-22 21:33:28 +05:30
|
|
|
|
"""
|
2021-12-14 21:35:00 +05:30
|
|
|
|
if not isinstance(other, Orientation):
|
2022-02-04 21:27:25 +05:30
|
|
|
|
return NotImplemented
|
2021-06-06 23:19:29 +05:30
|
|
|
|
matching_type = self.family == other.family and \
|
|
|
|
|
self.lattice == other.lattice and \
|
|
|
|
|
self.parameters == other.parameters
|
2021-04-10 11:59:42 +05:30
|
|
|
|
return np.logical_and(matching_type,super(self.__class__,self.reduced).__eq__(other.reduced))
|
2021-01-04 02:19:01 +05:30
|
|
|
|
|
2022-01-26 19:39:09 +05:30
|
|
|
|
def __ne__(self,
|
|
|
|
|
other: object) -> bool:
|
2021-01-04 02:19:01 +05:30
|
|
|
|
"""
|
2022-07-08 21:36:41 +05:30
|
|
|
|
Return self!=other.
|
|
|
|
|
|
2022-07-08 21:37:07 +05:30
|
|
|
|
Test inequality of other.
|
2021-01-04 02:19:01 +05:30
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
|
|
|
|
other : Orientation
|
|
|
|
|
Orientation to check for equality.
|
|
|
|
|
|
|
|
|
|
"""
|
2022-02-08 19:17:23 +05:30
|
|
|
|
return np.logical_not(self==other) if isinstance(other, Orientation) else NotImplemented
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
|
2022-02-16 02:38:12 +05:30
|
|
|
|
def isclose(self: MyType,
|
|
|
|
|
other: MyType,
|
2021-12-14 21:35:00 +05:30
|
|
|
|
rtol: float = 1e-5,
|
|
|
|
|
atol: float = 1e-8,
|
|
|
|
|
equal_nan: bool = True) -> bool:
|
2021-04-05 21:54:03 +05:30
|
|
|
|
"""
|
|
|
|
|
Report where values are approximately equal to corresponding ones of other Orientation.
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
|
|
|
|
other : Orientation
|
|
|
|
|
Orientation to compare against.
|
|
|
|
|
rtol : float, optional
|
|
|
|
|
Relative tolerance of equality.
|
|
|
|
|
atol : float, optional
|
|
|
|
|
Absolute tolerance of equality.
|
|
|
|
|
equal_nan : bool, optional
|
|
|
|
|
Consider matching NaN values as equal. Defaults to True.
|
|
|
|
|
|
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-02-13 05:54:02 +05:30
|
|
|
|
mask : numpy.ndarray of bool, shape (self.shape)
|
2021-04-05 21:54:03 +05:30
|
|
|
|
Mask indicating where corresponding orientations are close.
|
|
|
|
|
|
|
|
|
|
"""
|
2021-06-06 23:19:29 +05:30
|
|
|
|
matching_type = self.family == other.family and \
|
|
|
|
|
self.lattice == other.lattice and \
|
|
|
|
|
self.parameters == other.parameters
|
2021-04-10 11:59:42 +05:30
|
|
|
|
return np.logical_and(matching_type,super(self.__class__,self.reduced).isclose(other.reduced))
|
2021-04-05 21:54:03 +05:30
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2022-02-16 02:38:12 +05:30
|
|
|
|
def allclose(self: MyType,
|
|
|
|
|
other: MyType,
|
2021-12-14 21:35:00 +05:30
|
|
|
|
rtol: float = 1e-5,
|
|
|
|
|
atol: float = 1e-8,
|
2022-02-06 21:41:18 +05:30
|
|
|
|
equal_nan: bool = True) -> bool:
|
2021-04-05 21:54:03 +05:30
|
|
|
|
"""
|
|
|
|
|
Test whether all values are approximately equal to corresponding ones of other Orientation.
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
|
|
|
|
other : Orientation
|
|
|
|
|
Orientation to compare against.
|
|
|
|
|
rtol : float, optional
|
|
|
|
|
Relative tolerance of equality.
|
|
|
|
|
atol : float, optional
|
|
|
|
|
Absolute tolerance of equality.
|
|
|
|
|
equal_nan : bool, optional
|
|
|
|
|
Consider matching NaN values as equal. Defaults to True.
|
|
|
|
|
|
|
|
|
|
Returns
|
|
|
|
|
-------
|
|
|
|
|
answer : bool
|
|
|
|
|
Whether all values are close between both orientations.
|
|
|
|
|
|
|
|
|
|
"""
|
2022-02-06 21:41:18 +05:30
|
|
|
|
return bool(np.all(self.isclose(other,rtol,atol,equal_nan)))
|
2021-04-05 21:54:03 +05:30
|
|
|
|
|
|
|
|
|
|
2022-02-13 05:54:02 +05:30
|
|
|
|
def __mul__(self: MyType,
|
|
|
|
|
other: Union[Rotation, 'Orientation']) -> MyType:
|
2020-11-10 01:50:56 +05:30
|
|
|
|
"""
|
2022-08-11 00:21:50 +05:30
|
|
|
|
Return self*other.
|
|
|
|
|
|
|
|
|
|
Compose with other.
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-08-11 00:21:50 +05:30
|
|
|
|
other : Rotation or Orientation, shape (self.shape)
|
2021-01-13 05:26:40 +05:30
|
|
|
|
Object for composition.
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
Returns
|
|
|
|
|
-------
|
2021-01-13 05:26:40 +05:30
|
|
|
|
composition : Orientation
|
|
|
|
|
Compound rotation self*other, i.e. first other then self rotation.
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
"""
|
2023-02-24 00:19:08 +05:30
|
|
|
|
if not isinstance(other, (Orientation,Rotation)):
|
2021-05-27 23:42:20 +05:30
|
|
|
|
raise TypeError('use "O@b", i.e. matmul, to apply Orientation "O" to object "b"')
|
2023-02-24 00:19:08 +05:30
|
|
|
|
return self.copy(Rotation(self.quaternion)*Rotation(other.quaternion))
|
2021-02-12 21:34:35 +05:30
|
|
|
|
|
2021-02-11 01:38:48 +05:30
|
|
|
|
|
2020-11-10 01:50:56 +05:30
|
|
|
|
@classmethod
|
2022-11-19 13:40:00 +05:30
|
|
|
|
@util.extend_docstring(Rotation.from_random,
|
2023-03-03 00:16:00 +05:30
|
|
|
|
adopted_parameters=Crystal.__init__)
|
2023-02-24 00:19:08 +05:30
|
|
|
|
@util.pass_on('rotation', Rotation.from_random, wrapped=__init__)
|
2022-02-02 15:41:59 +05:30
|
|
|
|
def from_random(cls, **kwargs) -> 'Orientation':
|
2023-02-24 00:19:08 +05:30
|
|
|
|
return cls(**kwargs)
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
@classmethod
|
2022-11-19 13:40:00 +05:30
|
|
|
|
@util.extend_docstring(Rotation.from_quaternion,
|
2023-03-03 00:16:00 +05:30
|
|
|
|
adopted_parameters=Crystal.__init__)
|
2023-02-24 00:19:08 +05:30
|
|
|
|
@util.pass_on('rotation', Rotation.from_quaternion, wrapped=__init__)
|
2022-02-02 15:41:59 +05:30
|
|
|
|
def from_quaternion(cls, **kwargs) -> 'Orientation':
|
2023-02-24 00:19:08 +05:30
|
|
|
|
return cls(**kwargs)
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
@classmethod
|
2022-11-19 13:40:00 +05:30
|
|
|
|
@util.extend_docstring(Rotation.from_Euler_angles,
|
2023-03-03 00:16:00 +05:30
|
|
|
|
adopted_parameters=Crystal.__init__)
|
2023-02-24 00:19:08 +05:30
|
|
|
|
@util.pass_on('rotation', Rotation.from_Euler_angles, wrapped=__init__)
|
2022-02-02 15:41:59 +05:30
|
|
|
|
def from_Euler_angles(cls, **kwargs) -> 'Orientation':
|
2023-02-24 00:19:08 +05:30
|
|
|
|
return cls(**kwargs)
|
2021-02-12 21:34:35 +05:30
|
|
|
|
|
2020-11-10 01:50:56 +05:30
|
|
|
|
@classmethod
|
2022-11-19 13:40:00 +05:30
|
|
|
|
@util.extend_docstring(Rotation.from_axis_angle,
|
2023-03-03 00:16:00 +05:30
|
|
|
|
adopted_parameters=Crystal.__init__)
|
2023-02-24 00:19:08 +05:30
|
|
|
|
@util.pass_on('rotation', Rotation.from_axis_angle, wrapped=__init__)
|
2022-02-02 15:41:59 +05:30
|
|
|
|
def from_axis_angle(cls, **kwargs) -> 'Orientation':
|
2023-02-24 00:19:08 +05:30
|
|
|
|
return cls(**kwargs)
|
2021-02-12 21:34:35 +05:30
|
|
|
|
|
2020-11-10 01:50:56 +05:30
|
|
|
|
@classmethod
|
2022-11-19 13:40:00 +05:30
|
|
|
|
@util.extend_docstring(Rotation.from_basis,
|
2023-03-03 00:16:00 +05:30
|
|
|
|
adopted_parameters=Crystal.__init__)
|
2023-02-24 00:19:08 +05:30
|
|
|
|
@util.pass_on('rotation', Rotation.from_basis, wrapped=__init__)
|
2022-02-02 15:41:59 +05:30
|
|
|
|
def from_basis(cls, **kwargs) -> 'Orientation':
|
2023-02-24 00:19:08 +05:30
|
|
|
|
return cls(**kwargs)
|
2021-02-12 21:34:35 +05:30
|
|
|
|
|
2020-11-10 01:50:56 +05:30
|
|
|
|
@classmethod
|
2022-11-19 13:40:00 +05:30
|
|
|
|
@util.extend_docstring(Rotation.from_matrix,
|
2023-03-03 00:16:00 +05:30
|
|
|
|
adopted_parameters=Crystal.__init__)
|
2023-02-24 00:19:08 +05:30
|
|
|
|
@util.pass_on('rotation', Rotation.from_matrix, wrapped=__init__)
|
2022-02-02 15:41:59 +05:30
|
|
|
|
def from_matrix(cls, **kwargs) -> 'Orientation':
|
2023-02-24 00:19:08 +05:30
|
|
|
|
return cls(**kwargs)
|
2021-02-12 21:34:35 +05:30
|
|
|
|
|
2020-11-10 01:50:56 +05:30
|
|
|
|
@classmethod
|
2022-11-19 13:40:00 +05:30
|
|
|
|
@util.extend_docstring(Rotation.from_Rodrigues_vector,
|
2023-03-03 00:16:00 +05:30
|
|
|
|
adopted_parameters=Crystal.__init__)
|
2023-02-24 00:19:08 +05:30
|
|
|
|
@util.pass_on('rotation', Rotation.from_Rodrigues_vector, wrapped=__init__)
|
2022-02-02 15:41:59 +05:30
|
|
|
|
def from_Rodrigues_vector(cls, **kwargs) -> 'Orientation':
|
2023-02-24 00:19:08 +05:30
|
|
|
|
return cls(**kwargs)
|
2021-02-12 21:34:35 +05:30
|
|
|
|
|
2020-11-10 01:50:56 +05:30
|
|
|
|
@classmethod
|
2022-11-19 13:40:00 +05:30
|
|
|
|
@util.extend_docstring(Rotation.from_homochoric,
|
2023-03-03 00:16:00 +05:30
|
|
|
|
adopted_parameters=Crystal.__init__)
|
2023-02-24 00:19:08 +05:30
|
|
|
|
@util.pass_on('rotation', Rotation.from_homochoric, wrapped=__init__)
|
2022-02-02 15:41:59 +05:30
|
|
|
|
def from_homochoric(cls, **kwargs) -> 'Orientation':
|
2023-02-24 00:19:08 +05:30
|
|
|
|
return cls(**kwargs)
|
2021-02-12 21:34:35 +05:30
|
|
|
|
|
2020-11-10 01:50:56 +05:30
|
|
|
|
@classmethod
|
2022-11-19 13:40:00 +05:30
|
|
|
|
@util.extend_docstring(Rotation.from_cubochoric,
|
2023-03-03 00:16:00 +05:30
|
|
|
|
adopted_parameters=Crystal.__init__)
|
2023-02-24 00:19:08 +05:30
|
|
|
|
@util.pass_on('rotation', Rotation.from_cubochoric, wrapped=__init__)
|
2022-02-02 15:41:59 +05:30
|
|
|
|
def from_cubochoric(cls, **kwargs) -> 'Orientation':
|
2023-02-24 00:19:08 +05:30
|
|
|
|
return cls(**kwargs)
|
2021-02-12 21:34:35 +05:30
|
|
|
|
|
2020-11-10 01:50:56 +05:30
|
|
|
|
@classmethod
|
2022-11-19 13:40:00 +05:30
|
|
|
|
@util.extend_docstring(Rotation.from_spherical_component,
|
2023-03-03 00:16:00 +05:30
|
|
|
|
adopted_parameters=Crystal.__init__)
|
2023-02-24 00:19:08 +05:30
|
|
|
|
@util.pass_on('rotation', Rotation.from_spherical_component, wrapped=__init__)
|
2022-02-02 15:41:59 +05:30
|
|
|
|
def from_spherical_component(cls, **kwargs) -> 'Orientation':
|
2023-02-24 00:19:08 +05:30
|
|
|
|
return cls(**kwargs)
|
2021-02-12 21:34:35 +05:30
|
|
|
|
|
2020-11-10 01:50:56 +05:30
|
|
|
|
@classmethod
|
2022-11-19 13:40:00 +05:30
|
|
|
|
@util.extend_docstring(Rotation.from_fiber_component,
|
2023-03-03 00:16:00 +05:30
|
|
|
|
adopted_parameters=Crystal.__init__)
|
2023-02-24 00:19:08 +05:30
|
|
|
|
@util.pass_on('rotation', Rotation.from_fiber_component, wrapped=__init__)
|
2022-02-02 15:41:59 +05:30
|
|
|
|
def from_fiber_component(cls, **kwargs) -> 'Orientation':
|
2023-02-24 00:19:08 +05:30
|
|
|
|
return cls(**kwargs)
|
2021-02-12 21:34:35 +05:30
|
|
|
|
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
@classmethod
|
2023-03-03 00:16:00 +05:30
|
|
|
|
@util.extend_docstring(adopted_parameters=Crystal.__init__)
|
2021-12-14 21:35:00 +05:30
|
|
|
|
def from_directions(cls,
|
2022-01-14 19:07:48 +05:30
|
|
|
|
uvw: FloatSequence,
|
|
|
|
|
hkl: FloatSequence,
|
2022-02-02 15:41:59 +05:30
|
|
|
|
**kwargs) -> 'Orientation':
|
2020-03-22 21:33:28 +05:30
|
|
|
|
"""
|
2020-11-10 01:50:56 +05:30
|
|
|
|
Initialize orientation object from two crystallographic directions.
|
2020-03-22 21:33:28 +05:30
|
|
|
|
|
2020-11-10 01:50:56 +05:30
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-02-11 22:40:29 +05:30
|
|
|
|
uvw : numpy.ndarray, shape (...,3)
|
|
|
|
|
Lattice direction aligned with lab frame x-direction.
|
|
|
|
|
hkl : numpy.ndarray, shape (...,3)
|
|
|
|
|
Lattice plane normal aligned with lab frame z-direction.
|
2020-03-22 21:33:28 +05:30
|
|
|
|
|
2022-11-19 13:40:00 +05:30
|
|
|
|
Returns
|
|
|
|
|
-------
|
|
|
|
|
new : damask.Orientation
|
|
|
|
|
|
2020-11-10 01:50:56 +05:30
|
|
|
|
"""
|
|
|
|
|
o = cls(**kwargs)
|
2021-06-06 23:19:29 +05:30
|
|
|
|
x = o.to_frame(uvw=uvw)
|
|
|
|
|
z = o.to_frame(hkl=hkl)
|
2020-11-10 01:50:56 +05:30
|
|
|
|
om = np.stack([x,np.cross(z,x),z],axis=-2)
|
2022-02-13 05:54:02 +05:30
|
|
|
|
return o.copy(Rotation.from_matrix(tensor.transpose(om/np.linalg.norm(om,axis=-1,keepdims=True))))
|
2020-03-22 21:33:28 +05:30
|
|
|
|
|
2019-10-23 03:01:27 +05:30
|
|
|
|
|
2020-06-19 02:23:04 +05:30
|
|
|
|
@property
|
2022-02-13 05:54:02 +05:30
|
|
|
|
def equivalent(self: MyType) -> MyType:
|
2020-06-19 02:23:04 +05:30
|
|
|
|
"""
|
2020-11-10 01:50:56 +05:30
|
|
|
|
Orientations that are symmetrically equivalent.
|
2020-06-19 02:23:04 +05:30
|
|
|
|
|
2020-11-10 01:50:56 +05:30
|
|
|
|
One dimension (length corresponds to number of symmetrically equivalent orientations)
|
2020-07-01 04:07:02 +05:30
|
|
|
|
is added to the left of the Rotation array.
|
2020-06-19 02:23:04 +05:30
|
|
|
|
|
|
|
|
|
"""
|
2021-06-06 23:19:29 +05:30
|
|
|
|
sym_ops = self.symmetry_operations
|
2021-06-03 13:32:49 +05:30
|
|
|
|
o = sym_ops.broadcast_to(sym_ops.shape+self.shape,mode='right')
|
2022-02-13 05:54:02 +05:30
|
|
|
|
return self.copy(o*Rotation(self.quaternion).broadcast_to(o.shape,mode='left'))
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
|
|
|
|
|
@property
|
2022-02-13 05:54:02 +05:30
|
|
|
|
def reduced(self: MyType) -> MyType:
|
2020-11-10 01:50:56 +05:30
|
|
|
|
"""Select symmetrically equivalent orientation that falls into fundamental zone according to symmetry."""
|
|
|
|
|
eq = self.equivalent
|
|
|
|
|
ok = eq.in_FZ
|
|
|
|
|
ok &= np.cumsum(ok,axis=0) == 1
|
|
|
|
|
loc = np.where(ok)
|
|
|
|
|
sort = 0 if len(loc) == 1 else np.lexsort(loc[:0:-1])
|
|
|
|
|
return eq[ok][sort].reshape(self.shape)
|
|
|
|
|
|
2022-02-06 21:41:18 +05:30
|
|
|
|
|
2020-11-10 01:50:56 +05:30
|
|
|
|
@property
|
2021-12-14 21:35:00 +05:30
|
|
|
|
def in_FZ(self) -> Union[np.bool_, np.ndarray]:
|
2020-11-10 01:50:56 +05:30
|
|
|
|
"""
|
|
|
|
|
Check whether orientation falls into fundamental zone of own symmetry.
|
|
|
|
|
|
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-02-11 22:40:29 +05:30
|
|
|
|
in : numpy.ndarray of bool, shape (self.shape)
|
2022-01-13 03:43:38 +05:30
|
|
|
|
Whether Rodrigues-Frank vector falls into fundamental zone.
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
Notes
|
|
|
|
|
-----
|
|
|
|
|
Fundamental zones in Rodrigues space are point-symmetric around origin.
|
|
|
|
|
|
|
|
|
|
References
|
|
|
|
|
----------
|
|
|
|
|
A. Heinz and P. Neumann, Acta Crystallographica Section A 47:780-789, 1991
|
|
|
|
|
https://doi.org/10.1107/S0108767391006864
|
|
|
|
|
|
|
|
|
|
"""
|
2021-02-10 04:04:51 +05:30
|
|
|
|
rho_abs = np.abs(self.as_Rodrigues_vector(compact=True))*(1.-1.e-9)
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
with np.errstate(invalid='ignore'):
|
|
|
|
|
# using '*'/prod for 'and'
|
|
|
|
|
if self.family == 'cubic':
|
|
|
|
|
return (np.prod(np.sqrt(2)-1. >= rho_abs,axis=-1) *
|
2021-02-10 04:04:51 +05:30
|
|
|
|
(1. >= np.sum(rho_abs,axis=-1))).astype(bool)
|
2023-02-21 20:57:06 +05:30
|
|
|
|
if self.family == 'hexagonal':
|
2020-11-10 01:50:56 +05:30
|
|
|
|
return (np.prod(1. >= rho_abs,axis=-1) *
|
|
|
|
|
(2. >= np.sqrt(3)*rho_abs[...,0] + rho_abs[...,1]) *
|
|
|
|
|
(2. >= np.sqrt(3)*rho_abs[...,1] + rho_abs[...,0]) *
|
2021-02-10 04:04:51 +05:30
|
|
|
|
(2. >= np.sqrt(3) + rho_abs[...,2])).astype(bool)
|
2023-02-21 20:57:06 +05:30
|
|
|
|
if self.family == 'tetragonal':
|
2020-11-10 01:50:56 +05:30
|
|
|
|
return (np.prod(1. >= rho_abs[...,:2],axis=-1) *
|
|
|
|
|
(np.sqrt(2) >= rho_abs[...,0] + rho_abs[...,1]) *
|
2021-02-10 04:04:51 +05:30
|
|
|
|
(np.sqrt(2) >= rho_abs[...,2] + 1.)).astype(bool)
|
2023-02-21 20:57:06 +05:30
|
|
|
|
if self.family == 'orthorhombic':
|
2021-02-10 04:04:51 +05:30
|
|
|
|
return (np.prod(1. >= rho_abs,axis=-1)).astype(bool)
|
2023-02-21 20:57:06 +05:30
|
|
|
|
if self.family == 'monoclinic':
|
2022-06-10 02:08:13 +05:30
|
|
|
|
return np.logical_or( 1. >= rho_abs[...,1],
|
|
|
|
|
np.isnan(rho_abs[...,1]))
|
2023-02-21 20:57:06 +05:30
|
|
|
|
if self.family == 'triclinic':
|
2022-06-10 02:08:13 +05:30
|
|
|
|
return np.ones(rho_abs.shape[:-1]).astype(bool)
|
2023-02-21 20:57:06 +05:30
|
|
|
|
|
|
|
|
|
raise TypeError(f'unknown symmetry "{self.family}"')
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
|
|
|
|
|
@property
|
2021-12-14 21:35:00 +05:30
|
|
|
|
def in_disorientation_FZ(self) -> np.ndarray:
|
2020-11-10 01:50:56 +05:30
|
|
|
|
"""
|
|
|
|
|
Check whether orientation falls into fundamental zone of disorientations.
|
|
|
|
|
|
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-02-11 22:40:29 +05:30
|
|
|
|
in : numpy.ndarray of bool, shape (self.shape)
|
2022-01-13 03:43:38 +05:30
|
|
|
|
Whether Rodrigues-Frank vector falls into disorientation FZ.
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
References
|
|
|
|
|
----------
|
|
|
|
|
A. Heinz and P. Neumann, Acta Crystallographica Section A 47:780-789, 1991
|
|
|
|
|
https://doi.org/10.1107/S0108767391006864
|
|
|
|
|
|
|
|
|
|
"""
|
2021-02-10 04:04:51 +05:30
|
|
|
|
rho = self.as_Rodrigues_vector(compact=True)*(1.0-1.0e-9)
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
with np.errstate(invalid='ignore'):
|
|
|
|
|
if self.family == 'cubic':
|
|
|
|
|
return ((rho[...,0] >= rho[...,1]) &
|
|
|
|
|
(rho[...,1] >= rho[...,2]) &
|
2021-02-10 04:04:51 +05:30
|
|
|
|
(rho[...,2] >= 0)).astype(bool)
|
2023-02-21 20:57:06 +05:30
|
|
|
|
if self.family == 'hexagonal':
|
2020-11-10 01:50:56 +05:30
|
|
|
|
return ((rho[...,0] >= rho[...,1]*np.sqrt(3)) &
|
|
|
|
|
(rho[...,1] >= 0) &
|
2021-02-10 04:04:51 +05:30
|
|
|
|
(rho[...,2] >= 0)).astype(bool)
|
2023-02-21 20:57:06 +05:30
|
|
|
|
if self.family == 'tetragonal':
|
2020-11-10 01:50:56 +05:30
|
|
|
|
return ((rho[...,0] >= rho[...,1]) &
|
|
|
|
|
(rho[...,1] >= 0) &
|
2021-02-10 04:04:51 +05:30
|
|
|
|
(rho[...,2] >= 0)).astype(bool)
|
2023-02-21 20:57:06 +05:30
|
|
|
|
if self.family == 'orthorhombic':
|
2020-11-10 01:50:56 +05:30
|
|
|
|
return ((rho[...,0] >= 0) &
|
|
|
|
|
(rho[...,1] >= 0) &
|
2021-02-10 04:04:51 +05:30
|
|
|
|
(rho[...,2] >= 0)).astype(bool)
|
2023-02-21 20:57:06 +05:30
|
|
|
|
if self.family == 'monoclinic':
|
2020-11-10 01:50:56 +05:30
|
|
|
|
return ((rho[...,1] >= 0) &
|
2021-02-10 04:04:51 +05:30
|
|
|
|
(rho[...,2] >= 0)).astype(bool)
|
2023-02-21 20:57:06 +05:30
|
|
|
|
|
|
|
|
|
return np.ones_like(rho[...,0],dtype=bool)
|
2020-06-19 02:23:04 +05:30
|
|
|
|
|
2022-01-26 19:39:09 +05:30
|
|
|
|
def disorientation(self,
|
2022-02-04 21:27:25 +05:30
|
|
|
|
other: 'Orientation',
|
2022-02-04 14:27:42 +05:30
|
|
|
|
return_operators: bool = False) -> object:
|
2020-11-10 01:50:56 +05:30
|
|
|
|
"""
|
2023-02-21 20:57:06 +05:30
|
|
|
|
Calculate disorientation between self and given other orientation.
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
|
|
|
|
other : Orientation
|
|
|
|
|
Orientation to calculate disorientation for.
|
|
|
|
|
Shape of other blends with shape of own rotation array.
|
2023-02-21 20:57:06 +05:30
|
|
|
|
For example, shapes of (2,3) for own rotations
|
|
|
|
|
and (3,2) for other's result in (2,3,2) disorientations.
|
2020-11-10 01:50:56 +05:30
|
|
|
|
return_operators : bool, optional
|
2023-02-21 20:57:06 +05:30
|
|
|
|
Return index pair of symmetrically equivalent orientations
|
|
|
|
|
that result in disorientation axis falling into FZ.
|
2020-11-10 01:50:56 +05:30
|
|
|
|
Defaults to False.
|
|
|
|
|
|
|
|
|
|
Returns
|
|
|
|
|
-------
|
|
|
|
|
disorientation : Orientation
|
|
|
|
|
Disorientation between self and other.
|
2022-02-11 22:40:29 +05:30
|
|
|
|
operators : numpy.ndarray of int, shape (...,2), conditional
|
2020-11-10 01:50:56 +05:30
|
|
|
|
Index of symmetrically equivalent orientation that rotated vector to the SST.
|
|
|
|
|
|
|
|
|
|
Notes
|
|
|
|
|
-----
|
2022-11-19 13:40:00 +05:30
|
|
|
|
Requires same crystal family for both orientations.
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
2020-11-15 14:52:01 +05:30
|
|
|
|
Examples
|
|
|
|
|
--------
|
|
|
|
|
Disorientation between two specific orientations of hexagonal symmetry:
|
|
|
|
|
|
|
|
|
|
>>> import damask
|
2021-08-16 22:53:31 +05:30
|
|
|
|
>>> a = damask.Orientation.from_Euler_angles(phi=[123,32,21],degrees=True,family='hexagonal')
|
|
|
|
|
>>> b = damask.Orientation.from_Euler_angles(phi=[104,11,87],degrees=True,family='hexagonal')
|
2020-11-15 14:52:01 +05:30
|
|
|
|
>>> a.disorientation(b)
|
|
|
|
|
Crystal family hexagonal
|
|
|
|
|
Quaternion: (real=0.976, imag=<+0.189, +0.018, +0.103>)
|
|
|
|
|
Matrix:
|
|
|
|
|
[[ 0.97831006 0.20710935 0.00389135]
|
|
|
|
|
[-0.19363288 0.90765544 0.37238141]
|
|
|
|
|
[ 0.07359167 -0.36505797 0.92807163]]
|
|
|
|
|
Bunge Eulers / deg: (11.40, 21.86, 0.60)
|
|
|
|
|
|
2022-01-11 05:06:49 +05:30
|
|
|
|
Plot a sample from the Mackenzie distribution.
|
|
|
|
|
|
|
|
|
|
>>> import matplotlib.pyplot as plt
|
|
|
|
|
>>> import damask
|
|
|
|
|
>>> N = 10000
|
|
|
|
|
>>> a = damask.Orientation.from_random(shape=N,family='cubic')
|
|
|
|
|
>>> b = damask.Orientation.from_random(shape=N,family='cubic')
|
2023-02-21 20:57:06 +05:30
|
|
|
|
>>> n,omega = a.disorientation(b).as_axis_angle(degrees=True,pair=True)
|
|
|
|
|
>>> plt.hist(omega,25)
|
2022-01-11 05:06:49 +05:30
|
|
|
|
>>> plt.show()
|
|
|
|
|
|
2020-11-10 01:50:56 +05:30
|
|
|
|
"""
|
2022-11-19 13:40:00 +05:30
|
|
|
|
# For extension to cases with differing symmetry see
|
|
|
|
|
# https://doi.org/10.1107/S0021889808016373 and https://doi.org/10.1107/S0108767391006864
|
2020-11-10 01:50:56 +05:30
|
|
|
|
if self.family != other.family:
|
2020-11-20 03:06:19 +05:30
|
|
|
|
raise NotImplementedError('disorientation between different crystal families')
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
2022-02-03 20:41:09 +05:30
|
|
|
|
blend = util.shapeblender(self.shape,other.shape)
|
2021-12-14 21:35:00 +05:30
|
|
|
|
s = self.equivalent
|
2020-11-10 01:50:56 +05:30
|
|
|
|
o = other.equivalent
|
|
|
|
|
|
2022-02-13 15:11:10 +05:30
|
|
|
|
s_ = s.reshape((s.shape[0],1)+ self.shape).broadcast_to((s.shape[0],o.shape[0])+blend,mode='right')
|
|
|
|
|
o_ = o.reshape((1,o.shape[0])+other.shape).broadcast_to((s.shape[0],o.shape[0])+blend,mode='right')
|
2020-11-10 01:50:56 +05:30
|
|
|
|
r_ = s_.misorientation(o_)
|
|
|
|
|
_r = ~r_
|
|
|
|
|
|
|
|
|
|
forward = r_.in_FZ & r_.in_disorientation_FZ
|
|
|
|
|
reverse = _r.in_FZ & _r.in_disorientation_FZ
|
|
|
|
|
ok = forward | reverse
|
|
|
|
|
ok &= (np.cumsum(ok.reshape((-1,)+ok.shape[2:]),axis=0) == 1).reshape(ok.shape)
|
|
|
|
|
r = np.where(np.any(forward[...,np.newaxis],axis=(0,1),keepdims=True),
|
|
|
|
|
r_.quaternion,
|
|
|
|
|
_r.quaternion)
|
2022-02-04 14:27:42 +05:30
|
|
|
|
loc = np.where(ok)
|
|
|
|
|
sort = 0 if len(loc) == 2 else np.lexsort(loc[:1:-1])
|
|
|
|
|
quat = r[ok][sort].reshape(blend+(4,))
|
2022-02-06 21:41:18 +05:30
|
|
|
|
|
2020-11-10 01:50:56 +05:30
|
|
|
|
return (
|
|
|
|
|
(self.copy(rotation=quat),
|
|
|
|
|
(np.vstack(loc[:2]).T)[sort].reshape(blend+(2,)))
|
|
|
|
|
if return_operators else
|
|
|
|
|
self.copy(rotation=quat)
|
|
|
|
|
)
|
|
|
|
|
|
|
|
|
|
|
2022-01-26 19:39:09 +05:30
|
|
|
|
def average(self,
|
2022-11-23 02:56:15 +05:30
|
|
|
|
weights: Optional[FloatSequence] = None,
|
2022-02-11 22:40:29 +05:30
|
|
|
|
return_cloud: bool = False):
|
2020-11-10 01:50:56 +05:30
|
|
|
|
"""
|
|
|
|
|
Return orientation average over last dimension.
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-02-11 22:40:29 +05:30
|
|
|
|
weights : numpy.ndarray, shape (self.shape), optional
|
2020-11-10 01:50:56 +05:30
|
|
|
|
Relative weights of orientations.
|
2023-02-21 20:57:06 +05:30
|
|
|
|
Defaults to equal weights.
|
2020-11-10 01:50:56 +05:30
|
|
|
|
return_cloud : bool, optional
|
2022-12-06 04:59:03 +05:30
|
|
|
|
Return the specific (symmetrically equivalent) orientations that were averaged.
|
2020-11-10 01:50:56 +05:30
|
|
|
|
Defaults to False.
|
|
|
|
|
|
|
|
|
|
Returns
|
|
|
|
|
-------
|
|
|
|
|
average : Orientation
|
|
|
|
|
Weighted average of original Orientation field.
|
|
|
|
|
cloud : Orientations, conditional
|
2022-12-06 04:59:03 +05:30
|
|
|
|
Symmetrically equivalent version of each orientation that were actually used in averaging.
|
2015-03-28 13:13:49 +05:30
|
|
|
|
|
2020-11-10 01:50:56 +05:30
|
|
|
|
References
|
|
|
|
|
----------
|
2021-03-18 20:06:40 +05:30
|
|
|
|
J.C. Glez and J. Driver, Journal of Applied Crystallography 34:280-288, 2001
|
2020-11-10 01:50:56 +05:30
|
|
|
|
https://doi.org/10.1107/S0021889801003077
|
2020-02-21 03:33:37 +05:30
|
|
|
|
|
2020-11-10 01:50:56 +05:30
|
|
|
|
"""
|
2020-06-30 17:25:09 +05:30
|
|
|
|
eq = self.equivalent
|
2022-02-11 22:40:29 +05:30
|
|
|
|
m = eq.misorientation(self[...,0].reshape((1,)+self.shape[:-1]+(1,)) # type: ignore
|
|
|
|
|
.broadcast_to(eq.shape)).as_axis_angle()[...,3] # type: ignore
|
2020-11-10 01:50:56 +05:30
|
|
|
|
r = Rotation(np.squeeze(np.take_along_axis(eq.quaternion,
|
|
|
|
|
np.argmin(m,axis=0)[np.newaxis,...,np.newaxis],
|
|
|
|
|
axis=0),
|
|
|
|
|
axis=0))
|
2022-02-13 05:54:02 +05:30
|
|
|
|
return ((self.copy(Rotation(r).average(weights)),self.copy(Rotation(r))) if return_cloud else
|
|
|
|
|
self.copy(Rotation(r).average(weights))
|
2020-11-10 01:50:56 +05:30
|
|
|
|
)
|
|
|
|
|
|
|
|
|
|
|
2021-12-14 21:35:00 +05:30
|
|
|
|
def to_SST(self,
|
2022-02-11 22:40:29 +05:30
|
|
|
|
vector: FloatSequence,
|
2021-12-14 21:35:00 +05:30
|
|
|
|
proper: bool = False,
|
|
|
|
|
return_operators: bool = False) -> np.ndarray:
|
2020-11-10 01:50:56 +05:30
|
|
|
|
"""
|
2022-12-06 04:59:03 +05:30
|
|
|
|
Rotate lab frame vector to ensure it falls into (improper or proper) standard stereographic triangle of crystal symmetry.
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-02-11 22:40:29 +05:30
|
|
|
|
vector : numpy.ndarray, shape (...,3)
|
2020-11-10 01:50:56 +05:30
|
|
|
|
Lab frame vector to align with crystal frame direction.
|
2021-08-16 22:53:31 +05:30
|
|
|
|
Shape of vector blends with shape of own rotation array.
|
2022-02-11 22:40:29 +05:30
|
|
|
|
For example, a rotation array of shape (3,2) and a vector array of shape (2,4) result in (3,2,4) outputs.
|
2020-11-10 01:50:56 +05:30
|
|
|
|
proper : bool, optional
|
|
|
|
|
Consider only vectors with z >= 0, hence combine two neighboring SSTs.
|
|
|
|
|
Defaults to False.
|
|
|
|
|
return_operators : bool, optional
|
|
|
|
|
Return the symmetrically equivalent orientation that rotated vector to SST.
|
|
|
|
|
Defaults to False.
|
|
|
|
|
|
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-02-11 22:40:29 +05:30
|
|
|
|
vector_SST : numpy.ndarray, shape (...,3)
|
2020-11-10 01:50:56 +05:30
|
|
|
|
Rotated vector falling into SST.
|
2022-12-06 04:59:03 +05:30
|
|
|
|
operator : numpy.ndarray of int, shape (...), conditional
|
2020-11-10 01:50:56 +05:30
|
|
|
|
Index of symmetrically equivalent orientation that rotated vector to SST.
|
|
|
|
|
|
|
|
|
|
"""
|
2022-02-11 22:40:29 +05:30
|
|
|
|
vector_ = np.array(vector,float)
|
|
|
|
|
if vector_.shape[-1] != 3:
|
|
|
|
|
raise ValueError('input is not a field of three-dimensional vectors')
|
2020-11-10 01:50:56 +05:30
|
|
|
|
eq = self.equivalent
|
2022-02-11 22:40:29 +05:30
|
|
|
|
blend = util.shapeblender(eq.shape,vector_.shape[:-1])
|
2022-02-13 15:11:10 +05:30
|
|
|
|
poles = eq.broadcast_to(blend,mode='right') @ np.broadcast_to(vector_,blend+(3,))
|
2020-11-10 01:50:56 +05:30
|
|
|
|
ok = self.in_SST(poles,proper=proper)
|
|
|
|
|
ok &= np.cumsum(ok,axis=0) == 1
|
|
|
|
|
loc = np.where(ok)
|
|
|
|
|
sort = 0 if len(loc) == 1 else np.lexsort(loc[:0:-1])
|
|
|
|
|
return (
|
|
|
|
|
(poles[ok][sort].reshape(blend[1:]+(3,)), (np.vstack(loc[:1]).T)[sort].reshape(blend[1:]))
|
|
|
|
|
if return_operators else
|
|
|
|
|
poles[ok][sort].reshape(blend[1:]+(3,))
|
|
|
|
|
)
|
|
|
|
|
|
|
|
|
|
|
2022-01-26 19:39:09 +05:30
|
|
|
|
def in_SST(self,
|
2022-02-11 22:40:29 +05:30
|
|
|
|
vector: FloatSequence,
|
2022-01-26 19:39:09 +05:30
|
|
|
|
proper: bool = False) -> Union[np.bool_, np.ndarray]:
|
2021-06-02 20:55:24 +05:30
|
|
|
|
"""
|
|
|
|
|
Check whether given crystal frame vector falls into standard stereographic triangle of own symmetry.
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-02-11 22:40:29 +05:30
|
|
|
|
vector : numpy.ndarray, shape (...,3)
|
2021-06-02 20:55:24 +05:30
|
|
|
|
Vector to check.
|
|
|
|
|
proper : bool, optional
|
|
|
|
|
Consider only vectors with z >= 0, hence combine two neighboring SSTs.
|
|
|
|
|
Defaults to False.
|
|
|
|
|
|
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-01-13 03:43:38 +05:30
|
|
|
|
in : numpy.ndarray, shape (...)
|
|
|
|
|
Whether vector falls into SST.
|
2021-06-02 20:55:24 +05:30
|
|
|
|
|
|
|
|
|
"""
|
2022-02-11 22:40:29 +05:30
|
|
|
|
vector_ = np.array(vector,float)
|
|
|
|
|
if vector_.shape[-1] != 3:
|
2021-06-02 20:55:24 +05:30
|
|
|
|
raise ValueError('input is not a field of three-dimensional vectors')
|
|
|
|
|
|
2021-07-18 20:09:52 +05:30
|
|
|
|
if self.standard_triangle is None: # direct exit for no symmetry
|
2022-02-11 22:40:29 +05:30
|
|
|
|
return np.ones_like(vector_[...,0],bool)
|
2021-06-02 20:55:24 +05:30
|
|
|
|
|
|
|
|
|
if proper:
|
|
|
|
|
components_proper = np.around(np.einsum('...ji,...i',
|
2022-02-11 22:40:29 +05:30
|
|
|
|
np.broadcast_to(self.standard_triangle['proper'], vector_.shape+(3,)),
|
|
|
|
|
vector_), 12)
|
2021-06-02 20:55:24 +05:30
|
|
|
|
components_improper = np.around(np.einsum('...ji,...i',
|
2022-02-11 22:40:29 +05:30
|
|
|
|
np.broadcast_to(self.standard_triangle['improper'], vector_.shape+(3,)),
|
|
|
|
|
vector_), 12)
|
2021-06-02 20:55:24 +05:30
|
|
|
|
return np.all(components_proper >= 0.0,axis=-1) \
|
|
|
|
|
| np.all(components_improper >= 0.0,axis=-1)
|
|
|
|
|
else:
|
|
|
|
|
components = np.around(np.einsum('...ji,...i',
|
2022-02-11 22:40:29 +05:30
|
|
|
|
np.broadcast_to(self.standard_triangle['improper'], vector_.shape+(3,)),
|
|
|
|
|
np.block([vector_[...,:2],np.abs(vector_[...,2:3])])), 12)
|
2021-06-02 20:55:24 +05:30
|
|
|
|
|
|
|
|
|
return np.all(components >= 0.0,axis=-1)
|
|
|
|
|
|
|
|
|
|
|
2022-01-26 19:39:09 +05:30
|
|
|
|
def IPF_color(self,
|
2022-02-11 22:40:29 +05:30
|
|
|
|
vector: FloatSequence,
|
2022-01-26 19:39:09 +05:30
|
|
|
|
in_SST: bool = True,
|
|
|
|
|
proper: bool = False) -> np.ndarray:
|
2021-06-02 20:55:24 +05:30
|
|
|
|
"""
|
2022-12-06 04:59:03 +05:30
|
|
|
|
Map lab frame vector to RGB color within standard stereographic triangle of own symmetry.
|
2021-06-02 20:55:24 +05:30
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-02-11 22:40:29 +05:30
|
|
|
|
vector : numpy.ndarray, shape (...,3)
|
2022-12-06 04:59:03 +05:30
|
|
|
|
Lab frame vector to colorize.
|
2021-08-16 22:53:31 +05:30
|
|
|
|
Shape of vector blends with shape of own rotation array.
|
2022-02-11 22:40:29 +05:30
|
|
|
|
For example, a rotation array of shape (3,2) and a vector array of shape (2,4) result in (3,2,4) outputs.
|
2021-06-02 20:55:24 +05:30
|
|
|
|
in_SST : bool, optional
|
|
|
|
|
Consider symmetrically equivalent orientations such that poles are located in SST.
|
|
|
|
|
Defaults to True.
|
|
|
|
|
proper : bool, optional
|
|
|
|
|
Consider only vectors with z >= 0, hence combine two neighboring SSTs (with mirrored colors).
|
|
|
|
|
Defaults to False.
|
|
|
|
|
|
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-02-11 22:40:29 +05:30
|
|
|
|
rgb : numpy.ndarray, shape (...,3)
|
2021-06-02 20:55:24 +05:30
|
|
|
|
RGB array of IPF colors.
|
|
|
|
|
|
|
|
|
|
Examples
|
|
|
|
|
--------
|
2022-12-06 04:59:03 +05:30
|
|
|
|
Inverse pole figure color of the e_3 lab direction for a
|
|
|
|
|
crystal in "Cube" orientation with cubic symmetry:
|
2021-06-02 20:55:24 +05:30
|
|
|
|
|
2021-07-25 23:01:48 +05:30
|
|
|
|
>>> import damask
|
2021-08-16 22:53:31 +05:30
|
|
|
|
>>> o = damask.Orientation(family='cubic')
|
2021-06-02 20:55:24 +05:30
|
|
|
|
>>> o.IPF_color([0,0,1])
|
|
|
|
|
array([1., 0., 0.])
|
|
|
|
|
|
2022-12-03 03:25:22 +05:30
|
|
|
|
Sample standard triangle for hexagonal symmetry:
|
|
|
|
|
|
|
|
|
|
>>> import damask
|
|
|
|
|
>>> from matplotlib import pyplot as plt
|
2022-12-06 04:59:03 +05:30
|
|
|
|
>>> lab = [0,0,1]
|
2022-12-03 03:25:22 +05:30
|
|
|
|
>>> o = damask.Orientation.from_random(shape=500000,family='hexagonal')
|
2022-12-06 04:59:03 +05:30
|
|
|
|
>>> coord = damask.util.project_equal_area(o.to_SST(lab))
|
|
|
|
|
>>> color = o.IPF_color(lab)
|
2022-12-03 03:25:22 +05:30
|
|
|
|
>>> plt.scatter(coord[:,0],coord[:,1],color=color,s=.06)
|
|
|
|
|
>>> plt.axis('scaled')
|
|
|
|
|
>>> plt.show()
|
|
|
|
|
|
2021-06-02 20:55:24 +05:30
|
|
|
|
"""
|
|
|
|
|
if np.array(vector).shape[-1] != 3:
|
|
|
|
|
raise ValueError('input is not a field of three-dimensional vectors')
|
|
|
|
|
|
|
|
|
|
vector_ = self.to_SST(vector,proper) if in_SST else \
|
|
|
|
|
self @ np.broadcast_to(vector,self.shape+(3,))
|
|
|
|
|
|
2021-07-18 20:09:52 +05:30
|
|
|
|
if self.standard_triangle is None: # direct exit for no symmetry
|
2021-06-02 20:55:24 +05:30
|
|
|
|
return np.zeros_like(vector_)
|
|
|
|
|
|
|
|
|
|
if proper:
|
|
|
|
|
components_proper = np.around(np.einsum('...ji,...i',
|
2022-02-11 22:40:29 +05:30
|
|
|
|
np.broadcast_to(self.standard_triangle['proper'], vector_.shape+(3,)),
|
2021-06-02 20:55:24 +05:30
|
|
|
|
vector_), 12)
|
|
|
|
|
components_improper = np.around(np.einsum('...ji,...i',
|
2021-07-18 20:09:52 +05:30
|
|
|
|
np.broadcast_to(self.standard_triangle['improper'], vector_.shape+(3,)),
|
2021-06-02 20:55:24 +05:30
|
|
|
|
vector_), 12)
|
2021-12-14 21:35:00 +05:30
|
|
|
|
in_SST_ = np.all(components_proper >= 0.0,axis=-1) \
|
2022-02-03 20:41:09 +05:30
|
|
|
|
| np.all(components_improper >= 0.0,axis=-1)
|
2021-12-14 21:35:00 +05:30
|
|
|
|
components = np.where((in_SST_ & np.all(components_proper >= 0.0,axis=-1))[...,np.newaxis],
|
2021-06-02 20:55:24 +05:30
|
|
|
|
components_proper,components_improper)
|
|
|
|
|
else:
|
|
|
|
|
components = np.around(np.einsum('...ji,...i',
|
2021-07-18 20:09:52 +05:30
|
|
|
|
np.broadcast_to(self .standard_triangle['improper'], vector_.shape+(3,)),
|
2021-06-02 20:55:24 +05:30
|
|
|
|
np.block([vector_[...,:2],np.abs(vector_[...,2:3])])), 12)
|
|
|
|
|
|
2021-12-14 21:35:00 +05:30
|
|
|
|
in_SST_ = np.all(components >= 0.0,axis=-1)
|
2021-06-02 20:55:24 +05:30
|
|
|
|
|
|
|
|
|
with np.errstate(invalid='ignore',divide='ignore'):
|
2022-12-03 03:25:22 +05:30
|
|
|
|
rgb = (components/np.linalg.norm(components,axis=-1,keepdims=True))**(1./3.) # smoothen color ramps
|
2021-06-02 20:55:24 +05:30
|
|
|
|
rgb = np.clip(rgb,0.,1.) # clip intensity
|
|
|
|
|
rgb /= np.max(rgb,axis=-1,keepdims=True) # normalize to (HS)V = 1
|
2021-12-14 21:35:00 +05:30
|
|
|
|
rgb[np.broadcast_to(~in_SST_[...,np.newaxis],rgb.shape)] = 0.0
|
2021-06-02 20:55:24 +05:30
|
|
|
|
|
|
|
|
|
return rgb
|
|
|
|
|
|
|
|
|
|
|
2021-07-13 03:44:13 +05:30
|
|
|
|
####################################################################################################
|
2021-06-02 20:55:24 +05:30
|
|
|
|
# functions that require lattice, not just family
|
|
|
|
|
|
2021-12-14 21:35:00 +05:30
|
|
|
|
def to_pole(self, *,
|
2022-11-23 02:56:15 +05:30
|
|
|
|
uvw: Optional[FloatSequence] = None,
|
|
|
|
|
hkl: Optional[FloatSequence] = None,
|
2022-05-11 18:25:55 +05:30
|
|
|
|
with_symmetry: bool = False,
|
|
|
|
|
normalize: bool = True) -> np.ndarray:
|
2020-11-10 01:50:56 +05:30
|
|
|
|
"""
|
|
|
|
|
Calculate lab frame vector along lattice direction [uvw] or plane normal (hkl).
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-02-11 22:40:29 +05:30
|
|
|
|
uvw|hkl : numpy.ndarray, shape (...,3)
|
2020-11-10 01:50:56 +05:30
|
|
|
|
Miller indices of crystallographic direction or plane normal.
|
2021-08-16 22:53:31 +05:30
|
|
|
|
Shape of vector blends with shape of own rotation array.
|
2022-05-13 13:23:44 +05:30
|
|
|
|
For example, a rotation array of shape (3,2) and a vector
|
|
|
|
|
array of shape (2,4) result in (3,2,4) outputs.
|
2020-11-10 01:50:56 +05:30
|
|
|
|
with_symmetry : bool, optional
|
|
|
|
|
Calculate all N symmetrically equivalent vectors.
|
2022-05-11 18:25:55 +05:30
|
|
|
|
Defaults to False.
|
|
|
|
|
normalize : bool, optional
|
|
|
|
|
Normalize output vector.
|
|
|
|
|
Defaults to True.
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-02-11 22:40:29 +05:30
|
|
|
|
vector : numpy.ndarray, shape (...,3) or (...,N,3)
|
2022-05-13 13:23:44 +05:30
|
|
|
|
Lab frame vector (or vectors if with_symmetry) along
|
|
|
|
|
[uvw] direction or (hkl) plane normal.
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
"""
|
2021-06-06 23:19:29 +05:30
|
|
|
|
v = self.to_frame(uvw=uvw,hkl=hkl)
|
2021-08-16 22:53:31 +05:30
|
|
|
|
blend = util.shapeblender(self.shape,v.shape[:-1])
|
2022-05-11 18:25:55 +05:30
|
|
|
|
if normalize:
|
|
|
|
|
v /= np.linalg.norm(v,axis=-1,keepdims=len(v.shape)>1)
|
2021-06-02 12:58:27 +05:30
|
|
|
|
if with_symmetry:
|
2021-07-18 21:33:36 +05:30
|
|
|
|
sym_ops = self.symmetry_operations
|
2021-08-16 22:53:31 +05:30
|
|
|
|
shape = v.shape[:-1]+sym_ops.shape
|
|
|
|
|
blend += sym_ops.shape
|
|
|
|
|
v = sym_ops.broadcast_to(shape) \
|
|
|
|
|
@ np.broadcast_to(v.reshape(util.shapeshifter(v.shape,shape+(3,))),shape+(3,))
|
2023-09-20 03:33:55 +05:30
|
|
|
|
return ~(self.broadcast_to(blend))@np.broadcast_to(v,blend+(3,))
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
|
2021-12-14 21:35:00 +05:30
|
|
|
|
def Schmid(self, *,
|
2022-11-23 02:56:15 +05:30
|
|
|
|
N_slip: Optional[IntSequence] = None,
|
|
|
|
|
N_twin: Optional[IntSequence] = None) -> np.ndarray:
|
2020-11-10 01:50:56 +05:30
|
|
|
|
u"""
|
2021-08-04 21:15:25 +05:30
|
|
|
|
Calculate Schmid matrix P = d ⨂ n in the lab frame for selected deformation systems.
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-02-26 22:10:12 +05:30
|
|
|
|
N_slip|N_twin : '*' or sequence of int
|
2021-08-08 14:14:38 +05:30
|
|
|
|
Number of deformation systems per family of the deformation system.
|
2021-08-09 03:26:54 +05:30
|
|
|
|
Use '*' to select all.
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-02-11 22:40:29 +05:30
|
|
|
|
P : numpy.ndarray, shape (N,...,3,3)
|
2020-11-10 01:50:56 +05:30
|
|
|
|
Schmid matrix for each of the N deformation systems.
|
|
|
|
|
|
2020-11-15 14:52:01 +05:30
|
|
|
|
Examples
|
|
|
|
|
--------
|
2021-08-08 14:14:38 +05:30
|
|
|
|
Schmid matrix (in lab frame) of first octahedral slip system of a face-centered
|
2020-11-15 14:52:01 +05:30
|
|
|
|
cubic crystal in "Goss" orientation.
|
|
|
|
|
|
|
|
|
|
>>> import numpy as np
|
2023-02-21 20:57:06 +05:30
|
|
|
|
>>> import damask
|
2020-11-15 14:52:01 +05:30
|
|
|
|
>>> np.set_printoptions(3,suppress=True,floatmode='fixed')
|
2021-08-04 21:15:25 +05:30
|
|
|
|
>>> O = damask.Orientation.from_Euler_angles(phi=[0,45,0],degrees=True,lattice='cF')
|
2023-09-20 03:33:55 +05:30
|
|
|
|
>>> O.Schmid(N_slip=[12])[0]
|
2020-11-15 14:52:01 +05:30
|
|
|
|
array([[ 0.000, 0.000, 0.000],
|
|
|
|
|
[ 0.577, -0.000, 0.816],
|
|
|
|
|
[ 0.000, 0.000, 0.000]])
|
2020-11-15 17:36:26 +05:30
|
|
|
|
|
2020-11-10 01:50:56 +05:30
|
|
|
|
"""
|
2021-08-04 21:15:25 +05:30
|
|
|
|
if (N_slip is not None) ^ (N_twin is None):
|
2022-02-22 21:12:05 +05:30
|
|
|
|
raise KeyError('specify either "N_slip" or "N_twin"')
|
2021-08-04 21:15:25 +05:30
|
|
|
|
|
2021-08-09 03:26:54 +05:30
|
|
|
|
kinematics,active = (self.kinematics('slip'),N_slip) if N_twin is None else \
|
|
|
|
|
(self.kinematics('twin'),N_twin)
|
2021-08-08 14:14:38 +05:30
|
|
|
|
if active == '*': active = [len(a) for a in kinematics['direction']]
|
2021-08-04 21:15:25 +05:30
|
|
|
|
|
2022-02-03 20:41:09 +05:30
|
|
|
|
if not active:
|
2022-02-13 05:54:02 +05:30
|
|
|
|
raise ValueError('Schmid matrix not defined')
|
2021-08-04 21:15:25 +05:30
|
|
|
|
d = self.to_frame(uvw=np.vstack([kinematics['direction'][i][:n] for i,n in enumerate(active)]))
|
|
|
|
|
p = self.to_frame(hkl=np.vstack([kinematics['plane'][i][:n] for i,n in enumerate(active)]))
|
2021-07-25 23:01:48 +05:30
|
|
|
|
P = np.einsum('...i,...j',d/np.linalg.norm(d,axis=1,keepdims=True),
|
|
|
|
|
p/np.linalg.norm(p,axis=1,keepdims=True))
|
2020-11-10 01:50:56 +05:30
|
|
|
|
|
2021-07-26 00:02:24 +05:30
|
|
|
|
shape = P.shape[0:1]+self.shape+(3,3)
|
|
|
|
|
return ~self.broadcast_to(shape[:-2]) \
|
|
|
|
|
@ np.broadcast_to(P.reshape(util.shapeshifter(P.shape,shape)),shape)
|
2021-06-02 20:55:24 +05:30
|
|
|
|
|
|
|
|
|
|
2022-02-13 05:54:02 +05:30
|
|
|
|
def related(self: MyType,
|
|
|
|
|
model: str) -> MyType:
|
2021-06-02 20:55:24 +05:30
|
|
|
|
"""
|
2022-03-10 04:54:05 +05:30
|
|
|
|
All orientations related to self by given relationship model.
|
2021-06-02 20:55:24 +05:30
|
|
|
|
|
2022-02-26 22:10:12 +05:30
|
|
|
|
Parameters
|
|
|
|
|
----------
|
|
|
|
|
model : str
|
2022-03-10 04:54:05 +05:30
|
|
|
|
Orientation relationship model selected from self.orientation_relationships.
|
2022-02-26 22:10:12 +05:30
|
|
|
|
|
|
|
|
|
Returns
|
|
|
|
|
-------
|
2023-02-21 20:57:06 +05:30
|
|
|
|
rel : Orientation, shape (:,self.shape)
|
|
|
|
|
Orientations related to self according to the selected
|
|
|
|
|
model for the orientation relationship.
|
2022-02-26 22:10:12 +05:30
|
|
|
|
|
|
|
|
|
Examples
|
|
|
|
|
--------
|
2022-03-10 04:54:05 +05:30
|
|
|
|
Face-centered cubic orientations following from a
|
|
|
|
|
body-centered cubic crystal in "Cube" orientation according
|
|
|
|
|
to the Bain orientation relationship (cI -> cF).
|
2022-02-26 22:10:12 +05:30
|
|
|
|
|
|
|
|
|
>>> import numpy as np
|
|
|
|
|
>>> import damask
|
|
|
|
|
>>> np.set_printoptions(3,suppress=True,floatmode='fixed')
|
|
|
|
|
>>> damask.Orientation(lattice='cI').related('Bain')
|
|
|
|
|
Crystal family: cubic
|
|
|
|
|
Bravais lattice: cF
|
|
|
|
|
a=1 m, b=1 m, c=1 m
|
|
|
|
|
α=90°, β=90°, γ=90°
|
|
|
|
|
Quaternions of shape (3,)
|
|
|
|
|
[[0.924 0.383 0.000 0.000]
|
|
|
|
|
[0.924 0.000 0.383 0.000]
|
|
|
|
|
[0.924 0.000 0.000 0.383]]
|
2021-06-02 20:55:24 +05:30
|
|
|
|
|
|
|
|
|
"""
|
2021-06-06 23:19:29 +05:30
|
|
|
|
lattice,o = self.relation_operations(model)
|
2021-06-08 01:19:04 +05:30
|
|
|
|
target = Crystal(lattice=lattice)
|
2021-06-02 20:55:24 +05:30
|
|
|
|
o = o.broadcast_to(o.shape+self.shape,mode='right')
|
2021-06-03 13:48:27 +05:30
|
|
|
|
return Orientation(rotation=o*Rotation(self.quaternion).broadcast_to(o.shape,mode='left'),
|
|
|
|
|
lattice=lattice,
|
2021-06-06 23:19:29 +05:30
|
|
|
|
b = self.b if target.ratio['b'] is None else self.a*target.ratio['b'],
|
|
|
|
|
c = self.c if target.ratio['c'] is None else self.a*target.ratio['c'],
|
|
|
|
|
alpha = None if 'alpha' in target.immutable else self.alpha,
|
|
|
|
|
beta = None if 'beta' in target.immutable else self.beta,
|
|
|
|
|
gamma = None if 'gamma' in target.immutable else self.gamma,
|
2021-06-03 13:48:27 +05:30
|
|
|
|
)
|