-
Notifications
You must be signed in to change notification settings - Fork 191
First steps toward quadrature #146
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
Changes from all commits
Commits
Show all changes
12 commits
Select commit
Hold shift + click to select a range
97ce9d3
initial commit -- trapz
nshaffer 77a31ae
bugfix? redundant use stmt
nshaffer 32546f1
Merge branch 'master' of https://github.com/fortran-lang/stdlib into …
nshaffer 52386a6
Add simpson's rule API (not implemented)
nshaffer b1fb265
Add initial spec for trapz and trapz_weights
nshaffer 0b385e7
Add spec for simps
nshaffer 770671d
typos
nshaffer 8b6dcb7
minor wording changes
nshaffer f7d732d
add spec for simps_weights
nshaffer d1c91ec
replace f.p. equality checks with tolerance checks
nshaffer ccaf548
fix compile error
nshaffer 5e8294c
Merge branch 'master' into dev-quadrature
milancurcic File filter
Filter by extension
Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
There are no files selected for viewing
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,75 @@ | ||
#:include "common.fypp" | ||
|
||
module stdlib_experimental_quadrature | ||
use stdlib_experimental_kinds, only: sp, dp, qp | ||
|
||
implicit none | ||
|
||
private | ||
|
||
! array integration | ||
public :: trapz | ||
public :: trapz_weights | ||
public :: simps | ||
public :: simps_weights | ||
|
||
|
||
interface trapz | ||
#:for KIND in REAL_KINDS | ||
pure module function trapz_dx_${KIND}$(y, dx) result(integral) | ||
real(${KIND}$), dimension(:), intent(in) :: y | ||
real(${KIND}$), intent(in) :: dx | ||
real(${KIND}$) :: integral | ||
end function trapz_dx_${KIND}$ | ||
#:endfor | ||
#:for KIND in REAL_KINDS | ||
pure module function trapz_x_${KIND}$(y, x) result(integral) | ||
real(${KIND}$), dimension(:), intent(in) :: y | ||
real(${KIND}$), dimension(size(y)), intent(in) :: x | ||
real(${KIND}$) :: integral | ||
end function trapz_x_${KIND}$ | ||
#:endfor | ||
end interface trapz | ||
|
||
|
||
interface trapz_weights | ||
#:for KIND in REAL_KINDS | ||
pure module function trapz_weights_${KIND}$(x) result(w) | ||
real(${KIND}$), dimension(:), intent(in) :: x | ||
real(${KIND}$), dimension(size(x)) :: w | ||
end function trapz_weights_${KIND}$ | ||
#:endfor | ||
end interface trapz_weights | ||
|
||
|
||
interface simps | ||
#:for KIND in REAL_KINDS | ||
pure module function simps_dx_${KIND}$(y, dx, even) result(integral) | ||
real(${KIND}$), dimension(:), intent(in) :: y | ||
real(${KIND}$), intent(in) :: dx | ||
integer, intent(in), optional :: even | ||
real(${KIND}$) :: integral | ||
end function simps_dx_${KIND}$ | ||
#:endfor | ||
#:for KIND in REAL_KINDS | ||
pure module function simps_x_${KIND}$(y, x, even) result(integral) | ||
real(${KIND}$), dimension(:), intent(in) :: y | ||
real(${KIND}$), dimension(size(y)), intent(in) :: x | ||
integer, intent(in), optional :: even | ||
real(${KIND}$) :: integral | ||
end function simps_x_${KIND}$ | ||
#:endfor | ||
end interface simps | ||
|
||
|
||
interface simps_weights | ||
#:for KIND in REAL_KINDS | ||
pure module function simps_weights_${KIND}$(x, even) result(w) | ||
real(${KIND}$), dimension(:), intent(in) :: x | ||
real(${KIND}$), dimension(size(x)) :: w | ||
integer, intent(in), optional :: even | ||
end function simps_weights_${KIND}$ | ||
#:endfor | ||
end interface simps_weights | ||
|
||
end module stdlib_experimental_quadrature |
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,141 @@ | ||
# Numerical integration | ||
|
||
## Implemented | ||
|
||
* `trapz` | ||
* `trapz_weights` | ||
|
||
## `trapz` - integrate sampled values using trapezoidal rule | ||
|
||
Returns the trapezoidal rule integral of an array `y` representing discrete samples of a function. The integral is computed assuming either equidistant abscissas with spacing `dx` or arbitary abscissas `x`. | ||
|
||
### Syntax | ||
|
||
`result = trapz(y, x)` | ||
|
||
`result = trapz(y, dx)` | ||
|
||
### Arguments | ||
|
||
`y`: Shall be a rank-one array of type `real`. | ||
|
||
`x`: Shall be a rank-one array of type `real` having the same kind and size as `y`. | ||
|
||
`dx`: Shall be a scalar of type `real` having the same kind as `y`. | ||
|
||
### Return value | ||
|
||
The result is a scalar of type `real` having the same kind as `y`. | ||
|
||
If the size of `y` is zero or one, the result is zero. | ||
|
||
### Example | ||
|
||
```fortran | ||
program demo_trapz | ||
use stdlib_experimental_quadrature, only: trapz | ||
implicit none | ||
real :: x(5) = [0., 1., 2., 3., 4.] | ||
real :: y(5) = x**2 | ||
print *, trapz(y, x) | ||
! 22.0 | ||
print *, trapz(y, 0.5) | ||
! 11.0 | ||
end program | ||
``` | ||
|
||
## `trapz_weights` - trapezoidal rule weights for given abscissas | ||
|
||
Given an array of abscissas `x`, computes the array of weights `w` such that if `y` represented function values tabulated at `x`, then `sum(w*y)` produces a trapezoidal rule approximation to the integral. | ||
|
||
### Syntax | ||
|
||
`result = trapz_weights(x)` | ||
|
||
### Arguments | ||
|
||
`x`: Shall be a rank-one array of type `real`. | ||
|
||
### Return value | ||
|
||
The result is a `real` array with the same size and kind as `x`. | ||
|
||
If the size of `x` is one, then the sole element of the result is zero. | ||
|
||
### Example | ||
|
||
```fortran | ||
program demo_trapz_weights | ||
use stdlib_experimental_quadrature, only: trapz_weights | ||
implicit none | ||
real :: x(5) = [0., 1., 2., 3., 4.] | ||
real :: y(5) = x**2 | ||
real :: w(5) | ||
w = trapz_weight(x) | ||
print *, dot_product(w, y) | ||
! 22.0 | ||
end program | ||
|
||
``` | ||
|
||
# `simps` - integrate sampled values using Simpson's rule | ||
|
||
Returns the Simpson's rule integral of an array `y` representing discrete samples of a function. The integral is computed assuming either equidistant abscissas with spacing `dx` or arbitary abscissas `x`. | ||
|
||
Simpson's rule is defined for odd-length arrays only. For even-length arrays, an optional argument `even` may be used to specify at which index to replace Simpson's rule with Simpson's 3/8 rule. The 3/8 rule will be used for the array section `y(even:even+4)` and the ordinary Simpon's rule will be used elsewhere. | ||
|
||
### Syntax | ||
|
||
`result = simps(y, x [, even])` | ||
|
||
`result = simps(y, dx [, even])` | ||
|
||
### Arguments | ||
|
||
`y`: Shall be a rank-one array of type `real`. | ||
|
||
`x`: Shall be a rank-one array of type `real` having the same kind and size as `y`. | ||
|
||
`dx`: Shall be a scalar of type `real` having the same kind as `y`. | ||
|
||
`even`: (Optional) Shall be a scalar integer of default kind. Its default value is `1`. | ||
|
||
### Return value | ||
|
||
The result is a scalar of type `real` having the same kind as `y`. | ||
|
||
If the size of `y` is zero or one, the result is zero. | ||
|
||
If the size of `y` is two, the result is the same as if `trapz` had been called instead, regardless of the value of `even`. | ||
|
||
### Example | ||
|
||
TBD | ||
|
||
# `simps_weights` - Simpson's rule weights for given abscissas | ||
|
||
Given an array of abscissas `x`, computes the array of weights `w` such that if `y` represented function values tabulated at `x`, then `sum(w*y)` produces a Simpson's rule approximation to the integral. | ||
|
||
Simpson's rule is defined for odd-length arrays only. For even-length arrays, an optional argument `even` may be used to specify at which index to replace Simpson's rule with Simpson's 3/8 rule. The 3/8 rule will be used for the array section `x(even:even+4)` and the ordinary Simpon's rule will be used elsewhere. | ||
|
||
### Syntax | ||
|
||
`result = simps_weights(x [, even])` | ||
|
||
### Arguments | ||
|
||
`x`: Shall be a rank-one array of type `real`. | ||
|
||
`even`: (Optional) Shall be a scalar integer of default kind. Its default value is `1`. | ||
|
||
### Return value | ||
|
||
The result is a `real` array with the same size and kind as `x`. | ||
|
||
If the size of `x` is one, then the sole element of the result is zero. | ||
|
||
If the size of `x` is two, then the result is the same as if `trapz_weights` had been called instead. | ||
|
||
### Example | ||
|
||
TBD |
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,86 @@ | ||
#:include "common.fypp" | ||
|
||
submodule (stdlib_experimental_quadrature) stdlib_experimental_quadrature_trapz | ||
|
||
implicit none | ||
|
||
contains | ||
|
||
#:for KIND in REAL_KINDS | ||
|
||
pure module function trapz_dx_${KIND}$(y, dx) result(integral) | ||
real(${KIND}$), dimension(:), intent(in) :: y | ||
real(${KIND}$), intent(in) :: dx | ||
real(${KIND}$) :: integral | ||
|
||
integer :: n | ||
|
||
n = size(y) | ||
|
||
select case (n) | ||
case (0:1) | ||
integral = 0.0_${KIND}$ | ||
case (2) | ||
integral = 0.5_${KIND}$*dx*(y(1) + y(2)) | ||
case default | ||
integral = dx*(sum(y(2:n-1)) + 0.5_${KIND}$*(y(1) + y(n))) | ||
end select | ||
end function trapz_dx_${KIND}$ | ||
|
||
#:endfor | ||
#:for KIND in REAL_KINDS | ||
|
||
pure module function trapz_x_${KIND}$(y, x) result(integral) | ||
real(${KIND}$), dimension(:), intent(in) :: y | ||
real(${KIND}$), dimension(size(y)), intent(in) :: x | ||
real(${KIND}$) :: integral | ||
|
||
integer :: i | ||
integer :: n | ||
|
||
n = size(y) | ||
|
||
select case (n) | ||
case (0:1) | ||
integral = 0.0_${KIND}$ | ||
case (2) | ||
integral = 0.5_${KIND}$*(x(2) - x(1))*(y(1) + y(2)) | ||
case default | ||
integral = 0.0_${KIND}$ | ||
do i = 2, n | ||
integral = integral + (x(i) - x(i-1))*(y(i) + y(i-1)) | ||
end do | ||
integral = 0.5_${KIND}$*integral | ||
end select | ||
end function trapz_x_${KIND}$ | ||
|
||
#:endfor | ||
#:for KIND in REAL_KINDS | ||
|
||
pure module function trapz_weights_${KIND}$(x) result(w) | ||
real(${KIND}$), dimension(:), intent(in) :: x | ||
real(${KIND}$), dimension(size(x)) :: w | ||
|
||
integer :: i | ||
integer :: n | ||
|
||
n = size(x) | ||
|
||
select case (n) | ||
case (0) | ||
! no action needed | ||
case (1) | ||
w(1) = 0.0_${KIND}$ | ||
case (2) | ||
w = 0.5_${KIND}$*(x(2) - x(1)) | ||
case default | ||
w(1) = 0.5_${KIND}$*(x(2) - x(1)) | ||
w(n) = 0.5_${KIND}$*(x(n) - x(n-1)) | ||
do i = 2, size(x)-1 | ||
w(i) = 0.5_${KIND}$*(x(i+1) - x(i-1)) | ||
end do | ||
end select | ||
end function trapz_weights_${KIND}$ | ||
|
||
#:endfor | ||
end submodule stdlib_experimental_quadrature_trapz |
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,2 @@ | ||
ADDTEST(trapz) | ||
|
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,3 @@ | ||
PROGS_SRC = test_trapz.f90 | ||
|
||
include ../Makefile.manual.test.mk |
Oops, something went wrong.
Add this suggestion to a batch that can be applied as a single commit.
This suggestion is invalid because no changes were made to the code.
Suggestions cannot be applied while the pull request is closed.
Suggestions cannot be applied while viewing a subset of changes.
Only one suggestion per line can be applied in a batch.
Add this suggestion to a batch that can be applied as a single commit.
Applying suggestions on deleted lines is not supported.
You must change the existing code in this line in order to create a valid suggestion.
Outdated suggestions cannot be applied.
This suggestion has been applied or marked resolved.
Suggestions cannot be applied from pending reviews.
Suggestions cannot be applied on multi-line comments.
Suggestions cannot be applied while the pull request is queued to merge.
Suggestion cannot be applied right now. Please check back later.
Uh oh!
There was an error while loading. Please reload this page.