more vectorized

This commit is contained in:
f.basile 2020-06-28 19:03:06 +02:00
parent 6fa5ae6ebf
commit 352c4e95f1
2 changed files with 95 additions and 50 deletions

View File

@ -77,27 +77,19 @@ class Orientation: # ToDo: make subclass of lattice and Rotation
# self-->other: True, self<--other: False
def inFZ_vec(self):
"""
Check if orientations falls into Fundamental Zone.
self.rotation.as_Rodrigues() working fine
self.rotation.as_Rodrigues(vector=True) doesn't work for several rotations
i apply dirty fix
"""
"""Check if orientations fall into Fundamental Zone."""
if not self.rotation.shape:
return self.lattice.symmetry.inFZ(self.rotation.as_Rodrigues(vector=True))
else:
return [self.lattice.symmetry.inFZ(\
Rotation._qu2ro(self.rotation.as_quaternion())[l][...,:3]\
*Rotation._qu2ro(self.rotation.as_quaternion())[l][...,3])\
for l in range(self.rotation.shape[0])]
self.rotation.as_Rodrigues(vector=True)[l]) for l in range(self.rotation.shape[0])]
def inFZ(self):
return self.lattice.symmetry.inFZ(self.rotation.as_Rodrigues(vector=True))
@property
def equivalent(self):
def equivalent_vec(self):
"""
Return orientations which are symmetrically equivalent.
@ -105,12 +97,12 @@ class Orientation: # ToDo: make subclass of lattice and Rotation
is added to the left of the rotation array.
"""
s = self.lattice.symmetry.symmetry_operations
s = s.reshape(s.shape[:1]+(1,)*len(self.rotation.shape)+(4,))
s = self.lattice.symmetry.symmetry_operations #24 lines (sym) x 4 columns (quat)
s = s.reshape(s.shape[:1]+(1,)*len(self.rotation.shape)+(4,)) #reshape zo (24,1,4)
s = Rotation(np.broadcast_to(s,s.shape[:1]+self.rotation.quaternion.shape))
r = np.broadcast_to(self.rotation.quaternion,s.shape[:1]+self.rotation.quaternion.shape)
r = Rotation(r)
r = np.broadcast_to(self.rotation.quaternion,s.shape[:1]+self.rotation.quaternion.shape) #(24,NumRots,4)
r = Rotation(r) #(24, NumRot)
return self.__class__(s@r,self.lattice)
@ -127,14 +119,17 @@ class Orientation: # ToDo: make subclass of lattice and Rotation
def relatedOrientations_vec(self,model):
"""List of orientations related by the given orientation relationship."""
r = self.lattice.relationOperations(model)
if not self.rotation.shape:
return [self.__class__(o*self.rotation,r['lattice']) for o in r['rotations']]
else:
return np.reshape(\
[self.__class__(o*Rotation.from_quaternion(self.rotation.as_quaternion()[l])\
,r['lattice']) for o in r['rotations'] for l in range(self.rotation.shape[0])]
,(len(r['rotations']),self.rotation.shape[0]))
h = self.lattice.relationOperations(model)
rot= h['rotations']
op=np.array([o.as_quaternion() for o in rot])
s = op.reshape(op.shape[:1]+(1,)*len(self.rotation.shape)+(4,))
s = Rotation(np.broadcast_to(s,s.shape[:1]+self.rotation.quaternion.shape))
r = np.broadcast_to(self.rotation.quaternion,s.shape[:1]+self.rotation.quaternion.shape)
r = Rotation(r)
return self.__class__(s@r,h['lattice'])
def relatedOrientations(self,model):
@ -142,6 +137,19 @@ class Orientation: # ToDo: make subclass of lattice and Rotation
r = self.lattice.relationOperations(model)
return [self.__class__(o*self.rotation,r['lattice']) for o in r['rotations']]
@property
def reduced_vec(self):
"""Transform orientation to fall into fundamental zone according to symmetry."""
equi=self.equivalent_vec.rotation #24 x rot x 3(rodrigues vector)
r= 1 if not self.rotation.shape else equi.shape[1] #number of rotations
quat=np.empty( [r , 4])
for rot in range(r):
for sym in range(equi.shape[0]):
if self.lattice.symmetry.inFZ(equi.as_Rodrigues(vector=True)[sym,rot]) == True:
quat[rot]=equi.as_quaternion()[sym,rot]
return self.__class__(quat,self.lattice)
def reduced(self):
"""Transform orientation to fall into fundamental zone according to symmetry."""
@ -178,9 +186,9 @@ class Orientation: # ToDo: make subclass of lattice and Rotation
return color
def IPF_color(self,axis):
"""TSL color of inverse pole figure for given axis."""
eq = self.equivalent
def IPFcolor_vec(self,axis):
"""TSL color of inverse pole figure for given axis. Not for hex or triclinic lattices"""
eq = self.equivalent_vec
pole = eq.rotation @ np.broadcast_to(axis/np.linalg.norm(axis),eq.rotation.shape+(3,))
in_SST, color = self.lattice.symmetry.in_SST(pole,color=True)

View File

@ -5,32 +5,36 @@ from damask import Rotation
from damask import Orientation
from damask import Lattice
rot0= Rotation.from_random()
rot1= Rotation.from_random()
rot2= Rotation.from_random()
rot3= Rotation.from_random()
rot0= Rotation.from_random() ;
rot1= Rotation.from_random() ;
rot2= Rotation.from_random() ;
rot3= Rotation.from_random() ;
#disorientation
#fromaverage
#average
class TestOrientation_vec:
@pytest.mark.xfail
#@pytest.mark.xfail
@pytest.mark.parametrize('lattice',Lattice.lattices)
def test_equivalentOrientations_vec(self,lattice):
def test_equivalent_vec(self,lattice):
ori0=Orientation(rot0,lattice)
ori1=Orientation(rot1,lattice)
ori2=Orientation(rot2,lattice)
ori3=Orientation(rot3,lattice)
quat=np.array([rot0.as_quaternion(),rot1.as_quaternion(),rot2.as_quaternion(),rot3.as_quaternion()])
rot_vec=Rotation.from_quaternion(quat)
ori_vec=Orientation(rot_vec,lattice)
ori_vec=Orientation(quat,lattice)
for s in range(len(ori_vec.lattice.symmetry.symmetryOperations())):
assert all(ori_vec.equivalent_vec()[s,0].rotation.as_Eulers() == \
assert all(ori_vec.equivalent_vec.rotation.as_Eulers()[s,0] == \
ori0.equivalentOrientations()[s].rotation.as_Eulers())
assert all(ori_vec.equivalent_vec()[s,1].rotation.as_quaternion() == \
assert all(ori_vec.equivalent_vec.rotation.as_quaternion()[s,1] == \
ori1.equivalentOrientations()[s].rotation.as_quaternion())
assert all(ori_vec.equivalent_vec()[s,2].rotation.as_Rodrigues() == \
assert all(ori_vec.equivalent_vec.rotation.as_Rodrigues()[s,2] == \
ori2.equivalentOrientations()[s].rotation.as_Rodrigues())
assert all(ori_vec.equivalent_vec()[s,3].rotation.as_cubochoric() == \
assert all(ori_vec.equivalent_vec.rotation.as_cubochoric()[s,3] == \
ori3.equivalentOrientations()[s].rotation.as_cubochoric())
@pytest.mark.parametrize('lattice',Lattice.lattices)
@ -39,14 +43,11 @@ class TestOrientation_vec:
ori1=Orientation(rot1,lattice)
ori2=Orientation(rot2,lattice)
ori3=Orientation(rot3,lattice)
#ensure 1 of them is in FZ
ori4=ori0.reduced()
rot4=ori4.rotation
ori4=ori0.reduced() ; rot4=ori4.rotation #ensure 1 of them is in FZ
quat=np.array([rot0.as_quaternion(),rot1.as_quaternion(),\
rot2.as_quaternion(),rot3.as_quaternion(), rot4.as_quaternion()])
rot_vec=Rotation.from_quaternion(quat)
ori_vec=Orientation(rot_vec,lattice)
ori_vec=Orientation(quat,lattice)
assert ori_vec.inFZ_vec()[0] == ori0.inFZ()
assert ori_vec.inFZ_vec()[1] == ori1.inFZ()
@ -64,16 +65,52 @@ class TestOrientation_vec:
ori3=Orientation(rot3,lattice)
quat=np.array([rot0.as_quaternion(),rot1.as_quaternion(),rot2.as_quaternion(),rot3.as_quaternion()])
rot_vec=Rotation.from_quaternion(quat)
ori_vec=Orientation(rot_vec,lattice)
ori_vec=Orientation(quat,lattice)
for s in range(len(ori1.lattice.relationOperations(model)['rotations'])):
assert all(ori_vec.relatedOrientations_vec(model)[s,0].rotation.as_Eulers() == \
assert all(ori_vec.relatedOrientations_vec(model).rotation.as_Eulers()[s,0] == \
ori0.relatedOrientations(model)[s].rotation.as_Eulers())
assert all(ori_vec.relatedOrientations_vec(model)[s,1].rotation.as_quaternion() == \
assert all(ori_vec.relatedOrientations_vec(model).rotation.as_quaternion()[s,1] == \
ori1.relatedOrientations(model)[s].rotation.as_quaternion())
assert all(ori_vec.relatedOrientations_vec(model)[s,2].rotation.as_Rodrigues() == \
assert all(ori_vec.relatedOrientations_vec(model).rotation.as_Rodrigues()[s,2] == \
ori2.relatedOrientations(model)[s].rotation.as_Rodrigues())
assert all(ori_vec.relatedOrientations_vec(model)[s,3].rotation.as_cubochoric() == \
assert all(ori_vec.relatedOrientations_vec(model).rotation.as_cubochoric()[s,3] == \
ori3.relatedOrientations(model)[s].rotation.as_cubochoric())
@pytest.mark.parametrize('lattice',Lattice.lattices)
def test_reduced_vec(self,lattice):
ori0=Orientation(rot0,lattice)
ori1=Orientation(rot1,lattice)
ori2=Orientation(rot2,lattice)
ori3=Orientation(rot3,lattice)
#ensure 1 of them is in FZ
ori4=ori0.reduced()
rot4=ori4.rotation
quat=np.array([rot0.as_quaternion(),rot1.as_quaternion(),\
rot2.as_quaternion(),rot3.as_quaternion(), rot4.as_quaternion()])
ori_vec=Orientation(quat,lattice)
assert all(ori_vec.reduced_vec.rotation.as_Eulers()[0] == ori0.reduced().rotation.as_Eulers() )
assert all(ori_vec.reduced_vec.rotation.as_quaternion()[1] == ori1.reduced().rotation.as_quaternion() )
assert all(ori_vec.reduced_vec.rotation.as_Rodrigues()[2] == ori2.reduced().rotation.as_Rodrigues() )
assert all(ori_vec.reduced_vec.rotation.as_cubochoric()[3] == ori3.reduced().rotation.as_cubochoric() )
assert all(ori_vec.reduced_vec.rotation.as_axis_angle()[4] == ori4.reduced().rotation.as_axis_angle() )
@pytest.mark.parametrize('lattice',['bcc','fcc','bct'])
def test_IPFcolor_vec(self,lattice):
ori0=Orientation(rot0,lattice)
ori1=Orientation(rot1,lattice)
ori2=Orientation(rot2,lattice)
ori3=Orientation(rot3,lattice)
quat=np.array([rot0.as_quaternion(),rot1.as_quaternion(),\
rot2.as_quaternion(),rot3.as_quaternion()])
ori_vec=Orientation(quat,lattice)
assert np.allclose( ori_vec.IPFcolor_vec(np.array([0,0,1]))[0],ori0.IPFcolor(np.array([0,0,1])))
assert np.allclose( ori_vec.IPFcolor_vec(np.array([0,2,1]))[1],ori1.IPFcolor(np.array([0,2,1])))
assert np.allclose( ori_vec.IPFcolor_vec(np.array([0,3,1]))[2],ori2.IPFcolor(np.array([0,3,1])))
assert np.allclose( ori_vec.IPFcolor_vec(np.array([4,0,1]))[3],ori3.IPFcolor(np.array([4,0,1])))