following naming convention

This commit is contained in:
Martin Diehl 2022-01-29 15:30:59 +01:00
parent 487912cfb0
commit 762f93d724
10 changed files with 43 additions and 43 deletions

View File

@ -32,7 +32,7 @@ program DAMASK_grid
implicit none
type :: tLoadCase
type(rotation) :: rot !< rotation of BC
type(tRotation) :: rot !< rotation of BC
type(tBoundaryCondition) :: stress, & !< stress BC
deformation !< deformation BC (dot_F, F, or L)
real(pReal) :: t, & !< length of increment

View File

@ -339,7 +339,7 @@ subroutine grid_mechanical_FEM_forward(cutBack,guess,Delta_t,Delta_t_old,t_remai
type(tBoundaryCondition), intent(in) :: &
stress_BC, &
deformation_BC
type(rotation), intent(in) :: &
type(tRotation), intent(in) :: &
rotation_BC
PetscErrorCode :: err_PETSc
PetscScalar, pointer, dimension(:,:,:,:) :: &

View File

@ -310,7 +310,7 @@ subroutine grid_mechanical_spectral_basic_forward(cutBack,guess,Delta_t,Delta_t_
type(tBoundaryCondition), intent(in) :: &
stress_BC, &
deformation_BC
type(rotation), intent(in) :: &
type(tRotation), intent(in) :: &
rotation_BC
PetscErrorCode :: err_PETSc
PetscScalar, pointer, dimension(:,:,:,:) :: F

View File

@ -342,7 +342,7 @@ subroutine grid_mechanical_spectral_polarisation_forward(cutBack,guess,Delta_t,D
type(tBoundaryCondition), intent(in) :: &
stress_BC, &
deformation_BC
type(rotation), intent(in) :: &
type(tRotation), intent(in) :: &
rotation_BC
PetscErrorCode :: err_PETSc
PetscScalar, pointer, dimension(:,:,:,:) :: FandF_tau, F, F_tau

View File

@ -86,7 +86,7 @@ module spectral_utilities
type, public :: tSolutionParams
real(pReal), dimension(3,3) :: stress_BC
logical, dimension(3,3) :: stress_mask
type(rotation) :: rotation_BC
type(tRotation) :: rotation_BC
real(pReal) :: Delta_t
end type tSolutionParams
@ -666,7 +666,7 @@ function utilities_maskedCompliance(rot_BC,mask_stress,C)
real(pReal), dimension(3,3,3,3) :: utilities_maskedCompliance !< masked compliance
real(pReal), intent(in), dimension(3,3,3,3) :: C !< current average stiffness
type(rotation), intent(in) :: rot_BC !< rotation of load frame
type(tRotation), intent(in) :: rot_BC !< rotation of load frame
logical, intent(in), dimension(3,3) :: mask_stress !< mask of stress BC
integer :: i, j
@ -798,7 +798,7 @@ subroutine utilities_constitutiveResponse(P,P_av,C_volAvg,C_minmaxAvg,&
real(pReal), intent(out), dimension(3,3,cells(1),cells(2),cells3) :: P !< PK stress
real(pReal), intent(in), dimension(3,3,cells(1),cells(2),cells3) :: F !< deformation gradient target
real(pReal), intent(in) :: Delta_t !< loading time
type(rotation), intent(in), optional :: rotation_BC !< rotation of load frame
type(tRotation), intent(in), optional :: rotation_BC !< rotation of load frame
integer :: i

View File

@ -498,7 +498,7 @@ function lattice_C66_twin(Ntwin,C66,lattice,CoverA)
real(pReal), dimension(6,6,sum(Ntwin)) :: lattice_C66_twin
real(pReal), dimension(3,3,sum(Ntwin)):: coordinateSystem
type(rotation) :: R
type(tRotation) :: R
integer :: i
@ -537,7 +537,7 @@ function lattice_C66_trans(Ntrans,C_parent66,lattice_target, &
real(pReal), dimension(6,6) :: C_bar66, C_target_unrotated66
real(pReal), dimension(3,3,sum(Ntrans)) :: Q,S
type(rotation) :: R
type(tRotation) :: R
real(pReal) :: a_bcc, a_fcc, cOverA_trans
integer :: i
@ -599,7 +599,7 @@ function lattice_nonSchmidMatrix(Nslip,nonSchmidCoefficients,sense) result(nonSc
real(pReal), dimension(1:3,1:3,sum(Nslip)) :: coordinateSystem !< coordinate system of slip system
real(pReal), dimension(3) :: direction, normal, np
type(rotation) :: R
type(tRotation) :: R
integer :: i
@ -1976,7 +1976,7 @@ subroutine buildTransformationSystem(Q,S,Ntrans,cOverA,a_fcc,a_bcc)
a_bcc, & !< lattice parameter a for bcc target lattice
a_fcc !< lattice parameter a for fcc parent lattice
type(rotation) :: &
type(tRotation) :: &
R, & !< Pitsch rotation
B !< Rotation of fcc to Bain coordinate system
real(pReal), dimension(3,3) :: &

View File

@ -18,7 +18,7 @@ module material
private
type :: tRotationContainer
type(Rotation), dimension(:), allocatable :: data
type(tRotation), dimension(:), allocatable :: data
end type
type :: tTensorContainer
real(pReal), dimension(:,:,:), allocatable :: data

View File

@ -958,7 +958,7 @@ subroutine crystallite_results(group,ph)
!--------------------------------------------------------------------------------------------------
function to_quaternion(dataset)
type(rotation), dimension(:), intent(in) :: dataset
type(tRotation), dimension(:), intent(in) :: dataset
real(pReal), dimension(4,size(dataset,1)) :: to_quaternion
integer :: i

View File

@ -1395,7 +1395,7 @@ module subroutine plastic_nonlocal_updateCompatibility(orientation,ph,i,e)
nThresholdValues
logical, dimension(param(ph)%sum_N_sl) :: &
belowThreshold
type(rotation) :: mis
type(tRotation) :: mis
associate(prm => param(ph))

View File

@ -55,7 +55,7 @@ module rotations
real(pReal), parameter :: P = -1.0_pReal !< parameter for orientation conversion.
type, public :: rotation
type, public :: tRotation
real(pReal), dimension(4) :: q
contains
procedure, public :: asQuaternion
@ -78,7 +78,7 @@ module rotations
procedure, public :: rotStiffness
procedure, public :: misorientation
procedure, public :: standardize
end type rotation
end type tRotation
real(pReal), parameter :: &
PREF = sqrt(6.0_pReal/PI), &
@ -117,8 +117,8 @@ end subroutine rotations_init
!---------------------------------------------------------------------------------------------------
pure function asQuaternion(self)
class(rotation), intent(in) :: self
real(pReal), dimension(4) :: asQuaternion
class(tRotation), intent(in) :: self
real(pReal), dimension(4) :: asQuaternion
asQuaternion = self%q
@ -126,8 +126,8 @@ end function asQuaternion
!---------------------------------------------------------------------------------------------------
pure function asEulers(self)
class(rotation), intent(in) :: self
real(pReal), dimension(3) :: asEulers
class(tRotation), intent(in) :: self
real(pReal), dimension(3) :: asEulers
asEulers = qu2eu(self%q)
@ -135,8 +135,8 @@ end function asEulers
!---------------------------------------------------------------------------------------------------
pure function asAxisAngle(self)
class(rotation), intent(in) :: self
real(pReal), dimension(4) :: asAxisAngle
class(tRotation), intent(in) :: self
real(pReal), dimension(4) :: asAxisAngle
asAxisAngle = qu2ax(self%q)
@ -144,8 +144,8 @@ end function asAxisAngle
!---------------------------------------------------------------------------------------------------
pure function asMatrix(self)
class(rotation), intent(in) :: self
real(pReal), dimension(3,3) :: asMatrix
class(tRotation), intent(in) :: self
real(pReal), dimension(3,3) :: asMatrix
asMatrix = qu2om(self%q)
@ -153,8 +153,8 @@ end function asMatrix
!---------------------------------------------------------------------------------------------------
pure function asRodrigues(self)
class(rotation), intent(in) :: self
real(pReal), dimension(4) :: asRodrigues
class(tRotation), intent(in) :: self
real(pReal), dimension(4) :: asRodrigues
asRodrigues = qu2ro(self%q)
@ -162,8 +162,8 @@ end function asRodrigues
!---------------------------------------------------------------------------------------------------
pure function asHomochoric(self)
class(rotation), intent(in) :: self
real(pReal), dimension(3) :: asHomochoric
class(tRotation), intent(in) :: self
real(pReal), dimension(3) :: asHomochoric
asHomochoric = qu2ho(self%q)
@ -174,7 +174,7 @@ end function asHomochoric
!---------------------------------------------------------------------------------------------------
subroutine fromQuaternion(self,qu)
class(rotation), intent(out) :: self
class(tRotation), intent(out) :: self
real(pReal), dimension(4), intent(in) :: qu
if (dNeq(norm2(qu),1.0_pReal,1.0e-8_pReal)) call IO_error(402,ext_msg='fromQuaternion')
@ -185,7 +185,7 @@ end subroutine fromQuaternion
!---------------------------------------------------------------------------------------------------
subroutine fromEulers(self,eu,degrees)
class(rotation), intent(out) :: self
class(tRotation), intent(out) :: self
real(pReal), dimension(3), intent(in) :: eu
logical, intent(in), optional :: degrees
@ -206,7 +206,7 @@ end subroutine fromEulers
!---------------------------------------------------------------------------------------------------
subroutine fromAxisAngle(self,ax,degrees,P)
class(rotation), intent(out) :: self
class(tRotation), intent(out) :: self
real(pReal), dimension(4), intent(in) :: ax
logical, intent(in), optional :: degrees
integer, intent(in), optional :: P
@ -236,7 +236,7 @@ end subroutine fromAxisAngle
!---------------------------------------------------------------------------------------------------
subroutine fromMatrix(self,om)
class(rotation), intent(out) :: self
class(tRotation), intent(out) :: self
real(pReal), dimension(3,3), intent(in) :: om
if (dNeq(math_det33(om),1.0_pReal,tol=1.0e-5_pReal)) &
@ -253,10 +253,10 @@ end subroutine fromMatrix
!---------------------------------------------------------------------------------------------------
pure elemental function rotRot__(self,R) result(rRot)
type(rotation) :: rRot
class(rotation), intent(in) :: self,R
type(tRotation) :: rRot
class(tRotation), intent(in) :: self,R
rRot = rotation(multiply_quaternion(self%q,R%q))
rRot = tRotation(multiply_quaternion(self%q,R%q))
call rRot%standardize()
end function rotRot__
@ -267,7 +267,7 @@ end function rotRot__
!---------------------------------------------------------------------------------------------------
pure elemental subroutine standardize(self)
class(rotation), intent(inout) :: self
class(tRotation), intent(inout) :: self
if (sign(1.0_pReal,self%q(1)) < 0.0_pReal) self%q = - self%q
@ -281,7 +281,7 @@ end subroutine standardize
pure function rotVector(self,v,active) result(vRot)
real(pReal), dimension(3) :: vRot
class(rotation), intent(in) :: self
class(tRotation), intent(in) :: self
real(pReal), intent(in), dimension(3) :: v
logical, intent(in), optional :: active
@ -317,7 +317,7 @@ end function rotVector
pure function rotTensor2(self,T,active) result(tRot)
real(pReal), dimension(3,3) :: tRot
class(rotation), intent(in) :: self
class(tRotation), intent(in) :: self
real(pReal), intent(in), dimension(3,3) :: T
logical, intent(in), optional :: active
@ -346,7 +346,7 @@ end function rotTensor2
pure function rotTensor4(self,T,active) result(tRot)
real(pReal), dimension(3,3,3,3) :: tRot
class(rotation), intent(in) :: self
class(tRotation), intent(in) :: self
real(pReal), intent(in), dimension(3,3,3,3) :: T
logical, intent(in), optional :: active
@ -378,7 +378,7 @@ end function rotTensor4
pure function rotStiffness(self,C,active) result(cRot)
real(pReal), dimension(6,6) :: cRot
class(rotation), intent(in) :: self
class(tRotation), intent(in) :: self
real(pReal), intent(in), dimension(6,6) :: C
logical, intent(in), optional :: active
@ -415,8 +415,8 @@ end function rotStiffness
!---------------------------------------------------------------------------------------------------
pure elemental function misorientation(self,other)
type(rotation) :: misorientation
class(rotation), intent(in) :: self, other
type(tRotation) :: misorientation
class(tRotation), intent(in) :: self, other
misorientation%q = multiply_quaternion(other%q, conjugate_quaternion(self%q))
@ -1415,7 +1415,7 @@ end function conjugate_quaternion
!--------------------------------------------------------------------------------------------------
subroutine selfTest()
type(rotation) :: R
type(tRotation) :: R
real(pReal), dimension(4) :: qu, ax, ro
real(pReal), dimension(3) :: x, eu, ho, v3
real(pReal), dimension(3,3) :: om, t33