2012-08-09 16:31:53 +05:30
!--------------------------------------------------------------------------------------------------
2013-09-14 16:29:35 +05:30
!> @author Franz Roters, Max-Planck-Institut für Eisenforschung GmbH
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
!> @author Christoph Kords, Max-Planck-Institut für Eisenforschung GmbH
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
2014-05-08 20:25:19 +05:30
!> @author Luv Sharma, Max-Planck-Institut für Eisenforschung GmbH
2016-03-21 03:50:58 +05:30
!> @brief setting precision for real and int type
2012-08-09 16:31:53 +05:30
!--------------------------------------------------------------------------------------------------
2012-03-06 20:22:48 +05:30
module prec
2019-03-06 19:54:42 +05:30
use , intrinsic :: IEEE_arithmetic , only : &
IEEE_selected_real_kind
2019-03-07 23:59:48 +05:30
2019-03-06 19:54:42 +05:30
implicit none
2019-03-07 23:59:48 +05:30
private
2019-03-06 19:54:42 +05:30
! https://software.intel.com/en-us/blogs/2017/03/27/doctor-fortran-in-it-takes-all-kinds
2019-03-27 02:13:51 +05:30
#ifdef Abaqus
integer , parameter , public :: pReal = selected_real_kind ( 15 , 307 ) !< number with 15 significant digits, up to 1e+-307 (typically 64 bit)
#else
2019-03-07 23:07:22 +05:30
integer , parameter , public :: pReal = IEEE_selected_real_kind ( 15 , 307 ) !< number with 15 significant digits, up to 1e+-307 (typically 64 bit)
2019-03-27 02:13:51 +05:30
#endif
2019-03-06 20:22:52 +05:30
#if(INT==8)
2019-03-06 19:54:42 +05:30
integer , parameter , public :: pInt = selected_int_kind ( 18 ) !< number with at least up to +-1e18 (typically 64 bit)
2012-08-31 01:56:28 +05:30
#else
2019-03-06 20:22:52 +05:30
integer , parameter , public :: pInt = selected_int_kind ( 9 ) !< number with at least up to +-1e9 (typically 32 bit)
2012-08-28 21:38:17 +05:30
#endif
2019-03-06 19:54:42 +05:30
integer , parameter , public :: pLongInt = selected_int_kind ( 18 ) !< number with at least up to +-1e18 (typically 64 bit)
integer , parameter , public :: pStringLen = 256 !< default string length
real ( pReal ) , parameter , public :: tol_math_check = 1.0e-8_pReal !< tolerance for internal math self-checks (rotation)
type , public :: group_float !< variable length datatype used for storage of state
real ( pReal ) , dimension ( : ) , pointer :: p
end type group_float
type , public :: group_int
2019-03-07 23:59:48 +05:30
integer , dimension ( : ) , pointer :: p
2019-03-06 19:54:42 +05:30
end type group_int
! http://stackoverflow.com/questions/3948210/can-i-have-a-pointer-to-an-item-in-an-allocatable-array
type , public :: tState
2019-03-07 23:59:48 +05:30
integer :: &
sizeState = 0 , & !< size of state
sizeDotState = 0 , & !< size of dot state, i.e. state(1:sizeDot) follows time evolution by dotState rates
offsetDeltaState = 0 , & !< index offset of delta state
sizeDeltaState = 0 , & !< size of delta state, i.e. state(offset+1:offset+sizeDelta) follows time evolution by deltaState increments
sizePostResults = 0 !< size of output data
2019-03-06 19:54:42 +05:30
real ( pReal ) , pointer , dimension ( : ) , contiguous :: &
atolState
2019-03-07 23:59:48 +05:30
real ( pReal ) , pointer , dimension ( : , : ) , contiguous :: & ! a pointer is needed here because we might point to state/doState. However, they will never point to something, but are rather allocated and, hence, contiguous
2019-03-06 19:54:42 +05:30
state0 , &
state , & !< state
dotState , & !< rate of state change
deltaState !< increment of state change
real ( pReal ) , allocatable , dimension ( : , : ) :: &
partionedState0 , &
subState0 , &
previousDotState , & !< state rate of previous xxxx
previousDotState2 , & !< state rate two xxxx ago
RK4dotState
real ( pReal ) , allocatable , dimension ( : , : , : ) :: &
RKCK45dotState
end type
type , extends ( tState ) , public :: tPlasticState
2019-03-07 23:59:48 +05:30
integer :: &
nSlip = 0 , &
nTwin = 0 , &
nTrans = 0
logical :: &
2019-03-06 19:54:42 +05:30
nonlocal = . false .
real ( pReal ) , pointer , dimension ( : , : ) :: &
slipRate , & !< slip rate
accumulatedSlip !< accumulated plastic slip
end type
type , public :: tSourceState
type ( tState ) , dimension ( : ) , allocatable :: p !< tState for each active source mechanism in a phase
end type
2019-03-07 23:59:48 +05:30
2019-03-06 19:54:42 +05:30
type , public :: tHomogMapping
2019-03-07 23:59:48 +05:30
integer , pointer , dimension ( : , : ) :: p
end type
2019-03-06 19:54:42 +05:30
2019-03-07 23:59:48 +05:30
real ( pReal ) , private , parameter :: PREAL_EPSILON = epsilon ( 0.0_pReal ) !< minimum positive number such that 1.0 + EPSILON /= 1.0.
2019-03-07 15:32:27 +05:30
real ( pReal ) , private , parameter :: PREAL_MIN = tiny ( 0.0_pReal ) !< smallest normalized floating point number
2019-03-06 19:54:42 +05:30
public :: &
prec_init , &
dEq , &
dEq0 , &
cEq , &
dNeq , &
dNeq0 , &
cNeq
2019-03-07 23:59:48 +05:30
2012-03-06 20:22:48 +05:30
contains
2012-08-28 21:38:17 +05:30
2015-04-13 15:32:52 +05:30
2012-08-09 16:31:53 +05:30
!--------------------------------------------------------------------------------------------------
2017-05-04 04:02:44 +05:30
!> @brief reporting precision
2012-08-09 16:31:53 +05:30
!--------------------------------------------------------------------------------------------------
2012-03-06 20:22:48 +05:30
subroutine prec_init
2015-04-13 15:32:52 +05:30
2019-03-06 19:54:42 +05:30
implicit none
2019-03-07 23:59:48 +05:30
integer , allocatable , dimension ( : ) :: realloc_lhs_test
2019-03-06 19:54:42 +05:30
external :: &
quit
write ( 6 , '(/,a)' ) ' <<<+- prec init -+>>>'
2014-10-10 18:38:34 +05:30
2019-03-07 23:59:48 +05:30
write ( 6 , '(a,i3)' ) ' Size of integer in bit: ' , bit_size ( 0 )
write ( 6 , '(a,i19)' ) ' Maximum value: ' , huge ( 0 )
2019-03-06 19:54:42 +05:30
write ( 6 , '(/,a,i3)' ) ' Size of float in bit: ' , storage_size ( 0.0_pReal )
write ( 6 , '(a,e10.3)' ) ' Maximum value: ' , huge ( 0.0_pReal )
write ( 6 , '(a,e10.3)' ) ' Minimum value: ' , tiny ( 0.0_pReal )
write ( 6 , '(a,i3)' ) ' Decimal precision: ' , precision ( 0.0_pReal )
2014-04-15 14:50:38 +05:30
2019-03-07 23:59:48 +05:30
realloc_lhs_test = [ 1 , 2 ]
if ( realloc_lhs_test ( 2 ) / = 2 ) call quit ( 9000 )
2012-03-06 20:22:48 +05:30
end subroutine prec_init
2009-08-31 20:39:15 +05:30
2015-04-13 15:32:52 +05:30
2016-03-06 02:07:22 +05:30
!--------------------------------------------------------------------------------------------------
2016-05-29 14:15:03 +05:30
!> @brief equality comparison for float with double precision
2016-03-21 03:50:58 +05:30
! replaces "==" but for certain (relative) tolerance. Counterpart to dNeq
2016-10-18 23:20:05 +05:30
! https://randomascii.wordpress.com/2012/02/25/comparing-floating-point-numbers-2012-edition/
2017-11-21 13:48:03 +05:30
! AlmostEqualRelative
2016-03-06 02:07:22 +05:30
!--------------------------------------------------------------------------------------------------
logical elemental pure function dEq ( a , b , tol )
2016-05-25 11:22:56 +05:30
2019-03-06 20:27:42 +05:30
implicit none
real ( pReal ) , intent ( in ) :: a , b
real ( pReal ) , intent ( in ) , optional :: tol
real ( pReal ) :: eps
2019-03-07 23:59:48 +05:30
2019-03-06 20:27:42 +05:30
if ( present ( tol ) ) then
eps = tol
else
2019-03-07 15:32:27 +05:30
eps = PREAL_EPSILON * maxval ( abs ( [ a , b ] ) )
2019-03-06 20:27:42 +05:30
endif
dEq = merge ( . True . , . False . , abs ( a - b ) < eps )
2016-05-25 11:22:56 +05:30
2016-03-06 02:07:22 +05:30
end function dEq
!--------------------------------------------------------------------------------------------------
2016-05-29 14:15:03 +05:30
!> @brief inequality comparison for float with double precision
2016-03-21 03:50:58 +05:30
! replaces "!=" but for certain (relative) tolerance. Counterpart to dEq
2016-10-18 23:20:05 +05:30
! https://randomascii.wordpress.com/2012/02/25/comparing-floating-point-numbers-2012-edition/
2017-11-21 13:48:03 +05:30
! AlmostEqualRelative NOT
2016-03-06 02:07:22 +05:30
!--------------------------------------------------------------------------------------------------
logical elemental pure function dNeq ( a , b , tol )
2016-05-25 11:22:56 +05:30
2019-03-06 20:27:42 +05:30
implicit none
real ( pReal ) , intent ( in ) :: a , b
real ( pReal ) , intent ( in ) , optional :: tol
real ( pReal ) :: eps
2019-03-07 23:59:48 +05:30
2019-03-06 20:27:42 +05:30
if ( present ( tol ) ) then
eps = tol
else
2019-03-07 15:32:27 +05:30
eps = PREAL_EPSILON * maxval ( abs ( [ a , b ] ) )
2019-03-06 20:27:42 +05:30
endif
dNeq = merge ( . False . , . True . , abs ( a - b ) < = eps )
2016-05-25 11:22:56 +05:30
2016-03-06 02:07:22 +05:30
end function dNeq
2016-06-30 14:00:40 +05:30
!--------------------------------------------------------------------------------------------------
2016-10-18 23:20:05 +05:30
!> @brief equality to 0 comparison for float with double precision
2017-11-21 13:48:03 +05:30
! replaces "==0" but everything not representable as a normal number is treated as 0. Counterpart to dNeq0
! https://de.mathworks.com/help/matlab/ref/realmin.html
! https://docs.oracle.com/cd/E19957-01/806-3568/ncg_math.html
2016-06-30 14:00:40 +05:30
!--------------------------------------------------------------------------------------------------
logical elemental pure function dEq0 ( a , tol )
2019-03-06 20:27:42 +05:30
implicit none
real ( pReal ) , intent ( in ) :: a
real ( pReal ) , intent ( in ) , optional :: tol
real ( pReal ) :: eps
2019-03-07 23:59:48 +05:30
2019-03-06 20:27:42 +05:30
if ( present ( tol ) ) then
eps = tol
else
2019-03-07 15:32:27 +05:30
eps = PREAL_MIN * 1 0.0_pReal
2019-03-06 20:27:42 +05:30
endif
dEq0 = merge ( . True . , . False . , abs ( a ) < eps )
2016-06-30 14:00:40 +05:30
end function dEq0
!--------------------------------------------------------------------------------------------------
2016-10-18 23:20:05 +05:30
!> @brief inequality to 0 comparison for float with double precision
2017-11-21 13:48:03 +05:30
! replaces "!=0" but everything not representable as a normal number is treated as 0. Counterpart to dEq0
! https://de.mathworks.com/help/matlab/ref/realmin.html
! https://docs.oracle.com/cd/E19957-01/806-3568/ncg_math.html
2016-06-30 14:00:40 +05:30
!--------------------------------------------------------------------------------------------------
logical elemental pure function dNeq0 ( a , tol )
2019-03-06 20:27:42 +05:30
implicit none
real ( pReal ) , intent ( in ) :: a
real ( pReal ) , intent ( in ) , optional :: tol
real ( pReal ) :: eps
2019-03-07 23:59:48 +05:30
2019-03-06 20:27:42 +05:30
if ( present ( tol ) ) then
eps = tol
else
2019-03-07 15:32:27 +05:30
eps = PREAL_MIN * 1 0.0_pReal
2019-03-06 20:27:42 +05:30
endif
dNeq0 = merge ( . False . , . True . , abs ( a ) < = eps )
2016-06-30 14:00:40 +05:30
end function dNeq0
2016-05-29 14:15:03 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief equality comparison for complex with double precision
! replaces "==" but for certain (relative) tolerance. Counterpart to cNeq
2016-10-18 23:20:05 +05:30
! https://randomascii.wordpress.com/2012/02/25/comparing-floating-point-numbers-2012-edition/
2016-05-29 14:15:03 +05:30
! probably a component wise comparison would be more accurate than the comparsion of the absolute
! value
!--------------------------------------------------------------------------------------------------
logical elemental pure function cEq ( a , b , tol )
2019-03-06 20:27:42 +05:30
implicit none
complex ( pReal ) , intent ( in ) :: a , b
real ( pReal ) , intent ( in ) , optional :: tol
real ( pReal ) :: eps
2019-03-07 23:59:48 +05:30
2019-03-06 20:27:42 +05:30
if ( present ( tol ) ) then
eps = tol
else
2019-03-07 15:32:27 +05:30
eps = PREAL_EPSILON * maxval ( abs ( [ a , b ] ) )
2019-03-06 20:27:42 +05:30
endif
cEq = merge ( . True . , . False . , abs ( a - b ) < eps )
2016-05-29 14:15:03 +05:30
end function cEq
!--------------------------------------------------------------------------------------------------
!> @brief inequality comparison for complex with double precision
! replaces "!=" but for certain (relative) tolerance. Counterpart to cEq
2016-10-18 23:20:05 +05:30
! https://randomascii.wordpress.com/2012/02/25/comparing-floating-point-numbers-2012-edition/
2016-05-29 14:15:03 +05:30
! probably a component wise comparison would be more accurate than the comparsion of the absolute
! value
!--------------------------------------------------------------------------------------------------
logical elemental pure function cNeq ( a , b , tol )
2019-03-06 20:27:42 +05:30
implicit none
complex ( pReal ) , intent ( in ) :: a , b
real ( pReal ) , intent ( in ) , optional :: tol
real ( pReal ) :: eps
2019-03-07 23:59:48 +05:30
2019-03-06 20:27:42 +05:30
if ( present ( tol ) ) then
eps = tol
else
2019-03-07 15:32:27 +05:30
eps = PREAL_EPSILON * maxval ( abs ( [ a , b ] ) )
2019-03-06 20:27:42 +05:30
endif
cNeq = merge ( . False . , . True . , abs ( a - b ) < = eps )
2016-05-29 14:15:03 +05:30
end function cNeq
2012-03-06 20:22:48 +05:30
end module prec