tau is often easier to use than pi

https://tauday.com/tau-manifesto
This commit is contained in:
Martin Diehl 2022-01-29 14:39:32 +01:00
parent 4a84c42112
commit dd0f2cfa3c
4 changed files with 39 additions and 37 deletions

View File

@ -936,36 +936,36 @@ pure function utilities_getFreqDerivative(k_s)
select case (spectral_derivative_ID) select case (spectral_derivative_ID)
case (DERIVATIVE_CONTINUOUS_ID) case (DERIVATIVE_CONTINUOUS_ID)
utilities_getFreqDerivative = cmplx(0.0_pReal, 2.0_pReal*PI*real(k_s,pReal)/geomSize,pReal) utilities_getFreqDerivative = cmplx(0.0_pReal, TAU*real(k_s,pReal)/geomSize,pReal)
case (DERIVATIVE_CENTRAL_DIFF_ID) case (DERIVATIVE_CENTRAL_DIFF_ID)
utilities_getFreqDerivative = cmplx(0.0_pReal, sin(2.0_pReal*PI*real(k_s,pReal)/real(grid,pReal)), pReal)/ & utilities_getFreqDerivative = cmplx(0.0_pReal, sin(TAU*real(k_s,pReal)/real(grid,pReal)), pReal)/ &
cmplx(2.0_pReal*geomSize/real(grid,pReal), 0.0_pReal, pReal) cmplx(2.0_pReal*geomSize/real(grid,pReal), 0.0_pReal, pReal)
case (DERIVATIVE_FWBW_DIFF_ID) case (DERIVATIVE_FWBW_DIFF_ID)
utilities_getFreqDerivative(1) = & utilities_getFreqDerivative(1) = &
cmplx(cos(2.0_pReal*PI*real(k_s(1),pReal)/real(grid(1),pReal)) - 1.0_pReal, & cmplx(cos(TAU*real(k_s(1),pReal)/real(grid(1),pReal)) - 1.0_pReal, &
sin(2.0_pReal*PI*real(k_s(1),pReal)/real(grid(1),pReal)), pReal)* & sin(TAU*real(k_s(1),pReal)/real(grid(1),pReal)), pReal)* &
cmplx(cos(2.0_pReal*PI*real(k_s(2),pReal)/real(grid(2),pReal)) + 1.0_pReal, & cmplx(cos(TAU*real(k_s(2),pReal)/real(grid(2),pReal)) + 1.0_pReal, &
sin(2.0_pReal*PI*real(k_s(2),pReal)/real(grid(2),pReal)), pReal)* & sin(TAU*real(k_s(2),pReal)/real(grid(2),pReal)), pReal)* &
cmplx(cos(2.0_pReal*PI*real(k_s(3),pReal)/real(grid(3),pReal)) + 1.0_pReal, & cmplx(cos(TAU*real(k_s(3),pReal)/real(grid(3),pReal)) + 1.0_pReal, &
sin(2.0_pReal*PI*real(k_s(3),pReal)/real(grid(3),pReal)), pReal)/ & sin(TAU*real(k_s(3),pReal)/real(grid(3),pReal)), pReal)/ &
cmplx(4.0_pReal*geomSize(1)/real(grid(1),pReal), 0.0_pReal, pReal) cmplx(4.0_pReal*geomSize(1)/real(grid(1),pReal), 0.0_pReal, pReal)
utilities_getFreqDerivative(2) = & utilities_getFreqDerivative(2) = &
cmplx(cos(2.0_pReal*PI*real(k_s(1),pReal)/real(grid(1),pReal)) + 1.0_pReal, & cmplx(cos(TAU*real(k_s(1),pReal)/real(grid(1),pReal)) + 1.0_pReal, &
sin(2.0_pReal*PI*real(k_s(1),pReal)/real(grid(1),pReal)), pReal)* & sin(TAU*real(k_s(1),pReal)/real(grid(1),pReal)), pReal)* &
cmplx(cos(2.0_pReal*PI*real(k_s(2),pReal)/real(grid(2),pReal)) - 1.0_pReal, & cmplx(cos(TAU*real(k_s(2),pReal)/real(grid(2),pReal)) - 1.0_pReal, &
sin(2.0_pReal*PI*real(k_s(2),pReal)/real(grid(2),pReal)), pReal)* & sin(TAU*real(k_s(2),pReal)/real(grid(2),pReal)), pReal)* &
cmplx(cos(2.0_pReal*PI*real(k_s(3),pReal)/real(grid(3),pReal)) + 1.0_pReal, & cmplx(cos(TAU*real(k_s(3),pReal)/real(grid(3),pReal)) + 1.0_pReal, &
sin(2.0_pReal*PI*real(k_s(3),pReal)/real(grid(3),pReal)), pReal)/ & sin(TAU*real(k_s(3),pReal)/real(grid(3),pReal)), pReal)/ &
cmplx(4.0_pReal*geomSize(2)/real(grid(2),pReal), 0.0_pReal, pReal) cmplx(4.0_pReal*geomSize(2)/real(grid(2),pReal), 0.0_pReal, pReal)
utilities_getFreqDerivative(3) = & utilities_getFreqDerivative(3) = &
cmplx(cos(2.0_pReal*PI*real(k_s(1),pReal)/real(grid(1),pReal)) + 1.0_pReal, & cmplx(cos(TAU*real(k_s(1),pReal)/real(grid(1),pReal)) + 1.0_pReal, &
sin(2.0_pReal*PI*real(k_s(1),pReal)/real(grid(1),pReal)), pReal)* & sin(TAU*real(k_s(1),pReal)/real(grid(1),pReal)), pReal)* &
cmplx(cos(2.0_pReal*PI*real(k_s(2),pReal)/real(grid(2),pReal)) + 1.0_pReal, & cmplx(cos(TAU*real(k_s(2),pReal)/real(grid(2),pReal)) + 1.0_pReal, &
sin(2.0_pReal*PI*real(k_s(2),pReal)/real(grid(2),pReal)), pReal)* & sin(TAU*real(k_s(2),pReal)/real(grid(2),pReal)), pReal)* &
cmplx(cos(2.0_pReal*PI*real(k_s(3),pReal)/real(grid(3),pReal)) - 1.0_pReal, & cmplx(cos(TAU*real(k_s(3),pReal)/real(grid(3),pReal)) - 1.0_pReal, &
sin(2.0_pReal*PI*real(k_s(3),pReal)/real(grid(3),pReal)), pReal)/ & sin(TAU*real(k_s(3),pReal)/real(grid(3),pReal)), pReal)/ &
cmplx(4.0_pReal*geomSize(3)/real(grid(3),pReal), 0.0_pReal, pReal) cmplx(4.0_pReal*geomSize(3)/real(grid(3),pReal), 0.0_pReal, pReal)
end select end select

View File

@ -21,10 +21,12 @@ module math
config config
#endif #endif
real(pReal), parameter :: PI = acos(-1.0_pReal) !< ratio of a circle's circumference to its diameter real(pReal), parameter :: &
real(pReal), parameter :: INDEG = 180.0_pReal/PI !< conversion from radian to degree PI = acos(-1.0_pReal), & !< ratio of a circle's circumference to its diameter
real(pReal), parameter :: INRAD = PI/180.0_pReal !< conversion from degree to radian TAU = 2.0_pReal*PI, & !< ratio of a circle's circumference to its radius
complex(pReal), parameter :: TWOPIIMG = cmplx(0.0_pReal,2.0_pReal*PI) !< Re(0.0), Im(2xPi) INDEG = 360.0_pReal/TAU, & !< conversion from radian to degree
INRAD = TAU/360.0_pReal !< conversion from degree to radian
complex(pReal), parameter :: TWOPIIMG = cmplx(0.0_pReal,TAU) !< Re(0.0), Im(Tau)
real(pReal), dimension(3,3), parameter :: & real(pReal), dimension(3,3), parameter :: &
math_I3 = reshape([& math_I3 = reshape([&
@ -984,7 +986,7 @@ impure elemental subroutine math_normal(x,mu,sigma)
end if end if
call random_number(rnd) call random_number(rnd)
x = mu_ + sigma_ * sqrt(-2.0_pReal*log(1.0_pReal-rnd(1)))*cos(2.0_pReal*PI*(1.0_pReal - rnd(2))) x = mu_ + sigma_ * sqrt(-2.0_pReal*log(1.0_pReal-rnd(1)))*cos(TAU*(1.0_pReal - rnd(2)))
end subroutine math_normal end subroutine math_normal
@ -1088,7 +1090,7 @@ pure function math_rotationalPart(F) result(R)
if (dNeq0(x)) then if (dNeq0(x)) then
Phi = acos(math_clip((I_C(1)**3 -4.5_pReal*I_C(1)*I_C(2) +13.5_pReal*I_C(3))/x,-1.0_pReal,1.0_pReal)) Phi = acos(math_clip((I_C(1)**3 -4.5_pReal*I_C(1)*I_C(2) +13.5_pReal*I_C(3))/x,-1.0_pReal,1.0_pReal))
lambda = I_C(1) +(2.0_pReal * sqrt(math_clip(I_C(1)**2-3.0_pReal*I_C(2),0.0_pReal))) & lambda = I_C(1) +(2.0_pReal * sqrt(math_clip(I_C(1)**2-3.0_pReal*I_C(2),0.0_pReal))) &
*cos((Phi-2.0_pReal * PI*[1.0_pReal,2.0_pReal,3.0_pReal])/3.0_pReal) *cos((Phi-TAU*[1.0_pReal,2.0_pReal,3.0_pReal])/3.0_pReal)
lambda = sqrt(math_clip(lambda,0.0_pReal)/3.0_pReal) lambda = sqrt(math_clip(lambda,0.0_pReal)/3.0_pReal)
else else
lambda = sqrt(I_C(1)/3.0_pReal) lambda = sqrt(I_C(1)/3.0_pReal)
@ -1154,8 +1156,8 @@ pure function math_eigvalsh33(m)
phi=acos(math_clip(-Q/rho*0.5_pReal,-1.0_pReal,1.0_pReal)) phi=acos(math_clip(-Q/rho*0.5_pReal,-1.0_pReal,1.0_pReal))
math_eigvalsh33 = 2.0_pReal*rho**(1.0_pReal/3.0_pReal)* & math_eigvalsh33 = 2.0_pReal*rho**(1.0_pReal/3.0_pReal)* &
[cos( phi /3.0_pReal), & [cos( phi /3.0_pReal), &
cos((phi+2.0_pReal*PI)/3.0_pReal), & cos((phi+TAU)/3.0_pReal), &
cos((phi+4.0_pReal*PI)/3.0_pReal) & cos((phi+2.0_pReal*TAU)/3.0_pReal) &
] & ] &
+ I(1)/3.0_pReal + I(1)/3.0_pReal
endif endif

View File

@ -341,7 +341,7 @@ module subroutine dislotungsten_dotState(Mp,T,ph,en)
dot_rho_dip_formation = merge(2.0_pReal*(d_hat-prm%d_caron)*stt%rho_mob(:,en)*dot%gamma_sl(:,en)/prm%b_sl, & dot_rho_dip_formation = merge(2.0_pReal*(d_hat-prm%d_caron)*stt%rho_mob(:,en)*dot%gamma_sl(:,en)/prm%b_sl, &
0.0_pReal, & 0.0_pReal, &
prm%dipoleformation) prm%dipoleformation)
v_cl = (3.0_pReal*mu*prm%D_0*exp(-prm%Q_cl/(K_B*T))*prm%f_at/(2.0_pReal*PI*K_B*T)) & v_cl = (3.0_pReal*mu*prm%D_0*exp(-prm%Q_cl/(K_B*T))*prm%f_at/(TAU*K_B*T)) &
* (1.0_pReal/(d_hat+prm%d_caron)) * (1.0_pReal/(d_hat+prm%d_caron))
dot_rho_dip_climb = (4.0_pReal*v_cl*stt%rho_dip(:,en))/(d_hat-prm%d_caron) ! ToDo: Discuss with Franz: Stress dependency? dot_rho_dip_climb = (4.0_pReal*v_cl*stt%rho_dip(:,en))/(d_hat-prm%d_caron) ! ToDo: Discuss with Franz: Stress dependency?
end where end where

View File

@ -198,7 +198,7 @@ subroutine fromEulers(self,eu,degrees)
Eulers = merge(eu*INRAD,eu,degrees) Eulers = merge(eu*INRAD,eu,degrees)
endif endif
if (any(Eulers<0.0_pReal) .or. any(Eulers>2.0_pReal*PI) .or. Eulers(2) > PI) & if (any(Eulers<0.0_pReal) .or. any(Eulers>TAU) .or. Eulers(2) > PI) &
call IO_error(402,ext_msg='fromEulers') call IO_error(402,ext_msg='fromEulers')
self%q = eu2qu(Eulers) self%q = eu2qu(Eulers)
@ -480,7 +480,7 @@ pure function qu2eu(qu) result(eu)
atan2( 2.0_pReal*chi, q03-q12 ), & atan2( 2.0_pReal*chi, q03-q12 ), &
atan2(( P*qu(1)*qu(3)+qu(2)*qu(4))*chi, (-P*qu(1)*qu(2)+qu(3)*qu(4))*chi )] atan2(( P*qu(1)*qu(3)+qu(2)*qu(4))*chi, (-P*qu(1)*qu(2)+qu(3)*qu(4))*chi )]
endif degenerated endif degenerated
where(sign(1.0_pReal,eu)<0.0_pReal) eu = mod(eu+2.0_pReal*PI,[2.0_pReal*PI,PI,2.0_pReal*PI]) where(sign(1.0_pReal,eu)<0.0_pReal) eu = mod(eu+TAU,[TAU,PI,TAU])
end function qu2eu end function qu2eu
@ -628,7 +628,7 @@ pure function om2eu(om) result(eu)
eu = [atan2(om(1,2),om(1,1)), 0.5_pReal*PI*(1.0_pReal-om(3,3)),0.0_pReal ] eu = [atan2(om(1,2),om(1,1)), 0.5_pReal*PI*(1.0_pReal-om(3,3)),0.0_pReal ]
end if end if
where(abs(eu) < 1.e-8_pReal) eu = 0.0_pReal where(abs(eu) < 1.e-8_pReal) eu = 0.0_pReal
where(sign(1.0_pReal,eu)<0.0_pReal) eu = mod(eu+2.0_pReal*PI,[2.0_pReal*PI,PI,2.0_pReal*PI]) where(sign(1.0_pReal,eu)<0.0_pReal) eu = mod(eu+TAU,[TAU,PI,TAU])
end function om2eu end function om2eu
@ -1437,7 +1437,7 @@ subroutine selfTest()
elseif(i==2) then elseif(i==2) then
qu = eu2qu([0.0_pReal,0.0_pReal,0.0_pReal]) qu = eu2qu([0.0_pReal,0.0_pReal,0.0_pReal])
elseif(i==3) then elseif(i==3) then
qu = eu2qu([2.0_pReal*PI,PI,2.0_pReal*PI]) qu = eu2qu([TAU,PI,TAU])
elseif(i==4) then elseif(i==4) then
qu = [0.0_pReal,0.0_pReal,1.0_pReal,0.0_pReal] qu = [0.0_pReal,0.0_pReal,1.0_pReal,0.0_pReal]
elseif(i==5) then elseif(i==5) then
@ -1448,10 +1448,10 @@ subroutine selfTest()
call random_number(x) call random_number(x)
A = sqrt(x(3)) A = sqrt(x(3))
B = sqrt(1-0_pReal -x(3)) B = sqrt(1-0_pReal -x(3))
qu = [cos(2.0_pReal*PI*x(1))*A,& qu = [cos(TAU*x(1))*A,&
sin(2.0_pReal*PI*x(2))*B,& sin(TAU*x(2))*B,&
cos(2.0_pReal*PI*x(2))*B,& cos(TAU*x(2))*B,&
sin(2.0_pReal*PI*x(1))*A] sin(TAU*x(1))*A]
if(qu(1)<0.0_pReal) qu = qu * (-1.0_pReal) if(qu(1)<0.0_pReal) qu = qu * (-1.0_pReal)
endif endif