same names as in python

This commit is contained in:
Martin Diehl 2019-09-20 17:18:09 -07:00
parent 7e6265b128
commit ad83c8541d
6 changed files with 58 additions and 58 deletions

View File

@ -64,9 +64,9 @@ class Rotation:
def __repr__(self):
"""Orientation displayed as unit quaternion, rotation matrix, and Bunge-Euler angles."""
return '\n'.join([
'{}'.format(self.quaternion),
'Matrix:\n{}'.format( '\n'.join(['\t'.join(list(map(str,self.asMatrix()[i,:]))) for i in range(3)]) ),
'Bunge Eulers / deg: {}'.format('\t'.join(list(map(str,self.asEulers(degrees=True)))) ),
'{}'.format(self.quaternion),
'Matrix:\n{}'.format( '\n'.join(['\t'.join(list(map(str,self.asMatrix()[i,:]))) for i in range(3)]) ),
'Bunge Eulers / deg: {}'.format('\t'.join(list(map(str,self.asEulers(degrees=True)))) ),
])
def __mul__(self, other):

View File

@ -380,7 +380,7 @@ subroutine crystallite_init
do e = FEsolving_execElem(1),FEsolving_execElem(2)
myNcomponents = homogenization_Ngrains(material_homogenizationAt(e))
do i = FEsolving_execIP(1,e), FEsolving_execIP(2,e); do c = 1, myNcomponents
crystallite_Fp0(1:3,1:3,c,i,e) = material_orientation0(c,i,e)%asRotationMatrix() ! plastic def gradient reflects init orientation
crystallite_Fp0(1:3,1:3,c,i,e) = material_orientation0(c,i,e)%asMatrix() ! plastic def gradient reflects init orientation
crystallite_Fi0(1:3,1:3,c,i,e) = constitutive_initialFi(c,i,e)
crystallite_F0(1:3,1:3,c,i,e) = math_I3
crystallite_localPlasticity(c,i,e) = phase_localPlasticity(material_phaseAt(c,e))
@ -827,7 +827,7 @@ subroutine crystallite_orientations
do e = FEsolving_execElem(1),FEsolving_execElem(2)
do i = FEsolving_execIP(1,e),FEsolving_execIP(2,e)
do c = 1,homogenization_Ngrains(material_homogenizationAt(e))
call crystallite_orientation(c,i,e)%fromRotationMatrix(transpose(math_rotationalPart33(crystallite_Fe(1:3,1:3,c,i,e))))
call crystallite_orientation(c,i,e)%fromMatrix(transpose(math_rotationalPart33(crystallite_Fe(1:3,1:3,c,i,e))))
enddo; enddo; enddo
!$OMP END PARALLEL DO
@ -857,7 +857,7 @@ function crystallite_push33ToRef(ipc,ip,el, tensor33)
ip, &
ipc
T = matmul(material_orientation0(ipc,ip,el)%asRotationMatrix(), & ! ToDo: initial orientation correct?
T = matmul(material_orientation0(ipc,ip,el)%asMatrix(), & ! ToDo: initial orientation correct?
transpose(math_inv33(crystallite_subF(1:3,1:3,ipc,ip,el))))
crystallite_push33ToRef = matmul(transpose(T),matmul(tensor33,T))
@ -908,7 +908,7 @@ function crystallite_postResults(ipc, ip, el)
case (grainrotation_ID)
rot = material_orientation0(ipc,ip,el)%misorientation(crystallite_orientation(ipc,ip,el))
mySize = 4
crystallite_postResults(c+1:c+mySize) = rot%asAxisAnglePair()
crystallite_postResults(c+1:c+mySize) = rot%asAxisAngle()
crystallite_postResults(c+4) = inDeg * crystallite_postResults(c+4) ! angle in degree
! remark: tensor output is of the form 11,12,13, 21,22,23, 31,32,33

View File

@ -236,8 +236,8 @@ program DAMASK_spectral
do j = 1, 3
temp_valueVector(j) = IO_floatValue(line,chunkPos,i+k+j)
enddo
call R%fromEulerAngles(temp_valueVector(1:3),degrees=(l==1))
newLoadCase%rotation = R%asRotationMatrix()
call R%fromEulers(temp_valueVector(1:3),degrees=(l==1))
newLoadCase%rotation = R%asMatrix()
case('rotation','rot') ! assign values for the rotation matrix
temp_valueVector = 0.0_pReal
do j = 1, 9

View File

@ -939,8 +939,8 @@ function lattice_C66_twin(Ntwin,C66,structure,CoverA)
end select
do i = 1, sum(Ntwin)
call R%fromAxisAnglePair([coordinateSystem(1:3,2,i),PI],P=1) ! ToDo: Why always 180 deg?
lattice_C66_twin(1:6,1:6,i) = math_sym3333to66(math_rotate_forward3333(math_66toSym3333(C66),R%asRotationMatrix()))
call R%fromAxisAngle([coordinateSystem(1:3,2,i),PI],P=1) ! ToDo: Why always 180 deg?
lattice_C66_twin(1:6,1:6,i) = math_sym3333to66(math_rotate_forward3333(math_66toSym3333(C66),R%asMatrix()))
enddo
end function lattice_C66_twin
@ -1001,9 +1001,9 @@ function lattice_C66_trans(Ntrans,C_parent66,structure_target, &
call buildTransformationSystem(Q,S,Ntrans,CoverA_trans,a_fcc,a_bcc)
do i = 1, sum(Ntrans)
call R%fromRotationMatrix(Q(1:3,1:3,i))
call R%fromMatrix(Q(1:3,1:3,i))
lattice_C66_trans(1:6,1:6,i) &
= math_sym3333to66(math_rotate_forward3333(C_target_unrotated,R%asRotationMatrix()))
= math_sym3333to66(math_rotate_forward3333(C_target_unrotated,R%asMatrix()))
enddo
end function lattice_C66_trans
@ -1036,8 +1036,8 @@ function lattice_nonSchmidMatrix(Nslip,nonSchmidCoefficients,sense) result(nonSc
do i = 1,sum(Nslip)
direction = coordinateSystem(1:3,1,i)
normal = coordinateSystem(1:3,2,i)
call R%fromAxisAnglePair([direction,60.0_pReal],degrees=.true.,P=1)
np = R%rotVector(normal)
call R%fromAxisAngle([direction,60.0_pReal],degrees=.true.,P=1)
np = R%rotate(normal)
if (size(nonSchmidCoefficients)>0) nonSchmidMatrix(1:3,1:3,i) = nonSchmidMatrix(1:3,1:3,i) &
+ nonSchmidCoefficients(1) * math_outer(direction, np)
@ -2231,8 +2231,8 @@ subroutine buildTransformationSystem(Q,S,Ntrans,cOverA,a_fcc,a_bcc)
if (a_bcc > 0.0_pReal .and. dEq0(cOverA)) then ! fcc -> bcc transformation
do i = 1,sum(Ntrans)
call R%fromAxisAnglePair(LATTICE_FCCTOBCC_SYSTEMTRANS(:,i),degrees=.true.,P=1)
call B%fromAxisAnglePair(LATTICE_FCCTOBCC_BAINROT(:,i), degrees=.true.,P=1)
call R%fromAxisAngle(LATTICE_FCCTOBCC_SYSTEMTRANS(:,i),degrees=.true.,P=1)
call B%fromAxisAngle(LATTICE_FCCTOBCC_BAINROT(:,i), degrees=.true.,P=1)
x = real(LATTICE_FCCTOBCC_BAINVARIANT(1:3,i),pReal)
y = real(LATTICE_FCCTOBCC_BAINVARIANT(4:6,i),pReal)
z = real(LATTICE_FCCTOBCC_BAINVARIANT(7:9,i),pReal)
@ -2240,8 +2240,8 @@ subroutine buildTransformationSystem(Q,S,Ntrans,cOverA,a_fcc,a_bcc)
U = (a_bcc/a_fcc)*math_outer(x,x) &
+ (a_bcc/a_fcc)*math_outer(y,y) * sqrt(2.0_pReal) &
+ (a_bcc/a_fcc)*math_outer(z,z) * sqrt(2.0_pReal)
Q(1:3,1:3,i) = matmul(R%asRotationMatrix(),B%asRotationMatrix())
S(1:3,1:3,i) = matmul(R%asRotationMatrix(),U) - MATH_I3
Q(1:3,1:3,i) = matmul(R%asMatrix(),B%asMatrix())
S(1:3,1:3,i) = matmul(R%asMatrix(),U) - MATH_I3
enddo
elseif (cOverA > 0.0_pReal .and. dEq0(a_bcc)) then ! fcc -> hex transformation
ss = MATH_I3

View File

@ -718,7 +718,7 @@ subroutine material_parseTexture
Eulers(3) = IO_floatValue(strings(1),chunkPos,j+1)
end select
enddo
call texture_Eulers(t)%fromEulerAngles(Eulers,degrees=.true.)
call texture_Eulers(t)%fromEulers(Eulers,degrees=.true.)
if (config_texture(t)%keyExists('axes')) then
strings = config_texture(t)%getStrings('axes')
@ -741,7 +741,7 @@ subroutine material_parseTexture
end select
enddo
if(dNeq(math_det33(texture_transformation),1.0_pReal)) call IO_error(157,t)
call transformation%fromRotationMatrix(texture_transformation)
call transformation%fromMatrix(texture_transformation)
texture_Eulers(t) = texture_Eulers(t) * transformation
endif

View File

@ -59,14 +59,14 @@ module rotations
type(quaternion), private :: q
contains
procedure, public :: asQuaternion
procedure, public :: asEulerAngles
procedure, public :: asAxisAnglePair
procedure, public :: asRodriguesFrankVector
procedure, public :: asRotationMatrix
procedure, public :: asEulers
procedure, public :: asAxisAngle
procedure, public :: asRodrigues
procedure, public :: asMatrix
!------------------------------------------
procedure, public :: fromEulerAngles
procedure, public :: fromAxisAnglePair
procedure, public :: fromRotationMatrix
procedure, public :: fromEulers
procedure, public :: fromAxisAngle
procedure, public :: fromMatrix
!------------------------------------------
procedure, private :: rotRot__
generic, public :: operator(*) => rotRot__
@ -93,41 +93,41 @@ pure function asQuaternion(self)
end function asQuaternion
!---------------------------------------------------------------------------------------------------
pure function asEulerAngles(self)
pure function asEulers(self)
class(rotation), intent(in) :: self
real(pReal), dimension(3) :: asEulerAngles
real(pReal), dimension(3) :: asEulers
asEulerAngles = qu2eu(self%q%asArray())
asEulers = qu2eu(self%q%asArray())
end function asEulerAngles
end function asEulers
!---------------------------------------------------------------------------------------------------
pure function asAxisAnglePair(self)
pure function asAxisAngle(self)
class(rotation), intent(in) :: self
real(pReal), dimension(4) :: asAxisAnglePair
real(pReal), dimension(4) :: asAxisAngle
asAxisAnglePair = qu2ax(self%q%asArray())
asAxisAngle = qu2ax(self%q%asArray())
end function asAxisAnglePair
end function asAxisAngle
!---------------------------------------------------------------------------------------------------
pure function asRotationMatrix(self)
pure function asMatrix(self)
class(rotation), intent(in) :: self
real(pReal), dimension(3,3) :: asRotationMatrix
real(pReal), dimension(3,3) :: asMatrix
asRotationMatrix = qu2om(self%q%asArray())
asMatrix = qu2om(self%q%asArray())
end function asRotationMatrix
end function asMatrix
!---------------------------------------------------------------------------------------------------
pure function asRodriguesFrankVector(self)
pure function asRodrigues(self)
class(rotation), intent(in) :: self
real(pReal), dimension(4) :: asRodriguesFrankVector
real(pReal), dimension(4) :: asRodrigues
asRodriguesFrankVector = qu2ro(self%q%asArray())
asRodrigues = qu2ro(self%q%asArray())
end function asRodriguesFrankVector
end function asRodrigues
!---------------------------------------------------------------------------------------------------
pure function asHomochoric(self)
@ -141,7 +141,7 @@ end function asHomochoric
!---------------------------------------------------------------------------------------------------
! Initialize rotation from different representations
!---------------------------------------------------------------------------------------------------
subroutine fromEulerAngles(self,eu,degrees)
subroutine fromEulers(self,eu,degrees)
class(rotation), intent(out) :: self
real(pReal), dimension(3), intent(in) :: eu
@ -156,13 +156,13 @@ subroutine fromEulerAngles(self,eu,degrees)
endif
if (any(Eulers<0.0_pReal) .or. any(Eulers>2.0_pReal*PI) .or. Eulers(2) > PI) &
call IO_error(402,ext_msg='fromEulerAngles')
call IO_error(402,ext_msg='fromEulers')
self%q = eu2qu(Eulers)
end subroutine fromEulerAngles
end subroutine fromEulers
!---------------------------------------------------------------------------------------------------
subroutine fromAxisAnglePair(self,ax,degrees,P)
subroutine fromAxisAngle(self,ax,degrees,P)
class(rotation), intent(out) :: self
real(pReal), dimension(4), intent(in) :: ax
@ -182,27 +182,27 @@ subroutine fromAxisAnglePair(self,ax,degrees,P)
axis = ax(1:3)
else
axis = ax(1:3) * merge(-1.0_pReal,1.0_pReal,P == 1)
if(abs(P) /= 1) call IO_error(402,ext_msg='fromAxisAnglePair (P)')
if(abs(P) /= 1) call IO_error(402,ext_msg='fromAxisAngle (P)')
endif
if(dNeq(norm2(axis),1.0_pReal) .or. angle < 0.0_pReal .or. angle > PI) &
call IO_error(402,ext_msg='fromAxisAnglePair')
call IO_error(402,ext_msg='fromAxisAngle')
self%q = ax2qu([axis,angle])
end subroutine fromAxisAnglePair
end subroutine fromAxisAngle
!---------------------------------------------------------------------------------------------------
subroutine fromRotationMatrix(self,om)
subroutine fromMatrix(self,om)
class(rotation), intent(out) :: self
real(pReal), dimension(3,3), intent(in) :: om
if (dNeq(math_det33(om),1.0_pReal,tol=1.0e-5_pReal)) &
call IO_error(402,ext_msg='fromRotationMatrix')
call IO_error(402,ext_msg='fromMatrix')
self%q = om2qu(om)
end subroutine fromRotationMatrix
end subroutine fromMatrix
!---------------------------------------------------------------------------------------------------
@ -277,9 +277,9 @@ pure function rotTensor2(self,T,active) result(tRot)
endif
if (passive) then
tRot = matmul(matmul(self%asRotationMatrix(),T),transpose(self%asRotationMatrix()))
tRot = matmul(matmul(self%asMatrix(),T),transpose(self%asMatrix()))
else
tRot = matmul(matmul(transpose(self%asRotationMatrix()),T),self%asRotationMatrix())
tRot = matmul(matmul(transpose(self%asMatrix()),T),self%asMatrix())
endif
end function rotTensor2
@ -301,9 +301,9 @@ pure function rotTensor4(self,T,active) result(tRot)
integer :: i,j,k,l,m,n,o,p
if (present(active)) then
R = merge(transpose(self%asRotationMatrix()),self%asRotationMatrix(),active)
R = merge(transpose(self%asMatrix()),self%asMatrix(),active)
else
R = self%asRotationMatrix()
R = self%asMatrix()
endif
tRot = 0.0_pReal