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
2018-08-23 02:58:47 +05:30
! ToDo: use, intrinsic :: iso_fortran_env, only : I8 => int64, WP => real64
2012-03-06 20:22:48 +05:30
implicit none
private
2016-03-21 03:50:58 +05:30
#if (FLOAT==8)
2012-08-29 00:40:54 +05:30
integer , parameter , public :: pReal = 8 !< floating point double precision (was selected_real_kind(15,300), number with 15 significant digits, up to 1e+-300)
2012-08-31 01:56:28 +05:30
#else
2013-09-18 19:29:42 +05:30
NO SUITABLE PRECISION FOR REAL SELECTED , STOPPING COMPILATION
2012-02-13 19:38:07 +05:30
#endif
2007-03-20 19:25:22 +05:30
2012-08-28 21:38:17 +05:30
#if (INT==4)
integer , parameter , public :: pInt = 4 !< integer representation 32 bit (was selected_int_kind(9), number with at least up to +- 1e9)
#elif (INT==8)
integer , parameter , public :: pInt = 8 !< integer representation 64 bit (was selected_int_kind(12), number with at least up to +- 1e12)
2012-08-31 01:56:28 +05:30
#else
2013-09-18 19:29:42 +05:30
NO SUITABLE PRECISION FOR INTEGER SELECTED , STOPPING COMPILATION
2012-08-28 21:38:17 +05:30
#endif
2018-08-22 12:44:16 +05:30
integer , parameter , public :: pStringLen = 256 !< default string lenth
2012-08-28 21:38:17 +05:30
integer , parameter , public :: pLongInt = 8 !< integer representation 64 bit (was selected_int_kind(12), number with at least up to +- 1e12)
2013-09-14 16:29:35 +05:30
real ( pReal ) , parameter , public :: tol_math_check = 1.0e-8_pReal !< tolerance for internal math self-checks (rotation)
2012-08-28 21:38:17 +05:30
2015-08-31 22:00:04 +05:30
integer ( pInt ) , allocatable , dimension ( : ) :: realloc_lhs_test
2018-12-08 12:32:55 +05:30
type , public :: group_float !< variable length datatype used for storage of state
2015-05-28 22:32:23 +05:30
real ( pReal ) , dimension ( : ) , pointer :: p
2018-08-20 19:39:40 +05:30
end type group_float
2015-04-13 15:32:52 +05:30
2018-08-04 05:09:14 +05:30
type , public :: group_int
2015-05-28 22:32:23 +05:30
integer ( pInt ) , dimension ( : ) , pointer :: p
2018-08-04 05:09:14 +05:30
end type group_int
2013-02-25 18:43:52 +05:30
2014-04-15 14:50:38 +05:30
!http://stackoverflow.com/questions/3948210/can-i-have-a-pointer-to-an-item-in-an-allocatable-array
type , public :: tState
2014-12-09 17:42:53 +05:30
integer ( pInt ) :: &
2017-09-30 03:14:10 +05:30
sizeState = 0_pInt , & !< size of state
sizeDotState = 0_pInt , & !< size of dot state, i.e. state(1:sizeDot) follows time evolution by dotState rates
2018-08-14 03:57:51 +05:30
offsetDeltaState = 0_pInt , & !< index offset of delta state
sizeDeltaState = 0_pInt , & !< size of delta state, i.e. state(offset+1:offset+sizeDelta) follows time evolution by deltaState increments
2017-09-30 03:14:10 +05:30
sizePostResults = 0_pInt !< size of output data
2016-01-22 06:38:36 +05:30
real ( pReal ) , pointer , dimension ( : ) , contiguous :: &
2014-12-09 17:42:53 +05:30
atolState
2015-04-13 15:32:52 +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
2017-09-30 03:14:10 +05:30
state0 , &
2014-12-09 17:42:53 +05:30
state , & !< state
2017-09-30 03:14:10 +05:30
dotState , & !< rate of state change
deltaState !< increment of state change
2015-10-30 21:18:30 +05:30
real ( pReal ) , allocatable , dimension ( : , : ) :: &
2014-12-09 17:42:53 +05:30
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
integer ( pInt ) :: &
nSlip = 0_pInt , &
nTwin = 0_pInt , &
nTrans = 0_pInt
2015-05-28 22:32:23 +05:30
logical :: &
2016-07-27 17:59:16 +05:30
nonlocal = . false .
2018-05-16 03:39:33 +05:30
real ( pReal ) , pointer , dimension ( : , : ) :: &
2014-12-09 17:42:53 +05:30
slipRate , & !< slip rate
accumulatedSlip !< accumulated plastic slip
2014-04-15 14:50:38 +05:30
end type
2014-09-19 23:29:06 +05:30
2015-05-28 22:32:23 +05:30
type , public :: tSourceState
type ( tState ) , dimension ( : ) , allocatable :: p !< tState for each active source mechanism in a phase
end type
type , public :: tHomogMapping
integer ( pInt ) , pointer , dimension ( : , : ) :: p
end type
type , public :: tPhaseMapping
integer ( pInt ) , pointer , dimension ( : , : , : ) :: p
2014-08-21 23:18:20 +05:30
end type
2014-04-15 14:50:38 +05:30
2013-02-11 15:14:17 +05:30
public :: &
2015-04-13 15:32:52 +05:30
prec_init , &
2016-03-27 00:25:44 +05:30
dEq , &
2016-06-30 14:41:35 +05:30
dEq0 , &
2016-05-29 14:15:03 +05:30
cEq , &
dNeq , &
2016-06-30 14:41:35 +05:30
dNeq0 , &
2016-05-29 14:15:03 +05:30
cNeq
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
2018-02-02 17:06:09 +05:30
#if defined(__GFORTRAN__) || __INTEL_COMPILER >= 1800
2017-10-05 20:05:34 +05:30
use , intrinsic :: iso_fortran_env , only : &
compiler_version , &
compiler_options
#endif
2015-04-13 15:32:52 +05:30
2012-03-06 20:22:48 +05:30
implicit none
2013-02-13 23:24:56 +05:30
external :: &
2016-06-29 17:10:54 +05:30
quit
2014-10-10 18:38:34 +05:30
2016-06-29 17:10:54 +05:30
write ( 6 , '(/,a)' ) ' <<<+- prec init -+>>>'
2012-02-01 00:48:55 +05:30
#include "compilation_info.f90"
2016-06-29 17:10:54 +05:30
write ( 6 , '(a,i3)' ) ' Bytes for pReal: ' , pReal
write ( 6 , '(a,i3)' ) ' Bytes for pInt: ' , pInt
write ( 6 , '(a,i3)' ) ' Bytes for pLongInt: ' , pLongInt
2014-04-15 14:50:38 +05:30
2015-08-31 22:00:04 +05:30
realloc_lhs_test = [ 1_pInt , 2_pInt ]
if ( realloc_lhs_test ( 2 ) / = 2_pInt ) 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
implicit none
2016-05-29 14:15:03 +05:30
real ( pReal ) , intent ( in ) :: a , b
2016-05-25 11:22:56 +05:30
real ( pReal ) , intent ( in ) , optional :: tol
2016-05-29 14:15:03 +05:30
real ( pReal ) , parameter :: eps = 2.220446049250313E-16 ! DBL_EPSILON in C
2016-05-25 11:22:56 +05:30
2018-08-14 03:57:51 +05:30
dEq = merge ( . True . , . False . , abs ( a - b ) < = merge ( tol , eps , present ( tol ) ) * maxval ( abs ( [ a , b ] ) ) )
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
implicit none
real ( pReal ) , intent ( in ) :: a , b
real ( pReal ) , intent ( in ) , optional :: tol
2016-05-29 14:15:03 +05:30
real ( pReal ) , parameter :: eps = 2.220446049250313E-16 ! DBL_EPSILON in C
2016-05-25 11:22:56 +05:30
2018-08-14 03:57:51 +05:30
dNeq = merge ( . False . , . True . , abs ( a - b ) < = merge ( tol , eps , present ( tol ) ) * maxval ( abs ( [ a , b ] ) ) )
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 )
implicit none
real ( pReal ) , intent ( in ) :: a
real ( pReal ) , intent ( in ) , optional :: tol
2017-11-21 13:48:03 +05:30
real ( pReal ) , parameter :: eps = 2.2250738585072014E-308 ! smallest non-denormalized number
2016-06-30 14:00:40 +05:30
2018-08-14 03:57:51 +05:30
dEq0 = merge ( . True . , . False . , abs ( a ) < = merge ( tol , eps , present ( tol ) ) )
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 )
implicit none
real ( pReal ) , intent ( in ) :: a
real ( pReal ) , intent ( in ) , optional :: tol
2017-11-21 13:48:03 +05:30
real ( pReal ) , parameter :: eps = 2.2250738585072014E-308 ! smallest non-denormalized number
2016-06-30 14:00:40 +05:30
2018-08-14 03:57:51 +05:30
dNeq0 = merge ( . False . , . True . , abs ( a ) < = merge ( tol , eps , present ( tol ) ) )
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 )
implicit none
complex ( pReal ) , intent ( in ) :: a , b
real ( pReal ) , intent ( in ) , optional :: tol
real ( pReal ) , parameter :: eps = 2.220446049250313E-16 ! DBL_EPSILON in C
2018-08-14 03:57:51 +05:30
cEq = merge ( . True . , . False . , abs ( a - b ) < = merge ( tol , eps , present ( tol ) ) * maxval ( abs ( [ a , b ] ) ) )
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 )
implicit none
complex ( pReal ) , intent ( in ) :: a , b
real ( pReal ) , intent ( in ) , optional :: tol
real ( pReal ) , parameter :: eps = 2.220446049250313E-16 ! DBL_EPSILON in C
2018-08-14 03:57:51 +05:30
cNeq = merge ( . False . , . True . , abs ( a - b ) < = merge ( tol , eps , present ( tol ) ) * maxval ( abs ( [ a , b ] ) ) )
2016-05-29 14:15:03 +05:30
end function cNeq
2012-03-06 20:22:48 +05:30
end module prec