Skip to content

Fix erroneous gaussian quadrature points in gauss_legendre #660

New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Merged
merged 2 commits into from
Jun 13, 2022
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 0 additions & 2 deletions src/stdlib_quadrature_gauss.f90
Original file line number Diff line number Diff line change
Expand Up @@ -56,8 +56,6 @@ pure module subroutine gauss_legendre_fp64 (x, w, interval)
if (present(interval)) then
associate ( a => interval(1) , b => interval(2) )
x = 0.5_dp*(b-a)*x+0.5_dp*(b+a)
x(1) = interval(1)
x(size(x)) = interval(2)
w = 0.5_dp*(b-a)*w
end associate
end if
Expand Down
22 changes: 21 additions & 1 deletion src/tests/quadrature/test_gauss.f90
Original file line number Diff line number Diff line change
Expand Up @@ -21,7 +21,8 @@ subroutine collect_gauss(testsuite)
new_unittest("gauss-lobatto-analytic", test_gauss_lobatto_analytic), &
new_unittest("gauss-lobatto-5", test_gauss_lobatto_5), &
new_unittest("gauss-lobatto-32", test_gauss_lobatto_32), &
new_unittest("gauss-lobatto-64", test_gauss_lobatto_64) &
new_unittest("gauss-lobatto-64", test_gauss_lobatto_64), &
new_unittest("gauss-github-issue-619", test_fix_github_issue619) &
]
end subroutine

Expand All @@ -48,6 +49,25 @@ subroutine test_gauss_analytic(error)

end subroutine

subroutine test_fix_github_issue619(error)
!> See github issue https://github.com/fortran-lang/stdlib/issues/619
type(error_type), allocatable, intent(out) :: error
integer :: i

! test the values of nodes and weights
i = 5
block
real(dp), dimension(i) :: x1,w1,x2,w2
call gauss_legendre(x1,w1)
call gauss_legendre(x2,w2,interval=[-1._dp, 1._dp])

call check(error, all(abs(x1-x2) < 2*epsilon(x1(1))))
if (allocated(error)) return
call check(error, all(abs(w1-w2) < 2*epsilon(w1(1))))
end block

end subroutine

subroutine test_gauss_5(error)
!> Error handling
type(error_type), allocatable, intent(out) :: error
Expand Down