to_pole now blends; corrected help texts
This commit is contained in:
parent
af6a99921f
commit
b754617c76
|
@ -524,8 +524,8 @@ class Orientation(Rotation,Crystal):
|
||||||
Disorientation between two specific orientations of hexagonal symmetry:
|
Disorientation between two specific orientations of hexagonal symmetry:
|
||||||
|
|
||||||
>>> import damask
|
>>> import damask
|
||||||
>>> a = damask.Orientation.from_Eulers(phi=[123,32,21],degrees=True,lattice='hexagonal')
|
>>> a = damask.Orientation.from_Euler_angles(phi=[123,32,21],degrees=True,family='hexagonal')
|
||||||
>>> b = damask.Orientation.from_Eulers(phi=[104,11,87],degrees=True,lattice='hexagonal')
|
>>> b = damask.Orientation.from_Euler_angles(phi=[104,11,87],degrees=True,family='hexagonal')
|
||||||
>>> a.disorientation(b)
|
>>> a.disorientation(b)
|
||||||
Crystal family hexagonal
|
Crystal family hexagonal
|
||||||
Quaternion: (real=0.976, imag=<+0.189, +0.018, +0.103>)
|
Quaternion: (real=0.976, imag=<+0.189, +0.018, +0.103>)
|
||||||
|
@ -616,7 +616,7 @@ class Orientation(Rotation,Crystal):
|
||||||
----------
|
----------
|
||||||
vector : numpy.ndarray of shape (...,3)
|
vector : numpy.ndarray of shape (...,3)
|
||||||
Lab frame vector to align with crystal frame direction.
|
Lab frame vector to align with crystal frame direction.
|
||||||
Shape of other 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 (2,4) vector array 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.
|
||||||
|
@ -696,6 +696,8 @@ class Orientation(Rotation,Crystal):
|
||||||
----------
|
----------
|
||||||
vector : numpy.ndarray of shape (...,3)
|
vector : numpy.ndarray of shape (...,3)
|
||||||
Vector to colorize.
|
Vector to colorize.
|
||||||
|
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.
|
||||||
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.
|
||||||
|
@ -713,7 +715,7 @@ class Orientation(Rotation,Crystal):
|
||||||
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 direction for a crystal in "Cube" orientation with cubic symmetry:
|
||||||
|
|
||||||
>>> import damask
|
>>> import damask
|
||||||
>>> o = damask.Orientation(lattice='cubic')
|
>>> o = damask.Orientation(family='cubic')
|
||||||
>>> o.IPF_color([0,0,1])
|
>>> o.IPF_color([0,0,1])
|
||||||
array([1., 0., 0.])
|
array([1., 0., 0.])
|
||||||
|
|
||||||
|
@ -835,22 +837,27 @@ class Orientation(Rotation,Crystal):
|
||||||
----------
|
----------
|
||||||
uvw|hkl : numpy.ndarray of shape (...,3)
|
uvw|hkl : numpy.ndarray of 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.
|
||||||
|
For example, a rotation array of shape (3,2) and a (2,4) vector array 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 of 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.
|
||||||
|
|
||||||
"""
|
"""
|
||||||
v = self.to_frame(uvw=uvw,hkl=hkl)
|
v = self.to_frame(uvw=uvw,hkl=hkl)
|
||||||
|
blend = util.shapeblender(self.shape,v.shape[:-1])
|
||||||
if with_symmetry:
|
if with_symmetry:
|
||||||
sym_ops = self.symmetry_operations
|
sym_ops = self.symmetry_operations
|
||||||
v = sym_ops.broadcast_to(sym_ops.shape+v.shape[:-1],mode='right') \
|
shape = v.shape[:-1]+sym_ops.shape
|
||||||
@ np.broadcast_to(v,sym_ops.shape+v.shape)
|
blend += sym_ops.shape
|
||||||
return ~(self if self.shape+v.shape[:-1] == () else self.broadcast_to(self.shape+v.shape[:-1],mode='right')) \
|
v = sym_ops.broadcast_to(shape) \
|
||||||
@ np.broadcast_to(v,self.shape+v.shape)
|
@ np.broadcast_to(v.reshape(util.shapeshifter(v.shape,shape+(3,))),shape+(3,))
|
||||||
|
return ~(self.broadcast_to(blend)) \
|
||||||
|
@ np.broadcast_to(v,blend+(3,))
|
||||||
|
|
||||||
|
|
||||||
def Schmid(self,*,N_slip=None,N_twin=None):
|
def Schmid(self,*,N_slip=None,N_twin=None):
|
||||||
|
|
|
@ -77,6 +77,16 @@ class TestOrientation:
|
||||||
with pytest.raises(ValueError):
|
with pytest.raises(ValueError):
|
||||||
Orientation(**kwargs)
|
Orientation(**kwargs)
|
||||||
|
|
||||||
|
@pytest.mark.parametrize('invalid_family',[None,'fcc','bcc','hello'])
|
||||||
|
def test_invalid_family_init(self,invalid_family):
|
||||||
|
with pytest.raises(KeyError):
|
||||||
|
Orientation(family=invalid_family)
|
||||||
|
|
||||||
|
@pytest.mark.parametrize('invalid_lattice',[None,'fcc','bcc','hello'])
|
||||||
|
def test_invalid_family_init(self,invalid_lattice):
|
||||||
|
with pytest.raises(KeyError):
|
||||||
|
Orientation(lattice=invalid_lattice)
|
||||||
|
|
||||||
@pytest.mark.parametrize('kwargs',[
|
@pytest.mark.parametrize('kwargs',[
|
||||||
dict(lattice='aP',a=1.0,b=1.1,c=1.2,alpha=np.pi/4,beta=np.pi/3,gamma=np.pi/2),
|
dict(lattice='aP',a=1.0,b=1.1,c=1.2,alpha=np.pi/4,beta=np.pi/3,gamma=np.pi/2),
|
||||||
dict(lattice='mP',a=1.0,b=1.1,c=1.2, beta=np.pi/3 ),
|
dict(lattice='mP',a=1.0,b=1.1,c=1.2, beta=np.pi/3 ),
|
||||||
|
@ -203,6 +213,15 @@ class TestOrientation:
|
||||||
FZ = np.argmin(abs(eq.misorientation(i.broadcast_to(len(eq))).as_axis_angle(pair=True)[1]))
|
FZ = np.argmin(abs(eq.misorientation(i.broadcast_to(len(eq))).as_axis_angle(pair=True)[1]))
|
||||||
assert o.reduced == eq[FZ]
|
assert o.reduced == eq[FZ]
|
||||||
|
|
||||||
|
@pytest.mark.parametrize('family',crystal_families)
|
||||||
|
def test_reduced_corner_cases(self,family):
|
||||||
|
# test whether there is always a sym-eq rotation that falls into the FZ
|
||||||
|
N = np.random.randint(10,40)
|
||||||
|
size = np.ones(3)*np.pi**(2./3.)
|
||||||
|
grid = grid_filters.coordinates0_node([N+1,N+1,N+1],size,-size*.5)
|
||||||
|
evenly_distributed = Orientation.from_cubochoric(x=grid[:-2,:-2,:-2],family=family)
|
||||||
|
assert evenly_distributed.shape == evenly_distributed.reduced.shape
|
||||||
|
|
||||||
@pytest.mark.parametrize('family',crystal_families)
|
@pytest.mark.parametrize('family',crystal_families)
|
||||||
@pytest.mark.parametrize('N',[1,8,32])
|
@pytest.mark.parametrize('N',[1,8,32])
|
||||||
def test_disorientation(self,family,N):
|
def test_disorientation(self,family,N):
|
||||||
|
@ -221,81 +240,12 @@ class TestOrientation:
|
||||||
.misorientation(p[n].equivalent[ops[n][1]])
|
.misorientation(p[n].equivalent[ops[n][1]])
|
||||||
.as_quaternion())
|
.as_quaternion())
|
||||||
|
|
||||||
@pytest.mark.parametrize('family',crystal_families)
|
|
||||||
@pytest.mark.parametrize('a,b',[
|
|
||||||
((2,3,2),(2,3,2)),
|
|
||||||
((2,2),(4,4)),
|
|
||||||
((3,1),(1,3)),
|
|
||||||
(None,None),
|
|
||||||
])
|
|
||||||
def test_disorientation_blending(self,family,a,b):
|
|
||||||
o = Orientation.from_random(family=family,shape=a)
|
|
||||||
p = Orientation.from_random(family=family,shape=b)
|
|
||||||
blend = util.shapeblender(o.shape,p.shape)
|
|
||||||
for loc in np.random.randint(0,blend,(10,len(blend))):
|
|
||||||
assert o[tuple(loc[:len(o.shape)])].disorientation(p[tuple(loc[-len(p.shape):])]) \
|
|
||||||
.isclose(o.disorientation(p)[tuple(loc)])
|
|
||||||
|
|
||||||
@pytest.mark.parametrize('family',crystal_families)
|
@pytest.mark.parametrize('family',crystal_families)
|
||||||
def test_disorientation360(self,family):
|
def test_disorientation360(self,family):
|
||||||
o_1 = Orientation(Rotation(),family=family)
|
o_1 = Orientation(Rotation(),family=family)
|
||||||
o_2 = Orientation.from_Euler_angles(family=family,phi=[360,0,0],degrees=True)
|
o_2 = Orientation.from_Euler_angles(family=family,phi=[360,0,0],degrees=True)
|
||||||
assert np.allclose((o_1.disorientation(o_2)).as_matrix(),np.eye(3))
|
assert np.allclose((o_1.disorientation(o_2)).as_matrix(),np.eye(3))
|
||||||
|
|
||||||
@pytest.mark.parametrize('family',crystal_families)
|
|
||||||
@pytest.mark.parametrize('shape',[(1),(2,3),(4,3,2)])
|
|
||||||
def test_reduced_vectorization(self,family,shape):
|
|
||||||
o = Orientation.from_random(family=family,shape=shape)
|
|
||||||
for r, theO in zip(o.reduced.flatten(),o.flatten()):
|
|
||||||
assert r == theO.reduced
|
|
||||||
|
|
||||||
@pytest.mark.parametrize('family',crystal_families)
|
|
||||||
def test_reduced_corner_cases(self,family):
|
|
||||||
# test whether there is always a sym-eq rotation that falls into the FZ
|
|
||||||
N = np.random.randint(10,40)
|
|
||||||
size = np.ones(3)*np.pi**(2./3.)
|
|
||||||
grid = grid_filters.coordinates0_node([N+1,N+1,N+1],size,-size*.5)
|
|
||||||
evenly_distributed = Orientation.from_cubochoric(x=grid[:-2,:-2,:-2],family=family)
|
|
||||||
assert evenly_distributed.shape == evenly_distributed.reduced.shape
|
|
||||||
|
|
||||||
|
|
||||||
@pytest.mark.parametrize('family',crystal_families)
|
|
||||||
@pytest.mark.parametrize('shape',[(1),(2,3),(4,3,2)])
|
|
||||||
@pytest.mark.parametrize('vector',np.array([[1,0,0],[1,2,3],[-1,1,-1]]))
|
|
||||||
@pytest.mark.parametrize('proper',[True,False])
|
|
||||||
def test_to_SST_vectorization(self,family,shape,vector,proper):
|
|
||||||
o = Orientation.from_random(family=family,shape=shape)
|
|
||||||
for r, theO in zip(o.to_SST(vector=vector,proper=proper).reshape((-1,3)),o.flatten()):
|
|
||||||
assert np.allclose(r,theO.to_SST(vector=vector,proper=proper))
|
|
||||||
|
|
||||||
@pytest.mark.parametrize('family',crystal_families)
|
|
||||||
@pytest.mark.parametrize('shape',[(1),(2,3),(4,3,2)])
|
|
||||||
@pytest.mark.parametrize('vector',np.array([[1,0,0],[1,2,3],[-1,1,-1]]))
|
|
||||||
@pytest.mark.parametrize('proper',[True,False])
|
|
||||||
@pytest.mark.parametrize('in_SST',[True,False])
|
|
||||||
def test_IPF_color_vectorization(self,family,shape,vector,proper,in_SST):
|
|
||||||
o = Orientation.from_random(family=family,shape=shape)
|
|
||||||
for r, theO in zip(o.IPF_color(vector,in_SST=in_SST,proper=proper).reshape((-1,3)),o.flatten()):
|
|
||||||
assert np.allclose(r,theO.IPF_color(vector,in_SST=in_SST,proper=proper))
|
|
||||||
|
|
||||||
@pytest.mark.parametrize('family',crystal_families)
|
|
||||||
@pytest.mark.parametrize('a,b',[
|
|
||||||
((2,3,2),(2,3,2)),
|
|
||||||
((2,2),(4,4)),
|
|
||||||
((3,1),(1,3)),
|
|
||||||
(None,(3,)),
|
|
||||||
])
|
|
||||||
def test_to_SST_blending(self,family,a,b):
|
|
||||||
o = Orientation.from_random(family=family,shape=a)
|
|
||||||
v = np.random.random(b+(3,))
|
|
||||||
blend = util.shapeblender(o.shape,b)
|
|
||||||
for loc in np.random.randint(0,blend,(10,len(blend))):
|
|
||||||
print(f'{a}/{b} @ {loc}')
|
|
||||||
print(o[tuple(loc[:len(o.shape)])].to_SST(v[tuple(loc[-len(b):])]))
|
|
||||||
print(o.to_SST(v)[tuple(loc)])
|
|
||||||
assert np.allclose(o[tuple(loc[:len(o.shape)])].to_SST(v[tuple(loc[-len(b):])]),
|
|
||||||
o.to_SST(v)[tuple(loc)])
|
|
||||||
|
|
||||||
@pytest.mark.parametrize('color',[{'label':'red', 'RGB':[1,0,0],'direction':[0,0,1]},
|
@pytest.mark.parametrize('color',[{'label':'red', 'RGB':[1,0,0],'direction':[0,0,1]},
|
||||||
{'label':'green','RGB':[0,1,0],'direction':[0,1,1]},
|
{'label':'green','RGB':[0,1,0],'direction':[0,1,1]},
|
||||||
{'label':'blue', 'RGB':[0,0,1],'direction':[1,1,1]}])
|
{'label':'blue', 'RGB':[0,0,1],'direction':[1,1,1]}])
|
||||||
|
@ -314,33 +264,6 @@ class TestOrientation:
|
||||||
color = o.IPF_color(vector=direction,proper=proper)
|
color = o.IPF_color(vector=direction,proper=proper)
|
||||||
assert np.allclose(np.broadcast_to(color[0,...],color.shape),color)
|
assert np.allclose(np.broadcast_to(color[0,...],color.shape),color)
|
||||||
|
|
||||||
@pytest.mark.parametrize('family',crystal_families)
|
|
||||||
def test_in_FZ_vectorization(self,set_of_rodrigues,family):
|
|
||||||
result = Orientation.from_Rodrigues_vector(rho=set_of_rodrigues.reshape((-1,4,4)),family=family).in_FZ.reshape(-1)
|
|
||||||
for r,rho in zip(result,set_of_rodrigues[:len(result)]):
|
|
||||||
assert r == Orientation.from_Rodrigues_vector(rho=rho,family=family).in_FZ
|
|
||||||
|
|
||||||
@pytest.mark.parametrize('family',crystal_families)
|
|
||||||
def test_in_disorientation_FZ_vectorization(self,set_of_rodrigues,family):
|
|
||||||
result = Orientation.from_Rodrigues_vector(rho=set_of_rodrigues.reshape((-1,4,4)),
|
|
||||||
family=family).in_disorientation_FZ.reshape(-1)
|
|
||||||
for r,rho in zip(result,set_of_rodrigues[:len(result)]):
|
|
||||||
assert r == Orientation.from_Rodrigues_vector(rho=rho,family=family).in_disorientation_FZ
|
|
||||||
|
|
||||||
@pytest.mark.parametrize('proper',[True,False])
|
|
||||||
@pytest.mark.parametrize('family',crystal_families)
|
|
||||||
def test_in_SST_vectorization(self,family,proper):
|
|
||||||
vecs = np.random.rand(20,4,3)
|
|
||||||
result = Orientation(family=family).in_SST(vecs,proper).flatten()
|
|
||||||
for r,v in zip(result,vecs.reshape((-1,3))):
|
|
||||||
assert np.all(r == Orientation(family=family).in_SST(v,proper))
|
|
||||||
|
|
||||||
@pytest.mark.parametrize('invalid_family',[None,'fcc','bcc','hello'])
|
|
||||||
def test_invalid_lattice_init(self,invalid_family):
|
|
||||||
with pytest.raises(KeyError):
|
|
||||||
Orientation(family=invalid_family)
|
|
||||||
|
|
||||||
|
|
||||||
@pytest.mark.parametrize('relation',[None,'Peter','Paul'])
|
@pytest.mark.parametrize('relation',[None,'Peter','Paul'])
|
||||||
def test_unknown_relation(self,relation):
|
def test_unknown_relation(self,relation):
|
||||||
with pytest.raises(KeyError):
|
with pytest.raises(KeyError):
|
||||||
|
@ -371,12 +294,6 @@ class TestOrientation:
|
||||||
o = Orientation(family='cubic') # noqa
|
o = Orientation(family='cubic') # noqa
|
||||||
with pytest.raises(ValueError):
|
with pytest.raises(ValueError):
|
||||||
eval(f'o.{function}(np.ones(4))')
|
eval(f'o.{function}(np.ones(4))')
|
||||||
@pytest.mark.parametrize('model',['Bain','KS','GT','GT_prime','NW','Pitsch'])
|
|
||||||
@pytest.mark.parametrize('lattice',['cF','cI'])
|
|
||||||
def test_relationship_vectorize(self,set_of_quaternions,lattice,model):
|
|
||||||
r = Orientation(rotation=set_of_quaternions[:200].reshape((50,4,4)),lattice=lattice).related(model)
|
|
||||||
for i in range(200):
|
|
||||||
assert (r.reshape((-1,200))[:,i] == Orientation(set_of_quaternions[i],lattice=lattice).related(model)).all()
|
|
||||||
|
|
||||||
@pytest.mark.parametrize('model',['Bain','KS','GT','GT_prime','NW','Pitsch'])
|
@pytest.mark.parametrize('model',['Bain','KS','GT','GT_prime','NW','Pitsch'])
|
||||||
@pytest.mark.parametrize('lattice',['cF','cI'])
|
@pytest.mark.parametrize('lattice',['cF','cI'])
|
||||||
|
@ -411,7 +328,6 @@ class TestOrientation:
|
||||||
)
|
)
|
||||||
assert np.allclose(o.to_frame(uvw=np.eye(3)),basis), 'Lattice basis disagrees with initialization'
|
assert np.allclose(o.to_frame(uvw=np.eye(3)),basis), 'Lattice basis disagrees with initialization'
|
||||||
|
|
||||||
|
|
||||||
@pytest.mark.parametrize('lattice,a,b,c,alpha,beta,gamma',
|
@pytest.mark.parametrize('lattice,a,b,c,alpha,beta,gamma',
|
||||||
[
|
[
|
||||||
('aP',0.5,2.0,3.0,0.8,0.5,1.2),
|
('aP',0.5,2.0,3.0,0.8,0.5,1.2),
|
||||||
|
@ -421,7 +337,6 @@ class TestOrientation:
|
||||||
('hP',1.0,1.0,1.6,np.pi/2,np.pi/2,2*np.pi/3),
|
('hP',1.0,1.0,1.6,np.pi/2,np.pi/2,2*np.pi/3),
|
||||||
('cF',1.0,1.0,1.0,np.pi/2,np.pi/2,np.pi/2),
|
('cF',1.0,1.0,1.0,np.pi/2,np.pi/2,np.pi/2),
|
||||||
])
|
])
|
||||||
|
|
||||||
@pytest.mark.parametrize('kw',['uvw','hkl'])
|
@pytest.mark.parametrize('kw',['uvw','hkl'])
|
||||||
@pytest.mark.parametrize('with_symmetry',[False,True])
|
@pytest.mark.parametrize('with_symmetry',[False,True])
|
||||||
@pytest.mark.parametrize('shape',[None,1,(12,24)])
|
@pytest.mark.parametrize('shape',[None,1,(12,24)])
|
||||||
|
@ -436,7 +351,7 @@ class TestOrientation:
|
||||||
a=a,b=b,c=c,
|
a=a,b=b,c=c,
|
||||||
alpha=alpha,beta=beta,gamma=gamma)
|
alpha=alpha,beta=beta,gamma=gamma)
|
||||||
assert o.to_pole(**{kw:vector,'with_symmetry':with_symmetry}).shape \
|
assert o.to_pole(**{kw:vector,'with_symmetry':with_symmetry}).shape \
|
||||||
== o.shape + (o.symmetry_operations.shape if with_symmetry else ()) + vector.shape
|
== o.shape + vector.shape[:-1] + (o.symmetry_operations.shape if with_symmetry else ()) + vector.shape[-1:]
|
||||||
|
|
||||||
@pytest.mark.parametrize('lattice',['hP','cI','cF']) #tI not included yet
|
@pytest.mark.parametrize('lattice',['hP','cI','cF']) #tI not included yet
|
||||||
def test_Schmid(self,update,ref_path,lattice):
|
def test_Schmid(self,update,ref_path,lattice):
|
||||||
|
@ -449,11 +364,144 @@ class TestOrientation:
|
||||||
table.save(reference)
|
table.save(reference)
|
||||||
assert np.allclose(P,Table.load(reference).get('Schmid'))
|
assert np.allclose(P,Table.load(reference).get('Schmid'))
|
||||||
|
|
||||||
|
### vectorization tests ###
|
||||||
|
|
||||||
@pytest.mark.parametrize('lattice',['hP','cI','cF']) # tI not included yet
|
@pytest.mark.parametrize('lattice',['hP','cI','cF']) # tI not included yet
|
||||||
def test_Schmid_vectorize(self,lattice):
|
def test_Schmid_vectorization(self,lattice):
|
||||||
O = Orientation.from_random(shape=4,lattice=lattice) # noqa
|
O = Orientation.from_random(shape=4,lattice=lattice) # noqa
|
||||||
for mode in ['slip','twin']:
|
for mode in ['slip','twin']:
|
||||||
Ps = O.Schmid(N_slip='*') if mode == 'slip' else O.Schmid(N_twin='*')
|
Ps = O.Schmid(N_slip='*') if mode == 'slip' else O.Schmid(N_twin='*')
|
||||||
for i in range(4):
|
for i in range(4):
|
||||||
P = O[i].Schmid(N_slip='*') if mode == 'slip' else O[i].Schmid(N_twin='*')
|
P = O[i].Schmid(N_slip='*') if mode == 'slip' else O[i].Schmid(N_twin='*')
|
||||||
assert np.allclose(P,Ps[:,i])
|
assert np.allclose(P,Ps[:,i])
|
||||||
|
|
||||||
|
@pytest.mark.parametrize('family',crystal_families)
|
||||||
|
@pytest.mark.parametrize('shape',[(1),(2,3),(4,3,2)])
|
||||||
|
def test_reduced_vectorization(self,family,shape):
|
||||||
|
o = Orientation.from_random(family=family,shape=shape)
|
||||||
|
for r, theO in zip(o.reduced.flatten(),o.flatten()):
|
||||||
|
assert r == theO.reduced
|
||||||
|
|
||||||
|
|
||||||
|
@pytest.mark.parametrize('family',crystal_families)
|
||||||
|
@pytest.mark.parametrize('shape',[(1),(2,3),(4,3,2)])
|
||||||
|
@pytest.mark.parametrize('vector',np.array([[1,0,0],[1,2,3],[-1,1,-1]]))
|
||||||
|
@pytest.mark.parametrize('proper',[True,False])
|
||||||
|
def test_to_SST_vectorization(self,family,shape,vector,proper):
|
||||||
|
o = Orientation.from_random(family=family,shape=shape)
|
||||||
|
for r, theO in zip(o.to_SST(vector=vector,proper=proper).reshape((-1,3)),o.flatten()):
|
||||||
|
assert np.allclose(r,theO.to_SST(vector=vector,proper=proper))
|
||||||
|
|
||||||
|
@pytest.mark.parametrize('proper',[True,False])
|
||||||
|
@pytest.mark.parametrize('family',crystal_families)
|
||||||
|
def test_in_SST_vectorization(self,family,proper):
|
||||||
|
vecs = np.random.rand(20,4,3)
|
||||||
|
result = Orientation(family=family).in_SST(vecs,proper).flatten()
|
||||||
|
for r,v in zip(result,vecs.reshape((-1,3))):
|
||||||
|
assert np.all(r == Orientation(family=family).in_SST(v,proper))
|
||||||
|
|
||||||
|
@pytest.mark.parametrize('family',crystal_families)
|
||||||
|
@pytest.mark.parametrize('shape',[(1),(2,3),(4,3,2)])
|
||||||
|
@pytest.mark.parametrize('vector',np.array([[1,0,0],[1,2,3],[-1,1,-1]]))
|
||||||
|
@pytest.mark.parametrize('proper',[True,False])
|
||||||
|
@pytest.mark.parametrize('in_SST',[True,False])
|
||||||
|
def test_IPF_color_vectorization(self,family,shape,vector,proper,in_SST):
|
||||||
|
o = Orientation.from_random(family=family,shape=shape)
|
||||||
|
for r, theO in zip(o.IPF_color(vector,in_SST=in_SST,proper=proper).reshape((-1,3)),o.flatten()):
|
||||||
|
assert np.allclose(r,theO.IPF_color(vector,in_SST=in_SST,proper=proper))
|
||||||
|
|
||||||
|
@pytest.mark.parametrize('family',crystal_families)
|
||||||
|
def test_in_FZ_vectorization(self,set_of_rodrigues,family):
|
||||||
|
result = Orientation.from_Rodrigues_vector(rho=set_of_rodrigues.reshape((-1,4,4)),family=family).in_FZ.reshape(-1)
|
||||||
|
for r,rho in zip(result,set_of_rodrigues[:len(result)]):
|
||||||
|
assert r == Orientation.from_Rodrigues_vector(rho=rho,family=family).in_FZ
|
||||||
|
|
||||||
|
@pytest.mark.parametrize('family',crystal_families)
|
||||||
|
def test_in_disorientation_FZ_vectorization(self,set_of_rodrigues,family):
|
||||||
|
result = Orientation.from_Rodrigues_vector(rho=set_of_rodrigues.reshape((-1,4,4)),
|
||||||
|
family=family).in_disorientation_FZ.reshape(-1)
|
||||||
|
for r,rho in zip(result,set_of_rodrigues[:len(result)]):
|
||||||
|
assert r == Orientation.from_Rodrigues_vector(rho=rho,family=family).in_disorientation_FZ
|
||||||
|
|
||||||
|
@pytest.mark.parametrize('model',['Bain','KS','GT','GT_prime','NW','Pitsch'])
|
||||||
|
@pytest.mark.parametrize('lattice',['cF','cI'])
|
||||||
|
def test_relationship_vectorization(self,set_of_quaternions,lattice,model):
|
||||||
|
r = Orientation(rotation=set_of_quaternions[:200].reshape((50,4,4)),lattice=lattice).related(model)
|
||||||
|
for i in range(200):
|
||||||
|
assert (r.reshape((-1,200))[:,i] == Orientation(set_of_quaternions[i],lattice=lattice).related(model)).all()
|
||||||
|
|
||||||
|
### blending tests ###
|
||||||
|
|
||||||
|
@pytest.mark.parametrize('family',crystal_families)
|
||||||
|
@pytest.mark.parametrize('left,right',[
|
||||||
|
((2,3,2),(2,3,2)),
|
||||||
|
((2,2),(4,4)),
|
||||||
|
((3,1),(1,3)),
|
||||||
|
(None,None),
|
||||||
|
])
|
||||||
|
def test_disorientation_blending(self,family,left,right):
|
||||||
|
o = Orientation.from_random(family=family,shape=left)
|
||||||
|
p = Orientation.from_random(family=family,shape=right)
|
||||||
|
blend = util.shapeblender(o.shape,p.shape)
|
||||||
|
for loc in np.random.randint(0,blend,(10,len(blend))):
|
||||||
|
# print(f'{a}/{b} @ {loc}')
|
||||||
|
# print(o[tuple(loc[:len(o.shape)])].disorientation(p[tuple(loc[-len(p.shape):])]))
|
||||||
|
# print(o.disorientation(p)[tuple(loc)])
|
||||||
|
assert o[tuple(loc[:len(o.shape)])].disorientation(p[tuple(loc[-len(p.shape):])]) \
|
||||||
|
.isclose(o.disorientation(p)[tuple(loc)])
|
||||||
|
|
||||||
|
@pytest.mark.parametrize('family',crystal_families)
|
||||||
|
@pytest.mark.parametrize('left,right',[
|
||||||
|
((2,3,2),(2,3,2)),
|
||||||
|
((2,2),(4,4)),
|
||||||
|
((3,1),(1,3)),
|
||||||
|
(None,(3,)),
|
||||||
|
])
|
||||||
|
def test_IPF_color_blending(self,family,left,right):
|
||||||
|
o = Orientation.from_random(family=family,shape=left)
|
||||||
|
v = np.random.random(right+(3,))
|
||||||
|
blend = util.shapeblender(o.shape,v.shape[:-1])
|
||||||
|
for loc in np.random.randint(0,blend,(10,len(blend))):
|
||||||
|
assert np.allclose(o[tuple(loc[:len(o.shape)])].IPF_color(v[tuple(loc[-len(v.shape[:-1]):])]),
|
||||||
|
o.IPF_color(v)[tuple(loc)])
|
||||||
|
|
||||||
|
@pytest.mark.parametrize('family',crystal_families)
|
||||||
|
@pytest.mark.parametrize('left,right',[
|
||||||
|
((2,3,2),(2,3,2)),
|
||||||
|
((2,2),(4,4)),
|
||||||
|
((3,1),(1,3)),
|
||||||
|
(None,(3,)),
|
||||||
|
])
|
||||||
|
def test_to_SST_blending(self,family,left,right):
|
||||||
|
o = Orientation.from_random(family=family,shape=left)
|
||||||
|
v = np.random.random(right+(3,))
|
||||||
|
blend = util.shapeblender(o.shape,v.shape[:-1])
|
||||||
|
for loc in np.random.randint(0,blend,(10,len(blend))):
|
||||||
|
assert np.allclose(o[tuple(loc[:len(o.shape)])].to_SST(v[tuple(loc[-len(v.shape[:-1]):])]),
|
||||||
|
o.to_SST(v)[tuple(loc)])
|
||||||
|
|
||||||
|
@pytest.mark.parametrize('lattice,a,b,c,alpha,beta,gamma',
|
||||||
|
[
|
||||||
|
('aP',0.5,2.0,3.0,0.8,0.5,1.2),
|
||||||
|
('mP',1.0,2.0,3.0,np.pi/2,0.5,np.pi/2),
|
||||||
|
('oI',0.5,1.5,3.0,np.pi/2,np.pi/2,np.pi/2),
|
||||||
|
('tP',0.5,0.5,3.0,np.pi/2,np.pi/2,np.pi/2),
|
||||||
|
('hP',1.0,1.0,1.6,np.pi/2,np.pi/2,2*np.pi/3),
|
||||||
|
('cF',1.0,1.0,1.0,np.pi/2,np.pi/2,np.pi/2),
|
||||||
|
])
|
||||||
|
@pytest.mark.parametrize('left,right',[
|
||||||
|
((2,3,2),(2,3,2)),
|
||||||
|
((2,2),(4,4)),
|
||||||
|
((3,1),(1,3)),
|
||||||
|
(None,(3,)),
|
||||||
|
])
|
||||||
|
def test_to_pole_blending(self,lattice,a,b,c,alpha,beta,gamma,left,right):
|
||||||
|
o = Orientation.from_random(shape=left,
|
||||||
|
lattice=lattice,
|
||||||
|
a=a,b=b,c=c,
|
||||||
|
alpha=alpha,beta=beta,gamma=gamma)
|
||||||
|
v = np.random.random(right+(3,))
|
||||||
|
blend = util.shapeblender(o.shape,v.shape[:-1])
|
||||||
|
for loc in np.random.randint(0,blend,(10,len(blend))):
|
||||||
|
assert np.allclose(o[tuple(loc[:len(o.shape)])].to_pole(uvw=v[tuple(loc[-len(v.shape[:-1]):])]),
|
||||||
|
o.to_pole(uvw=v)[tuple(loc)])
|
||||||
|
|
Loading…
Reference in New Issue