-
Notifications
You must be signed in to change notification settings - Fork 11
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
Numerical Integrator added to SRC along with the __init__() file #391
Merged
Merged
Changes from all commits
Commits
Show all changes
10 commits
Select commit
Hold shift + click to select a range
eb21002
added numerical integrator and libray call to __innit
Witt-D bbe94ce
tidied up linting
Witt-D e8e5ae2
fixed lint
Witt-D 855c10e
added numerical integrator test
Witt-D 1284f36
typo
Witt-D 7c05364
adjusted docstring
Witt-D 4c873b2
adjusted sine test
Witt-D 6af56ff
put in alphabetical order
Witt-D ec18dbe
fixed lint
Witt-D e8ec956
err tol too low
Witt-D 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 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 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,63 @@ | ||
import numpy as np | ||
|
||
|
||
class NumericalIntegral(object): | ||
""" | ||
A class for numerically evaluating and tabulating some 1D integral. | ||
Args: | ||
lower_bound(float): lower bound of integral | ||
upper_bound(float): upper bound of integral | ||
num_points(float): number of points to tabulate integral at | ||
|
||
""" | ||
def __init__(self, lower_bound, upper_bound, num_points=500): | ||
|
||
# if upper_bound <= lower_bound: | ||
# raise ValueError('lower_bound must be lower than upper_bound') | ||
self.x = np.linspace(lower_bound, upper_bound, num_points) | ||
self.x_double = np.linspace(lower_bound, upper_bound, 2*num_points-1) | ||
self.lower_bound = lower_bound | ||
self.upper_bound = upper_bound | ||
self.num_points = num_points | ||
self.tabulated = False | ||
|
||
def tabulate(self, expression): | ||
""" | ||
Tabulate some integral expression using Simpson's rule. | ||
Args: | ||
expression (func): a function representing the integrand to be | ||
evaluated. should take a numpy array as an argument. | ||
""" | ||
|
||
self.cumulative = np.zeros_like(self.x) | ||
self.interval_areas = np.zeros(len(self.x)-1) | ||
# Evaluate expression in advance to make use of numpy optimisation | ||
# We evaluate at the tabulation points and the midpoints of the intervals | ||
f = expression(self.x_double) | ||
|
||
# Just do Simpson's rule for evaluating area of each interval | ||
self.interval_areas = ((self.x[1:] - self.x[:-1]) / 6.0 | ||
* (f[2::2] + 4.0 * f[1::2] + f[:-1:2])) | ||
|
||
# Add the interval areas together to create cumulative integral | ||
for i in range(self.num_points - 1): | ||
self.cumulative[i+1] = self.cumulative[i] + self.interval_areas[i] | ||
|
||
self.tabulated = True | ||
|
||
def evaluate_at(self, points): | ||
""" | ||
Evaluates the integral at some point using linear interpolation. | ||
Args: | ||
points (float or iter) the point value, or array of point values to | ||
evaluate the integral at. | ||
Return: | ||
returns the numerical approximation of the integral from lower | ||
bound to point(s) | ||
""" | ||
# Do linear interpolation from tabulated values | ||
if not self.tabulated: | ||
raise RuntimeError( | ||
'Integral must be tabulated before we can evaluate it at a point') | ||
|
||
return np.interp(points, self.x, self.cumulative) |
This file contains 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,34 @@ | ||
""" | ||
Tests the numerical integrator. | ||
""" | ||
from gusto import NumericalIntegral | ||
from numpy import sin, pi | ||
import pytest | ||
|
||
|
||
def quadratic(x): | ||
return x**2 | ||
|
||
|
||
def sine(x): | ||
return sin(x) | ||
|
||
|
||
@pytest.mark.parametrize("integrand_name", ["quadratic", "sine"]) | ||
def test_numerical_integrator(integrand_name): | ||
if integrand_name == "quadratic": | ||
integrand = quadratic | ||
upperbound = 3 | ||
answer = 9 | ||
elif integrand_name == "sine": | ||
integrand = sine | ||
upperbound = pi | ||
answer = 2 | ||
else: | ||
raise ValueError(f'{integrand_name} integrand not recognised') | ||
numerical_integral = NumericalIntegral(0, upperbound) | ||
numerical_integral.tabulate(integrand) | ||
area = numerical_integral.evaluate_at(upperbound) | ||
err_tol = 1e-10 | ||
assert abs(area-answer) < err_tol, \ | ||
f'numerical integrator is incorrect for {integrand_name} function' |
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.
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
This docstring should also have a
Returns:
statement, as the function routines something