DAMASK_EICMD/src/polynomials.f90

189 lines
5.8 KiB
Fortran
Raw Normal View History

2022-01-31 19:35:15 +05:30
!--------------------------------------------------------------------------------------------------
!> @author Martin Diehl, KU Leuven
!> @brief Polynomial representation for variable data
!--------------------------------------------------------------------------------------------------
module polynomials
use prec
use IO
use YAML_parse
use YAML_types
implicit none
private
type, public :: tPolynomial
real(pReal), dimension(:), allocatable :: coef
real(pReal) :: x_ref = huge(0.0_pReal)
2022-01-31 19:35:15 +05:30
contains
procedure, public :: at => eval
procedure, public :: der1_at => eval_der1
end type tPolynomial
interface polynomial
module procedure polynomial_from_dict
module procedure polynomial_from_coef
end interface polynomial
public :: &
polynomial, &
polynomials_init
contains
!--------------------------------------------------------------------------------------------------
!> @brief Run self-test.
!--------------------------------------------------------------------------------------------------
subroutine polynomials_init()
print'(/,1x,a)', '<<<+- polynomials init -+>>>'; flush(IO_STDOUT)
call selfTest()
end subroutine polynomials_init
!--------------------------------------------------------------------------------------------------
!> @brief Initialize a Polynomial from Coefficients.
!--------------------------------------------------------------------------------------------------
2022-05-08 01:09:43 +05:30
pure function polynomial_from_coef(coef,x_ref) result(p)
2022-01-31 19:35:15 +05:30
2022-05-08 01:09:43 +05:30
real(pReal), dimension(0:), intent(in) :: coef
2022-01-31 19:35:15 +05:30
real(pReal), intent(in) :: x_ref
type(tPolynomial) :: p
2022-05-08 01:09:43 +05:30
p%coef = coef
2022-01-31 19:35:15 +05:30
p%x_ref = x_ref
end function polynomial_from_coef
!--------------------------------------------------------------------------------------------------
!> @brief Initialize a Polynomial from a Dictionary with Coefficients.
!--------------------------------------------------------------------------------------------------
function polynomial_from_dict(dict,y,x) result(p)
type(tDict), intent(in) :: dict
character(len=*), intent(in) :: y, x
type(tPolynomial) :: p
real(pReal), dimension(:), allocatable :: coef
real(pReal) :: x_ref
integer :: i
2022-01-31 19:35:15 +05:30
allocate(coef(1),source=dict%get_asFloat(y))
if (dict%contains(y//','//x)) then
x_ref = dict%get_asFloat(x//'_ref')
coef = [coef,dict%get_asFloat(y//','//x)]
end if
if (dict%contains(y//','//x//'^2')) then
x_ref = dict%get_asFloat(x//'_ref')
coef = [coef,[(0.0_pReal,i=size(coef),2-1)],dict%get_asFloat(y//','//x//'^2')]
2022-01-31 19:35:15 +05:30
end if
p = Polynomial(coef,x_ref)
end function polynomial_from_dict
!--------------------------------------------------------------------------------------------------
!> @brief Evaluate a Polynomial.
!--------------------------------------------------------------------------------------------------
pure function eval(self,x) result(y)
class(tPolynomial), intent(in) :: self
real(pReal), intent(in) :: x
real(pReal) :: y
integer :: i
y = self%coef(0)
do i = 1, ubound(self%coef,1)
y = y + self%coef(i) * (x-self%x_ref)**i
enddo
end function eval
!--------------------------------------------------------------------------------------------------
!> @brief Evaluate a first derivative of Polynomial.
!--------------------------------------------------------------------------------------------------
pure function eval_der1(self,x) result(y)
class(tPolynomial), intent(in) :: self
real(pReal), intent(in) :: x
real(pReal) :: y
integer :: i
y = 0.0_pReal
do i = 1, ubound(self%coef,1)
y = y + real(i,pReal)*self%coef(i) * (x-self%x_ref)**(i-1)
enddo
end function eval_der1
!--------------------------------------------------------------------------------------------------
!> @brief Check correctness of polynomical functionality.
!--------------------------------------------------------------------------------------------------
subroutine selfTest()
2022-01-31 19:35:15 +05:30
type(tPolynomial) :: p1, p2
real(pReal), dimension(3) :: coef
integer :: i
2022-01-31 19:35:15 +05:30
real(pReal) :: x_ref, x
class(tNode), pointer :: dict
character(len=pStringLen), dimension(size(coef)) :: coef_s
2022-01-31 19:35:15 +05:30
character(len=pStringLen) :: x_ref_s, x_s, YAML_s
2022-01-31 19:35:15 +05:30
call random_number(coef)
call random_number(x_ref)
call random_number(x)
coef = coef*10_pReal -0.5_pReal
x_ref = x_ref*10_pReal -0.5_pReal
x = x*10_pReal -0.5_pReal
p1 = polynomial(coef,x_ref)
if (dNeq(p1%at(x_ref),coef(1))) error stop 'polynomial: @ref'
do i = 1, size(coef_s)
write(coef_s(i),*) coef(i)
end do
2022-01-31 19:35:15 +05:30
write(x_ref_s,*) x_ref
write(x_s,*) x
YAML_s = 'C: '//trim(adjustl(coef_s(1)))//IO_EOL//&
'C,T: '//trim(adjustl(coef_s(2)))//IO_EOL//&
'C,T^2: '//trim(adjustl(coef_s(3)))//IO_EOL//&
'T_ref: '//trim(adjustl(x_ref_s))//IO_EOL
Dict => YAML_parse_str(trim(YAML_s))
p2 = polynomial(dict%asDict(),'C','T')
if (dNeq(p1%at(x),p2%at(x),1.0e-10_pReal)) error stop 'polynomials: init'
2022-01-31 19:35:15 +05:30
YAML_s = 'C: 0.0'//IO_EOL//&
'C,T: '//trim(adjustl(coef_s(2)))//IO_EOL//&
'T_ref: '//trim(adjustl(x_ref_s))//IO_EOL
Dict => YAML_parse_str(trim(YAML_s))
p1 = polynomial(dict%asDict(),'C','T')
2022-01-31 19:35:15 +05:30
if (dNeq(p1%at(x_ref+x),-p1%at(x_ref-x),1.0e-10_pReal)) error stop 'polynomials: eval(odd)'
if (dNeq(p1%der1_at(x),p1%der1_at(5.0_pReal*x),1.0e-10_pReal)) error stop 'polynomials: eval_der(odd)'
YAML_s = 'C: 0.0'//IO_EOL//&
'C,T^2: '//trim(adjustl(coef_s(3)))//IO_EOL//&
'T_ref: '//trim(adjustl(x_ref_s))//IO_EOL
Dict => YAML_parse_str(trim(YAML_s))
p1 = polynomial(dict%asDict(),'C','T')
2022-01-31 19:35:15 +05:30
if (dNeq(p1%at(x_ref+x),p1%at(x_ref-x),1e-10_pReal)) error stop 'polynomials: eval(even)'
if (dNeq(p1%der1_at(x_ref+x),-p1%der1_at(x_ref-x),1e-10_pReal)) error stop 'polynomials: eval_der(even)'
end subroutine selfTest
end module polynomials