Skip to content

Commit b1fb265

Browse files
committed
Add initial spec for trapz and trapz_weights
1 parent 52386a6 commit b1fb265

File tree

1 file changed

+79
-0
lines changed

1 file changed

+79
-0
lines changed

src/stdlib_experimental_quadrature.md

Lines changed: 79 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,79 @@
1+
# Numerical integration
2+
3+
## Implemented
4+
5+
* `trapz`
6+
* `trapz_weights`
7+
8+
## `trapz` - integrate sampled values using trapezoidal rule
9+
10+
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`.
11+
12+
### Syntax
13+
14+
`result = trapz(y, x)`
15+
16+
`result = trapz(y, dx)`
17+
18+
### Arguments
19+
20+
`y`: Shall be a rank-one array of type `real`.
21+
22+
`x`: Shall be a rank-one array of type `real` having the same kind and size as `y`.
23+
24+
`dx`: Shall be a scalar of type `real` having the same kind as `y`.
25+
26+
### Return value
27+
28+
The result is a scalar of type `real` having the same kind as `y`.
29+
30+
If the size of `y` is zero or one, the result is zero.
31+
32+
### Example
33+
34+
```fortran
35+
program demo_trapz
36+
use stdlib_experimental_quadrature, only: trapz
37+
implicit none
38+
real :: x(5) = [0., 1., 2., 3., 4.]
39+
real :: y(5) = x**2
40+
print *, trapz(y, x)
41+
! 22.0
42+
print *, trapz(y, 0.5)
43+
! 11.0
44+
end program
45+
```
46+
47+
## `trapz_weights` - compute trapezoidal rule weights for given abscissas
48+
49+
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.
50+
51+
### Syntax
52+
53+
`result = trapz_weights(x)`
54+
55+
### Arguments
56+
57+
`x`: Shall be an array of type `real`.
58+
59+
### Return value
60+
61+
The result is a `real` array with the same size and kind as `x`.
62+
63+
If the size of `x` is one, then the only element of the result is zero.
64+
65+
### Example
66+
67+
```fortran
68+
program demo_trapz_weights
69+
use stdlib_experimental_quadrature, only: trapz_weights
70+
implicit none
71+
real :: x(5) = [0., 1., 2., 3., 4.]
72+
real :: y(5) = x**2
73+
real :: w(5)
74+
w = trapz_weight(x)
75+
print *, dot_product(w, y)
76+
! 22.0
77+
end program
78+
79+
```

0 commit comments

Comments
 (0)