added explicit method Rotation.fromBasis, which can treat real and reciprocal basis sets

This commit is contained in:
Philip Eisenlohr 2019-05-14 12:56:24 -04:00
parent 7c22c88df5
commit 4c7af713f1
1 changed files with 17 additions and 5 deletions

View File

@ -227,12 +227,17 @@ class Rotation:
return cls(ax2qu(ax))
@classmethod
def fromMatrix(cls,
matrix,
containsStretch = False): #ToDo: better name?
def fromBasis(cls,
basis,
orthonormal = True,
reciprocal = False,
):
om = matrix if isinstance(matrix, np.ndarray) else np.array(matrix).reshape((3,3)) # ToDo: Reshape here or require explicit?
if containsStretch:
om = basis if isinstance(basis, np.ndarray) else np.array(basis).reshape((3,3))
if reciprocal:
om = np.linalg.inv(om.T/np.pi) # transform reciprocal basis set
orthonormal = False # contains stretch
if not orthonormal:
(U,S,Vh) = np.linalg.svd(om) # singular value decomposition
om = np.dot(U,Vh)
if not np.isclose(np.linalg.det(om),1.0):
@ -244,6 +249,13 @@ class Rotation:
return cls(om2qu(om))
@classmethod
def fromMatrix(cls,
om,
):
return cls.fromBasis(om)
@classmethod
def fromRodrigues(cls,
rodrigues,