polishing of help and style; relax to FloatSequence type where appropriate but keep doc at np.ndarray

This commit is contained in:
Philip Eisenlohr 2022-02-11 12:10:29 -05:00
parent 8c6225794d
commit 0a52ae3b6f
4 changed files with 113 additions and 91 deletions

View File

@ -131,9 +131,8 @@ class Crystal():
Crystal to check for equality. Crystal to check for equality.
""" """
if not isinstance(other, Crystal): return NotImplemented if not isinstance(other, Crystal) else \
return NotImplemented self.lattice == other.lattice and \
return self.lattice == other.lattice and \
self.parameters == other.parameters and \ self.parameters == other.parameters and \
self.family == other.family self.family == other.family
@ -316,8 +315,8 @@ class Crystal():
self.lattice[-1],None),dtype=float) self.lattice[-1],None),dtype=float)
def to_lattice(self, *, def to_lattice(self, *,
direction: np.ndarray = None, direction: FloatSequence = None,
plane: np.ndarray = None) -> np.ndarray: plane: FloatSequence = None) -> np.ndarray:
""" """
Calculate lattice vector corresponding to crystal frame direction or plane normal. Calculate lattice vector corresponding to crystal frame direction or plane normal.

View File

@ -189,7 +189,7 @@ class Orientation(Rotation,Crystal):
Returns Returns
------- -------
mask : numpy.ndarray bool mask : numpy.ndarray of bool
Mask indicating where corresponding orientations are close. Mask indicating where corresponding orientations are close.
""" """
@ -372,10 +372,10 @@ class Orientation(Rotation,Crystal):
Parameters Parameters
---------- ----------
uvw : list, numpy.ndarray of shape (...,3) uvw : numpy.ndarray, shape (...,3)
lattice direction aligned with lab frame x-direction. Lattice direction aligned with lab frame x-direction.
hkl : list, numpy.ndarray of shape (...,3) hkl : numpy.ndarray, shape (...,3)
lattice plane normal aligned with lab frame z-direction. Lattice plane normal aligned with lab frame z-direction.
""" """
o = cls(**kwargs) o = cls(**kwargs)
@ -417,7 +417,7 @@ class Orientation(Rotation,Crystal):
Returns Returns
------- -------
in : numpy.ndarray of bool, quaternion.shape in : numpy.ndarray of bool, shape (self.shape)
Whether Rodrigues-Frank vector falls into fundamental zone. Whether Rodrigues-Frank vector falls into fundamental zone.
Notes Notes
@ -461,7 +461,7 @@ class Orientation(Rotation,Crystal):
Returns Returns
------- -------
in : numpy.ndarray of bool, quaternion.shape in : numpy.ndarray of bool, shape (self.shape)
Whether Rodrigues-Frank vector falls into disorientation FZ. Whether Rodrigues-Frank vector falls into disorientation FZ.
References References
@ -515,7 +515,7 @@ class Orientation(Rotation,Crystal):
------- -------
disorientation : Orientation disorientation : Orientation
Disorientation between self and other. Disorientation between self and other.
operators : numpy.ndarray int of shape (...,2), conditional operators : numpy.ndarray of int, shape (...,2), conditional
Index of symmetrically equivalent orientation that rotated vector to the SST. Index of symmetrically equivalent orientation that rotated vector to the SST.
Notes Notes
@ -583,14 +583,14 @@ class Orientation(Rotation,Crystal):
def average(self, def average(self,
weights = None, weights: FloatSequence = None,
return_cloud = False): return_cloud: bool = False):
""" """
Return orientation average over last dimension. Return orientation average over last dimension.
Parameters Parameters
---------- ----------
weights : numpy.ndarray, optional weights : numpy.ndarray, shape (self.shape), optional
Relative weights of orientations. Relative weights of orientations.
return_cloud : bool, optional return_cloud : bool, optional
Return the set of symmetrically equivalent orientations that was used in averaging. Return the set of symmetrically equivalent orientations that was used in averaging.
@ -610,8 +610,8 @@ class Orientation(Rotation,Crystal):
""" """
eq = self.equivalent eq = self.equivalent
m = eq.misorientation(self[...,0].reshape((1,)+self.shape[:-1]+(1,)) m = eq.misorientation(self[...,0].reshape((1,)+self.shape[:-1]+(1,)) # type: ignore
.broadcast_to(eq.shape)).as_axis_angle()[...,3] .broadcast_to(eq.shape)).as_axis_angle()[...,3] # type: ignore
r = Rotation(np.squeeze(np.take_along_axis(eq.quaternion, r = Rotation(np.squeeze(np.take_along_axis(eq.quaternion,
np.argmin(m,axis=0)[np.newaxis,...,np.newaxis], np.argmin(m,axis=0)[np.newaxis,...,np.newaxis],
axis=0), axis=0),
@ -625,7 +625,7 @@ class Orientation(Rotation,Crystal):
def to_SST(self, def to_SST(self,
vector: np.ndarray, vector: FloatSequence,
proper: bool = False, proper: bool = False,
return_operators: bool = False) -> np.ndarray: return_operators: bool = False) -> np.ndarray:
""" """
@ -633,10 +633,10 @@ class Orientation(Rotation,Crystal):
Parameters Parameters
---------- ----------
vector : numpy.ndarray of shape (...,3) vector : numpy.ndarray, shape (...,3)
Lab frame vector to align with crystal frame direction. Lab frame vector to align with crystal frame direction.
Shape of vector blends with shape of own rotation array. Shape of vector blends with shape of own rotation array.
For example, a rotation array of shape (3,2) and a (2,4) vector array result in (3,2,4) outputs. For example, a rotation array of shape (3,2) and a vector array of shape (2,4) result in (3,2,4) outputs.
proper : bool, optional proper : bool, optional
Consider only vectors with z >= 0, hence combine two neighboring SSTs. Consider only vectors with z >= 0, hence combine two neighboring SSTs.
Defaults to False. Defaults to False.
@ -646,15 +646,18 @@ class Orientation(Rotation,Crystal):
Returns Returns
------- -------
vector_SST : numpy.ndarray of shape (...,3) vector_SST : numpy.ndarray, shape (...,3)
Rotated vector falling into SST. Rotated vector falling into SST.
operators : numpy.ndarray int of shape (...), conditional operators : numpy.ndarray of int, shape (...), conditional
Index of symmetrically equivalent orientation that rotated vector to SST. Index of symmetrically equivalent orientation that rotated vector to SST.
""" """
vector_ = np.array(vector,float)
if vector_.shape[-1] != 3:
raise ValueError('input is not a field of three-dimensional vectors')
eq = self.equivalent eq = self.equivalent
blend = util.shapeblender(eq.shape,np.array(vector).shape[:-1]) blend = util.shapeblender(eq.shape,vector_.shape[:-1])
poles = eq.broadcast_to(blend,mode='right') @ np.broadcast_to(np.array(vector),blend+(3,)) #type: ignore poles = eq.broadcast_to(blend,mode='right') @ np.broadcast_to(vector_,blend+(3,)) #type: ignore
ok = self.in_SST(poles,proper=proper) ok = self.in_SST(poles,proper=proper)
ok &= np.cumsum(ok,axis=0) == 1 ok &= np.cumsum(ok,axis=0) == 1
loc = np.where(ok) loc = np.where(ok)
@ -667,14 +670,14 @@ class Orientation(Rotation,Crystal):
def in_SST(self, def in_SST(self,
vector: np.ndarray, vector: FloatSequence,
proper: bool = False) -> Union[np.bool_, np.ndarray]: proper: bool = False) -> Union[np.bool_, np.ndarray]:
""" """
Check whether given crystal frame vector falls into standard stereographic triangle of own symmetry. Check whether given crystal frame vector falls into standard stereographic triangle of own symmetry.
Parameters Parameters
---------- ----------
vector : numpy.ndarray of shape (...,3) vector : numpy.ndarray, shape (...,3)
Vector to check. Vector to check.
proper : bool, optional proper : bool, optional
Consider only vectors with z >= 0, hence combine two neighboring SSTs. Consider only vectors with z >= 0, hence combine two neighboring SSTs.
@ -686,31 +689,32 @@ class Orientation(Rotation,Crystal):
Whether vector falls into SST. Whether vector falls into SST.
""" """
if not isinstance(vector,np.ndarray) or vector.shape[-1] != 3: vector_ = np.array(vector,float)
if vector_.shape[-1] != 3:
raise ValueError('input is not a field of three-dimensional vectors') raise ValueError('input is not a field of three-dimensional vectors')
if self.standard_triangle is None: # direct exit for no symmetry if self.standard_triangle is None: # direct exit for no symmetry
return np.ones_like(vector[...,0],bool) return np.ones_like(vector_[...,0],bool)
if proper: if proper:
components_proper = np.around(np.einsum('...ji,...i', components_proper = np.around(np.einsum('...ji,...i',
np.broadcast_to(self.standard_triangle['proper'], vector.shape+(3,)), np.broadcast_to(self.standard_triangle['proper'], vector_.shape+(3,)),
vector), 12) vector_), 12)
components_improper = np.around(np.einsum('...ji,...i', components_improper = np.around(np.einsum('...ji,...i',
np.broadcast_to(self.standard_triangle['improper'], vector.shape+(3,)), np.broadcast_to(self.standard_triangle['improper'], vector_.shape+(3,)),
vector), 12) vector_), 12)
return np.all(components_proper >= 0.0,axis=-1) \ return np.all(components_proper >= 0.0,axis=-1) \
| np.all(components_improper >= 0.0,axis=-1) | np.all(components_improper >= 0.0,axis=-1)
else: else:
components = np.around(np.einsum('...ji,...i', components = np.around(np.einsum('...ji,...i',
np.broadcast_to(self.standard_triangle['improper'], vector.shape+(3,)), np.broadcast_to(self.standard_triangle['improper'], vector_.shape+(3,)),
np.block([vector[...,:2],np.abs(vector[...,2:3])])), 12) np.block([vector_[...,:2],np.abs(vector_[...,2:3])])), 12)
return np.all(components >= 0.0,axis=-1) return np.all(components >= 0.0,axis=-1)
def IPF_color(self, def IPF_color(self,
vector: np.ndarray, vector: FloatSequence,
in_SST: bool = True, in_SST: bool = True,
proper: bool = False) -> np.ndarray: proper: bool = False) -> np.ndarray:
""" """
@ -718,10 +722,10 @@ class Orientation(Rotation,Crystal):
Parameters Parameters
---------- ----------
vector : numpy.ndarray of shape (...,3) vector : numpy.ndarray, shape (...,3)
Vector to colorize. Vector to colorize.
Shape of vector blends with shape of own rotation array. Shape of vector blends with shape of own rotation array.
For example, a rotation array of shape (3,2) and a (2,4) vector array result in (3,2,4) outputs. For example, a rotation array of shape (3,2) and a vector array of shape (2,4) result in (3,2,4) outputs.
in_SST : bool, optional in_SST : bool, optional
Consider symmetrically equivalent orientations such that poles are located in SST. Consider symmetrically equivalent orientations such that poles are located in SST.
Defaults to True. Defaults to True.
@ -731,7 +735,7 @@ class Orientation(Rotation,Crystal):
Returns Returns
------- -------
rgb : numpy.ndarray of shape (...,3) rgb : numpy.ndarray, shape (...,3)
RGB array of IPF colors. RGB array of IPF colors.
Examples Examples
@ -755,7 +759,7 @@ class Orientation(Rotation,Crystal):
if proper: if proper:
components_proper = np.around(np.einsum('...ji,...i', components_proper = np.around(np.einsum('...ji,...i',
np.broadcast_to(self.standard_triangle['proper'], vector_.shape+(3,)), np.broadcast_to(self.standard_triangle['proper'], vector_.shape+(3,)),
vector_), 12) vector_), 12)
components_improper = np.around(np.einsum('...ji,...i', components_improper = np.around(np.einsum('...ji,...i',
np.broadcast_to(self.standard_triangle['improper'], vector_.shape+(3,)), np.broadcast_to(self.standard_triangle['improper'], vector_.shape+(3,)),
@ -862,16 +866,16 @@ class Orientation(Rotation,Crystal):
Parameters Parameters
---------- ----------
uvw|hkl : numpy.ndarray of shape (...,3) uvw|hkl : numpy.ndarray, shape (...,3)
Miller indices of crystallographic direction or plane normal. Miller indices of crystallographic direction or plane normal.
Shape of vector blends with shape of own rotation array. Shape of vector blends with shape of own rotation array.
For example, a rotation array of shape (3,2) and a (2,4) vector array result in (3,2,4) outputs. For example, a rotation array, shape (3,2) and a vector array of shape (2,4) result in (3,2,4) outputs.
with_symmetry : bool, optional with_symmetry : bool, optional
Calculate all N symmetrically equivalent vectors. Calculate all N symmetrically equivalent vectors.
Returns Returns
------- -------
vector : numpy.ndarray of shape (...,3) or (...,N,3) vector : numpy.ndarray, shape (...,3) or (...,N,3)
Lab frame vector (or vectors if with_symmetry) along [uvw] direction or (hkl) plane normal. Lab frame vector (or vectors if with_symmetry) along [uvw] direction or (hkl) plane normal.
""" """
@ -894,13 +898,13 @@ class Orientation(Rotation,Crystal):
Parameters Parameters
---------- ----------
N_slip|N_twin : iterable of int N_slip|N_twin : '*' or iterable of int
Number of deformation systems per family of the deformation system. Number of deformation systems per family of the deformation system.
Use '*' to select all. Use '*' to select all.
Returns Returns
------- -------
P : numpy.ndarray of shape (N,...,3,3) P : numpy.ndarray, shape (N,...,3,3)
Schmid matrix for each of the N deformation systems. Schmid matrix for each of the N deformation systems.
Examples Examples

View File

@ -108,12 +108,12 @@ class Rotation:
def __getitem__(self, def __getitem__(self,
item: Union[Tuple[int], int, bool, np.bool_, np.ndarray]): item: Union[Tuple[int], int, bool, np.bool_, np.ndarray]):
"""Return slice according to item.""" """Return slice according to item."""
return self.copy() \ return self.copy() if self.shape == () else \
if self.shape == () else \
self.copy(rotation=self.quaternion[item+(slice(None),)] if isinstance(item,tuple) else self.quaternion[item]) self.copy(rotation=self.quaternion[item+(slice(None),)] if isinstance(item,tuple) else self.quaternion[item])
def __eq__(self, other: object) -> bool: def __eq__(self,
other: object) -> bool:
""" """
Equal to other. Equal to other.
@ -123,9 +123,8 @@ class Rotation:
Rotation to check for equality. Rotation to check for equality.
""" """
if not isinstance(other, Rotation): return NotImplemented if not isinstance(other, Rotation) else \
return NotImplemented np.logical_or(np.all(self.quaternion == other.quaternion,axis=-1),
return np.logical_or(np.all(self.quaternion == other.quaternion,axis=-1),
np.all(self.quaternion == -1.0*other.quaternion,axis=-1)) np.all(self.quaternion == -1.0*other.quaternion,axis=-1))
@ -163,7 +162,7 @@ class Rotation:
Returns Returns
------- -------
mask : numpy.ndarray bool mask : numpy.ndarray of bool
Mask indicating where corresponding rotations are close. Mask indicating where corresponding rotations are close.
""" """
@ -228,13 +227,13 @@ class Rotation:
def __pow__(self: MyType, def __pow__(self: MyType,
exp: int) -> MyType: exp: Union[float, int]) -> MyType:
""" """
Perform the rotation 'exp' times. Perform the rotation 'exp' times.
Parameters Parameters
---------- ----------
exp : float exp : scalar
Exponent. Exponent.
""" """
@ -243,13 +242,13 @@ class Rotation:
return self.copy(rotation=Rotation(np.block([np.cos(exp*phi),np.sin(exp*phi)*p]))._standardize()) return self.copy(rotation=Rotation(np.block([np.cos(exp*phi),np.sin(exp*phi)*p]))._standardize())
def __ipow__(self: MyType, def __ipow__(self: MyType,
exp: int) -> MyType: exp: Union[float, int]) -> MyType:
""" """
Perform the rotation 'exp' times (in-place). Perform the rotation 'exp' times (in-place).
Parameters Parameters
---------- ----------
exp : float exp : scalar
Exponent. Exponent.
""" """
@ -263,7 +262,7 @@ class Rotation:
Parameters Parameters
---------- ----------
other : Rotation of shape (self.shape) other : Rotation, shape (self.shape)
Rotation for composition. Rotation for composition.
Returns Returns
@ -290,7 +289,7 @@ class Rotation:
Parameters Parameters
---------- ----------
other : Rotation of shape (self.shape) other : Rotation, shape (self.shape)
Rotation for composition. Rotation for composition.
""" """
@ -304,7 +303,7 @@ class Rotation:
Parameters Parameters
---------- ----------
other : damask.Rotation of shape (self.shape) other : damask.Rotation, shape (self.shape)
Rotation to invert for composition. Rotation to invert for composition.
Returns Returns
@ -325,7 +324,7 @@ class Rotation:
Parameters Parameters
---------- ----------
other : Rotation of shape (self.shape) other : Rotation, shape (self.shape)
Rotation to invert for composition. Rotation to invert for composition.
""" """
@ -339,12 +338,12 @@ class Rotation:
Parameters Parameters
---------- ----------
other : numpy.ndarray of shape (...,3), (...,3,3), or (...,3,3,3,3) other : numpy.ndarray, shape (...,3), (...,3,3), or (...,3,3,3,3)
Vector or tensor on which to apply the rotation. Vector or tensor on which to apply the rotation.
Returns Returns
------- -------
rotated : numpy.ndarray of shape (...,3), (...,3,3), or (...,3,3,3,3) rotated : numpy.ndarray, shape (...,3), (...,3,3), or (...,3,3,3,3)
Rotated vector or tensor, i.e. transformed to frame defined by rotation. Rotated vector or tensor, i.e. transformed to frame defined by rotation.
""" """
@ -401,6 +400,15 @@ class Rotation:
""" """
Flatten array. Flatten array.
Parameters
----------
order : {'C', 'F', 'A'}, optional
'C' flattens in row-major (C-style) order.
'F' flattens in column-major (Fortran-style) order.
'A' flattens in column-major order if object is Fortran contiguous in memory,
row-major order otherwise.
Defaults to 'C'.
Returns Returns
------- -------
flattened : damask.Rotation flattened : damask.Rotation
@ -416,6 +424,18 @@ class Rotation:
""" """
Reshape array. Reshape array.
Parameters
----------
shape : int or tuple of ints
The new shape should be compatible with the original shape.
If an integer is supplied, then the result will be a 1-D array of that length.
order : {'C', 'F', 'A'}, optional
'C' flattens in row-major (C-style) order.
'F' flattens in column-major (Fortran-style) order.
'A' flattens in column-major order if object is Fortran contiguous in memory,
row-major order otherwise.
Defaults to 'C'.
Returns Returns
------- -------
reshaped : damask.Rotation reshaped : damask.Rotation
@ -434,7 +454,7 @@ class Rotation:
Parameters Parameters
---------- ----------
shape : int, tuple shape : int or tuple of ints
Shape of broadcasted array. Shape of broadcasted array.
mode : str, optional mode : str, optional
Where to preferentially locate missing dimensions. Where to preferentially locate missing dimensions.
@ -458,7 +478,7 @@ class Rotation:
Parameters Parameters
---------- ----------
weights : FloatSequence, optional weights : numpy.ndarray, optional
Relative weight of each rotation. Relative weight of each rotation.
Returns Returns
@ -518,7 +538,7 @@ class Rotation:
Returns Returns
------- -------
q : numpy.ndarray of shape (...,4) q : numpy.ndarray, shape (...,4)
Unit quaternion (q_0, q_1, q_2, q_3) in positive real hemisphere, i.e. ǀqǀ = 1, q_0 0. Unit quaternion (q_0, q_1, q_2, q_3) in positive real hemisphere, i.e. ǀqǀ = 1, q_0 0.
""" """
@ -536,7 +556,7 @@ class Rotation:
Returns Returns
------- -------
phi : numpy.ndarray of shape (...,3) phi : numpy.ndarray, shape (...,3)
Bunge Euler angles (φ_1 [0,2π], ϕ [0,π], φ_2 [0,2π]) Bunge Euler angles (φ_1 [0,2π], ϕ [0,π], φ_2 [0,2π])
or (φ_1 [0,360], ϕ [0,180], φ_2 [0,360]) if degrees == True. or (φ_1 [0,360], ϕ [0,180], φ_2 [0,360]) if degrees == True.
@ -555,8 +575,7 @@ class Rotation:
""" """
eu = Rotation._qu2eu(self.quaternion) eu = Rotation._qu2eu(self.quaternion)
if degrees: eu = np.degrees(eu) return np.degrees(eu) if degrees else eu
return eu
def as_axis_angle(self, def as_axis_angle(self,
degrees: bool = False, degrees: bool = False,
@ -573,7 +592,7 @@ class Rotation:
Returns Returns
------- -------
axis_angle : numpy.ndarray of shape (...,4) or tuple ((...,3), (...)) if pair == True axis_angle : numpy.ndarray, shape (...,4) or tuple ((...,3), (...)) if pair == True
Axis and angle [n_1, n_2, n_3, ω] with ǀnǀ = 1 and ω [0,π] Axis and angle [n_1, n_2, n_3, ω] with ǀnǀ = 1 and ω [0,π]
or ω [0,180] if degrees == True. or ω [0,180] if degrees == True.
@ -597,7 +616,7 @@ class Rotation:
Returns Returns
------- -------
R : numpy.ndarray of shape (...,3,3) R : numpy.ndarray, shape (...,3,3)
Rotation matrix R with det(R) = 1, R.T R = I. Rotation matrix R with det(R) = 1, R.T R = I.
Examples Examples
@ -627,7 +646,7 @@ class Rotation:
Returns Returns
------- -------
rho : numpy.ndarray of shape (...,4) or (...,3) if compact == True rho : numpy.ndarray, shape (...,4) or (...,3) if compact == True
RodriguesFrank vector [n_1, n_2, n_3, tan(ω/2)] with ǀnǀ = 1 and ω [0,π] RodriguesFrank vector [n_1, n_2, n_3, tan(ω/2)] with ǀnǀ = 1 and ω [0,π]
or [n_1, n_2, n_3] with ǀnǀ = tan(ω/2) and ω [0,π] if compact == True. or [n_1, n_2, n_3] with ǀnǀ = tan(ω/2) and ω [0,π] if compact == True.
@ -654,7 +673,7 @@ class Rotation:
Returns Returns
------- -------
h : numpy.ndarray of shape (...,3) h : numpy.ndarray, shape (...,3)
Homochoric vector (h_1, h_2, h_3) with ǀhǀ < (3/4*π)^(1/3). Homochoric vector (h_1, h_2, h_3) with ǀhǀ < (3/4*π)^(1/3).
Examples Examples
@ -675,7 +694,7 @@ class Rotation:
Returns Returns
------- -------
x : numpy.ndarray of shape (...,3) x : numpy.ndarray, shape (...,3)
Cubochoric vector (x_1, x_2, x_3) with max(x_i) < 1/2*π^(2/3). Cubochoric vector (x_1, x_2, x_3) with max(x_i) < 1/2*π^(2/3).
Examples Examples
@ -702,7 +721,7 @@ class Rotation:
Parameters Parameters
---------- ----------
q : numpy.ndarray of shape (...,4) q : numpy.ndarray, shape (...,4)
Unit quaternion (q_0, q_1, q_2, q_3) in positive real hemisphere, i.e. ǀqǀ = 1, q_0 0. Unit quaternion (q_0, q_1, q_2, q_3) in positive real hemisphere, i.e. ǀqǀ = 1, q_0 0.
accept_homomorph : bool, optional accept_homomorph : bool, optional
Allow homomorphic variants, i.e. q_0 < 0 (negative real hemisphere). Allow homomorphic variants, i.e. q_0 < 0 (negative real hemisphere).
@ -736,7 +755,7 @@ class Rotation:
Parameters Parameters
---------- ----------
phi : numpy.ndarray of shape (...,3) phi : numpy.ndarray, shape (...,3)
Euler angles (φ_1 [0,2π], ϕ [0,π], φ_2 [0,2π]) Euler angles (φ_1 [0,2π], ϕ [0,π], φ_2 [0,2π])
or (φ_1 [0,360], ϕ [0,180], φ_2 [0,360]) if degrees == True. or (φ_1 [0,360], ϕ [0,180], φ_2 [0,360]) if degrees == True.
degrees : bool, optional degrees : bool, optional
@ -767,7 +786,7 @@ class Rotation:
Parameters Parameters
---------- ----------
axis_angle : numpy.ndarray of shape (...,4) axis_angle : numpy.ndarray, shape (...,4)
Axis and angle (n_1, n_2, n_3, ω) with ǀnǀ = 1 and ω [0,π] Axis and angle (n_1, n_2, n_3, ω) with ǀnǀ = 1 and ω [0,π]
or ω [0,180] if degrees == True. or ω [0,180] if degrees == True.
degrees : bool, optional degrees : bool, optional
@ -804,7 +823,7 @@ class Rotation:
Parameters Parameters
---------- ----------
basis : numpy.ndarray of shape (...,3,3) basis : numpy.ndarray, shape (...,3,3)
Three three-dimensional lattice basis vectors. Three three-dimensional lattice basis vectors.
orthonormal : bool, optional orthonormal : bool, optional
Basis is strictly orthonormal, i.e. is free of stretch components. Defaults to True. Basis is strictly orthonormal, i.e. is free of stretch components. Defaults to True.
@ -838,7 +857,7 @@ class Rotation:
Parameters Parameters
---------- ----------
R : numpy.ndarray of shape (...,3,3) R : numpy.ndarray, shape (...,3,3)
Rotation matrix with det(R) = 1, R.T R = I. Rotation matrix with det(R) = 1, R.T R = I.
""" """
@ -852,9 +871,9 @@ class Rotation:
Parameters Parameters
---------- ----------
a : numpy.ndarray of shape (...,2,3) a : numpy.ndarray, shape (...,2,3)
Two three-dimensional lattice vectors of first orthogonal basis. Two three-dimensional lattice vectors of first orthogonal basis.
b : numpy.ndarray of shape (...,2,3) b : numpy.ndarray, shape (...,2,3)
Corresponding three-dimensional lattice vectors of second basis. Corresponding three-dimensional lattice vectors of second basis.
""" """
@ -882,7 +901,7 @@ class Rotation:
Parameters Parameters
---------- ----------
rho : numpy.ndarray of shape (...,4) rho : numpy.ndarray, shape (...,4)
RodriguesFrank vector (n_1, n_2, n_3, tan(ω/2)) with ǀnǀ = 1 and ω [0,π]. RodriguesFrank vector (n_1, n_2, n_3, tan(ω/2)) with ǀnǀ = 1 and ω [0,π].
normalize : bool, optional normalize : bool, optional
Allow ǀnǀ 1. Defaults to False. Allow ǀnǀ 1. Defaults to False.
@ -913,7 +932,7 @@ class Rotation:
Parameters Parameters
---------- ----------
h : numpy.ndarray of shape (...,3) h : numpy.ndarray, shape (...,3)
Homochoric vector (h_1, h_2, h_3) with ǀhǀ < (3/4*π)^(1/3). Homochoric vector (h_1, h_2, h_3) with ǀhǀ < (3/4*π)^(1/3).
P : int {-1,1}, optional P : int {-1,1}, optional
Sign convention. Defaults to -1. Sign convention. Defaults to -1.
@ -940,7 +959,7 @@ class Rotation:
Parameters Parameters
---------- ----------
x : numpy.ndarray of shape (...,3) x : numpy.ndarray, shape (...,3)
Cubochoric vector (x_1, x_2, x_3) with max(x_i) < 1/2*π^(2/3). Cubochoric vector (x_1, x_2, x_3) with max(x_i) < 1/2*π^(2/3).
P : int {-1,1}, optional P : int {-1,1}, optional
Sign convention. Defaults to -1. Sign convention. Defaults to -1.
@ -1002,9 +1021,9 @@ class Rotation:
Parameters Parameters
---------- ----------
weights : numpy.ndarray of shape (n) weights : numpy.ndarray, shape (n)
Texture intensity values (probability density or volume fraction) at Euler space grid points. Texture intensity values (probability density or volume fraction) at Euler space grid points.
phi : numpy.ndarray of shape (n,3) phi : numpy.ndarray, shape (n,3)
Grid coordinates in Euler space at which weights are defined. Grid coordinates in Euler space at which weights are defined.
N : integer, optional N : integer, optional
Number of discrete orientations to be sampled from the given ODF. Number of discrete orientations to be sampled from the given ODF.
@ -1020,14 +1039,14 @@ class Rotation:
Returns Returns
------- -------
samples : damask.Rotation of shape (N) samples : damask.Rotation, shape (N)
Array of sampled rotations closely representing the input ODF. Array of sampled rotations that approximate the input ODF.
Notes Notes
----- -----
Due to the distortion of Euler space in the vicinity of ϕ = 0, probability densities, p, defined on Due to the distortion of Euler space in the vicinity of ϕ = 0, probability densities, p, defined on
grid points with ϕ = 0 will never result in reconstructed orientations as their dV/V = p dγ = p × 0. grid points with ϕ = 0 will never result in reconstructed orientations as their dV/V = p dγ = p × 0.
Hence, it is recommended to transform any such dataset to cell centers that avoid grid points at ϕ = 0. Hence, it is recommended to transform any such dataset to a cell-centered version, which avoids grid points at ϕ = 0.
References References
---------- ----------
@ -1095,9 +1114,9 @@ class Rotation:
Parameters Parameters
---------- ----------
alpha : numpy.ndarray of shape (2) alpha : numpy.ndarray, shape (2)
Polar coordinates (phi from x, theta from z) of fiber direction in crystal frame. Polar coordinates (phi from x, theta from z) of fiber direction in crystal frame.
beta : numpy.ndarray of shape (2) beta : numpy.ndarray, shape (2)
Polar coordinates (phi from x, theta from z) of fiber direction in sample frame. Polar coordinates (phi from x, theta from z) of fiber direction in sample frame.
sigma : float, optional sigma : float, optional
Standard deviation of (Gaussian) misorientation distribution. Standard deviation of (Gaussian) misorientation distribution.

View File

@ -185,7 +185,7 @@ class Table:
Returns Returns
------- -------
mask : numpy.ndarray bool mask : numpy.ndarray of bool
Mask indicating where corresponding table values are close. Mask indicating where corresponding table values are close.
""" """