Merge branch '241-same-quaternion-but-different-ipfs-calculated-by-damask-and-orix-python-package' into 'development'

lower exponent gives smoother gradients

Closes #241

See merge request damask/DAMASK!675
This commit is contained in:
Sharan Roongta 2022-12-06 14:15:35 +00:00
commit 4717536f1a
1 changed files with 21 additions and 8 deletions

View File

@ -627,7 +627,7 @@ class Orientation(Rotation,Crystal):
weights : numpy.ndarray, shape (self.shape), optional
Relative weights of orientations.
return_cloud : bool, optional
Return the set of symmetrically equivalent orientations that was used in averaging.
Return the specific (symmetrically equivalent) orientations that were averaged.
Defaults to False.
Returns
@ -635,7 +635,7 @@ class Orientation(Rotation,Crystal):
average : Orientation
Weighted average of original Orientation field.
cloud : Orientations, conditional
Set of symmetrically equivalent orientations that were used in averaging.
Symmetrically equivalent version of each orientation that were actually used in averaging.
References
----------
@ -660,7 +660,7 @@ class Orientation(Rotation,Crystal):
proper: bool = False,
return_operators: bool = False) -> np.ndarray:
"""
Rotate vector to ensure it falls into (improper or proper) standard stereographic triangle of crystal symmetry.
Rotate lab frame vector to ensure it falls into (improper or proper) standard stereographic triangle of crystal symmetry.
Parameters
----------
@ -679,7 +679,7 @@ class Orientation(Rotation,Crystal):
-------
vector_SST : numpy.ndarray, shape (...,3)
Rotated vector falling into SST.
operators : numpy.ndarray of int, shape (...), conditional
operator : numpy.ndarray of int, shape (...), conditional
Index of symmetrically equivalent orientation that rotated vector to SST.
"""
@ -749,12 +749,12 @@ class Orientation(Rotation,Crystal):
in_SST: bool = True,
proper: bool = False) -> np.ndarray:
"""
Map vector to RGB color within standard stereographic triangle of own symmetry.
Map lab frame vector to RGB color within standard stereographic triangle of own symmetry.
Parameters
----------
vector : numpy.ndarray, shape (...,3)
Vector to colorize.
Lab frame vector to colorize.
Shape of vector blends with shape of own rotation array.
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
@ -771,13 +771,26 @@ class Orientation(Rotation,Crystal):
Examples
--------
Inverse pole figure color of the e_3 direction for a crystal in "Cube" orientation with cubic symmetry:
Inverse pole figure color of the e_3 lab direction for a
crystal in "Cube" orientation with cubic symmetry:
>>> import damask
>>> o = damask.Orientation(family='cubic')
>>> o.IPF_color([0,0,1])
array([1., 0., 0.])
Sample standard triangle for hexagonal symmetry:
>>> import damask
>>> from matplotlib import pyplot as plt
>>> lab = [0,0,1]
>>> o = damask.Orientation.from_random(shape=500000,family='hexagonal')
>>> coord = damask.util.project_equal_area(o.to_SST(lab))
>>> color = o.IPF_color(lab)
>>> plt.scatter(coord[:,0],coord[:,1],color=color,s=.06)
>>> plt.axis('scaled')
>>> plt.show()
"""
if np.array(vector).shape[-1] != 3:
raise ValueError('input is not a field of three-dimensional vectors')
@ -807,7 +820,7 @@ class Orientation(Rotation,Crystal):
in_SST_ = np.all(components >= 0.0,axis=-1)
with np.errstate(invalid='ignore',divide='ignore'):
rgb = (components/np.linalg.norm(components,axis=-1,keepdims=True))**0.5 # smoothen color ramps
rgb = (components/np.linalg.norm(components,axis=-1,keepdims=True))**(1./3.) # smoothen color ramps
rgb = np.clip(rgb,0.,1.) # clip intensity
rgb /= np.max(rgb,axis=-1,keepdims=True) # normalize to (HS)V = 1
rgb[np.broadcast_to(~in_SST_[...,np.newaxis],rgb.shape)] = 0.0