Skip to content

Commit 0d0cb44

Browse files
authored
Merge pull request #70 from GComitini/mainline-complex-integration
IMPL: generic and complex integration (#35)
2 parents eede4db + 3c28e69 commit 0d0cb44

File tree

3 files changed

+178
-45
lines changed

3 files changed

+178
-45
lines changed

src/complex/integral.rs

+53
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,53 @@
1+
use crate::complex::C64;
2+
use crate::numerical::integral::{GKIntegrable, GLKIntegrable, NCIntegrable};
3+
use crate::structure::polynomial::{lagrange_polynomial, Calculus, Polynomial};
4+
5+
// Newton Cotes Quadrature for Complex Functions of one Real Variable
6+
impl NCIntegrable for C64 {
7+
type NodeY = (Vec<f64>, Vec<f64>);
8+
type NCPolynomial = (Polynomial, Polynomial);
9+
10+
fn compute_node_y<F>(f: F, node_x: &[f64]) -> Self::NodeY
11+
where
12+
F: Fn(f64) -> Self,
13+
{
14+
node_x
15+
.iter()
16+
.map(|x| {
17+
let z = f(*x);
18+
(z.re, z.im)
19+
})
20+
.unzip()
21+
}
22+
23+
fn compute_polynomial(node_x: &[f64], node_y: &Self::NodeY) -> Self::NCPolynomial {
24+
(
25+
lagrange_polynomial(node_x.to_vec(), node_y.0.to_vec()),
26+
lagrange_polynomial(node_x.to_vec(), node_y.1.to_vec()),
27+
)
28+
}
29+
30+
fn integrate_polynomial(p: &Self::NCPolynomial) -> Self::NCPolynomial {
31+
(p.0.integral(), p.1.integral())
32+
}
33+
34+
fn evaluate_polynomial(p: &Self::NCPolynomial, x: f64) -> Self {
35+
p.0.eval(x) + C64::I * p.1.eval(x)
36+
}
37+
}
38+
39+
// Gauss Lagrange and Kronrod Quadrature for Complex Functions of one Real Variable
40+
impl GLKIntegrable for C64 {
41+
const ZERO: Self = C64::ZERO;
42+
}
43+
44+
// Gauss Kronrod Quadrature for Complex Functions of one Real Variable
45+
impl GKIntegrable for C64 {
46+
fn is_finite(&self) -> bool {
47+
C64::is_finite(*self)
48+
}
49+
50+
fn gk_norm(&self) -> f64 {
51+
self.norm()
52+
}
53+
}

src/complex/mod.rs

+1
Original file line numberDiff line numberDiff line change
@@ -2,5 +2,6 @@ use num_complex::Complex;
22

33
pub type C64 = Complex<f64>;
44

5+
pub mod integral;
56
pub mod matrix;
67
pub mod vector;

src/numerical/integral.rs

+124-45
Original file line numberDiff line numberDiff line change
@@ -1,4 +1,4 @@
1-
use crate::structure::polynomial::{lagrange_polynomial, Calculus};
1+
use crate::structure::polynomial::{lagrange_polynomial, Calculus, Polynomial};
22
use crate::traits::fp::FPVector;
33
use crate::util::non_macro::seq;
44

@@ -132,6 +132,78 @@ impl Integral {
132132
}
133133
}
134134

135+
// Types that can be integrated using Newton Cotes Quadrature
136+
pub trait NCIntegrable: std::ops::Sub<Self, Output = Self> + Sized {
137+
type NodeY;
138+
type NCPolynomial;
139+
140+
fn compute_node_y<F>(f: F, node_x: &[f64]) -> Self::NodeY
141+
where
142+
F: Fn(f64) -> Self;
143+
144+
fn compute_polynomial(node_x: &[f64], node_y: &Self::NodeY) -> Self::NCPolynomial;
145+
146+
fn integrate_polynomial(p: &Self::NCPolynomial) -> Self::NCPolynomial;
147+
148+
fn evaluate_polynomial(p: &Self::NCPolynomial, x: f64) -> Self;
149+
}
150+
151+
impl NCIntegrable for f64 {
152+
type NodeY = Vec<f64>;
153+
type NCPolynomial = Polynomial;
154+
155+
fn compute_node_y<F>(f: F, node_x: &[f64]) -> Vec<f64>
156+
where
157+
F: Fn(f64) -> Self,
158+
{
159+
node_x.to_vec().fmap(|x| f(x))
160+
}
161+
162+
fn compute_polynomial(node_x: &[f64], node_y: &Vec<f64>) -> Polynomial {
163+
lagrange_polynomial(node_x.to_vec(), node_y.to_vec())
164+
}
165+
166+
fn integrate_polynomial(p: &Polynomial) -> Polynomial {
167+
p.integral()
168+
}
169+
170+
fn evaluate_polynomial(p: &Polynomial, x: f64) -> Self {
171+
p.eval(x)
172+
}
173+
}
174+
175+
// Types that can be integrated using Gauss Lagrange or Kronrod Quadrature
176+
pub trait GLKIntegrable:
177+
std::ops::Add<Self, Output = Self> + std::ops::Mul<f64, Output = Self> + Sized
178+
{
179+
const ZERO: Self;
180+
}
181+
182+
impl GLKIntegrable for f64 {
183+
const ZERO: Self = 0f64;
184+
}
185+
186+
// Types that can be integrated using Gauss Kronrod Quadrature
187+
pub trait GKIntegrable: GLKIntegrable + std::ops::Sub<Self, Output = Self> + Sized + Clone {
188+
fn is_finite(&self) -> bool;
189+
190+
fn gk_norm(&self) -> f64;
191+
192+
fn is_within_tol(&self, tol: f64) -> bool {
193+
self.gk_norm() < tol
194+
}
195+
}
196+
197+
impl GKIntegrable for f64 {
198+
fn is_finite(&self) -> bool {
199+
f64::is_finite(*self)
200+
}
201+
202+
fn gk_norm(&self) -> f64 {
203+
self.abs()
204+
}
205+
}
206+
135207
/// Numerical Integration
136208
///
137209
/// # Description
@@ -159,29 +231,31 @@ impl Integral {
159231
/// * `G20K41R`
160232
/// * `G25K51R`
161233
/// * `G30K61R`
162-
pub fn integrate<F>(f: F, (a, b): (f64, f64), method: Integral) -> f64
234+
pub fn integrate<F, Y>(f: F, (a, b): (f64, f64), method: Integral) -> Y
163235
where
164-
F: Fn(f64) -> f64 + Copy,
236+
F: Fn(f64) -> Y + Copy,
237+
Y: GKIntegrable + NCIntegrable,
165238
{
166239
match method {
167240
Integral::GaussLegendre(n) => gauss_legendre_quadrature(f, n, (a, b)),
168241
Integral::NewtonCotes(n) => newton_cotes_quadrature(f, n, (a, b)),
169-
method => gauss_kronrod_quadrature(f, (a,b), method),
242+
method => gauss_kronrod_quadrature(f, (a, b), method),
170243
}
171244
}
172245

173246
/// Newton Cotes Quadrature
174-
pub fn newton_cotes_quadrature<F>(f: F, n: usize, (a, b): (f64, f64)) -> f64
247+
pub fn newton_cotes_quadrature<F, Y>(f: F, n: usize, (a, b): (f64, f64)) -> Y
175248
where
176-
F: Fn(f64) -> f64,
249+
F: Fn(f64) -> Y,
250+
Y: NCIntegrable,
177251
{
178252
let h = (b - a) / (n as f64);
179253
let node_x = seq(a, b, h);
180-
let node_y = node_x.fmap(|x| f(x));
181-
let p = lagrange_polynomial(node_x, node_y);
182-
let q = p.integral();
254+
let node_y = Y::compute_node_y(f, &node_x);
255+
let p = Y::compute_polynomial(&node_x, &node_y);
256+
let q = Y::integrate_polynomial(&p);
183257

184-
q.eval(b) - q.eval(a)
258+
Y::evaluate_polynomial(&q, b) - Y::evaluate_polynomial(&q, a)
185259
}
186260

187261
/// Gauss Legendre Quadrature
@@ -195,11 +269,12 @@ where
195269
/// # Reference
196270
/// * A. N. Lowan et al. (1942)
197271
/// * [Keisan Online Calculator](https://keisan.casio.com/exec/system/1329114617)
198-
pub fn gauss_legendre_quadrature<F>(f: F, n: usize, (a, b): (f64, f64)) -> f64
272+
pub fn gauss_legendre_quadrature<F, Y>(f: F, n: usize, (a, b): (f64, f64)) -> Y
199273
where
200-
F: Fn(f64) -> f64,
274+
F: Fn(f64) -> Y,
275+
Y: GLKIntegrable,
201276
{
202-
(b - a) / 2f64 * unit_gauss_legendre_quadrature(|x| f(x * (b - a) / 2f64 + (a + b) / 2f64), n)
277+
unit_gauss_legendre_quadrature(|x| f(x * (b - a) / 2f64 + (a + b) / 2f64), n) * ((b - a) / 2f64)
203278
}
204279

205280
/// Gauss Kronrod Quadrature
@@ -215,16 +290,17 @@ where
215290
/// * arXiv: [1003.4629](https://arxiv.org/abs/1003.4629)
216291
/// * [Keisan Online Calculator](https://keisan.casio.com/exec/system/1329114617)
217292
#[allow(non_snake_case)]
218-
pub fn gauss_kronrod_quadrature<F, T, S>(f: F, (a, b): (T, S), method: Integral) -> f64
293+
pub fn gauss_kronrod_quadrature<F, T, S, Y>(f: F, (a, b): (T, S), method: Integral) -> Y
219294
where
220-
F: Fn(f64) -> f64 + Copy,
221-
T: Into<f64>,
222-
S: Into<f64>,
295+
F: Fn(f64) -> Y + Copy,
296+
T: Into<f64>,
297+
S: Into<f64>,
298+
Y: GKIntegrable,
223299
{
224300
let (g, k) = method.get_gauss_kronrod_order();
225301
let tol = method.get_tol();
226302
let max_iter = method.get_max_iter();
227-
let mut I = 0f64;
303+
let mut I = Y::ZERO;
228304
let mut S: Vec<(f64, f64, f64, u32)> = vec![];
229305
S.push((a.into(), b.into(), tol, max_iter));
230306

@@ -235,15 +311,15 @@ where
235311
let K = kronrod_quadrature(f, k as usize, (a, b));
236312
let c = (a + b) / 2f64;
237313
let tol_curr = if method.is_relative() {
238-
tol * G
314+
tol * G.gk_norm()
239315
} else {
240316
tol
241317
};
242-
if (G - K).abs() < tol_curr || a == b || max_iter == 0 {
243-
if ! G.is_finite() {
318+
if (G.clone() - K).is_within_tol(tol_curr) || a == b || max_iter == 0 {
319+
if !G.is_finite() {
244320
return G;
245321
}
246-
I += G;
322+
I = I + G;
247323
} else {
248324
S.push((a, c, tol / 2f64, max_iter - 1));
249325
S.push((c, b, tol / 2f64, max_iter - 1));
@@ -255,24 +331,26 @@ where
255331
I
256332
}
257333

258-
pub fn kronrod_quadrature<F>(f: F, n: usize, (a, b): (f64, f64)) -> f64
334+
pub fn kronrod_quadrature<F, Y>(f: F, n: usize, (a, b): (f64, f64)) -> Y
259335
where
260-
F: Fn(f64) -> f64,
336+
F: Fn(f64) -> Y,
337+
Y: GLKIntegrable,
261338
{
262-
(b - a) / 2f64 * unit_kronrod_quadrature(|x| f(x * (b-a) / 2f64 + (a + b) / 2f64), n)
339+
unit_kronrod_quadrature(|x| f(x * (b - a) / 2f64 + (a + b) / 2f64), n) * ((b - a) / 2f64)
263340
}
264341

265342
// =============================================================================
266343
// Gauss Legendre Backends
267344
// =============================================================================
268-
fn unit_gauss_legendre_quadrature<F>(f: F, n: usize) -> f64
345+
fn unit_gauss_legendre_quadrature<F, Y>(f: F, n: usize) -> Y
269346
where
270-
F: Fn(f64) -> f64,
347+
F: Fn(f64) -> Y,
348+
Y: GLKIntegrable,
271349
{
272350
let (a, x) = gauss_legendre_table(n);
273-
let mut s = 0f64;
351+
let mut s = Y::ZERO;
274352
for i in 0..a.len() {
275-
s += a[i] * f(x[i]);
353+
s = s + f(x[i]) * a[i];
276354
}
277355
s
278356
}
@@ -378,14 +456,15 @@ fn gauss_legendre_table(n: usize) -> (Vec<f64>, Vec<f64>) {
378456
// Gauss Kronrod Backends
379457
// =============================================================================
380458

381-
fn unit_kronrod_quadrature<F>(f: F, n: usize) -> f64
459+
fn unit_kronrod_quadrature<F, Y>(f: F, n: usize) -> Y
382460
where
383-
F: Fn(f64) -> f64,
461+
F: Fn(f64) -> Y,
462+
Y: GLKIntegrable,
384463
{
385-
let (a,x) = kronrod_table(n);
386-
let mut s = 0f64;
387-
for i in 0 .. a.len() {
388-
s += a[i] * f(x[i]);
464+
let (a, x) = kronrod_table(n);
465+
let mut s = Y::ZERO;
466+
for i in 0..a.len() {
467+
s = s + f(x[i]) * a[i];
389468
}
390469
s
391470
}
@@ -414,28 +493,28 @@ fn kronrod_table(n: usize) -> (Vec<f64>, Vec<f64>) {
414493

415494
match n % 2 {
416495
0 => {
417-
for i in 0 .. ref_node.len() {
496+
for i in 0..ref_node.len() {
418497
result_node[i] = ref_node[i];
419498
result_weight[i] = ref_weight[i];
420499
}
421500

422-
for i in ref_node.len() .. n {
423-
result_node[i] = -ref_node[n-i-1];
424-
result_weight[i] = ref_weight[n-i-1];
501+
for i in ref_node.len()..n {
502+
result_node[i] = -ref_node[n - i - 1];
503+
result_weight[i] = ref_weight[n - i - 1];
425504
}
426505
}
427506
1 => {
428-
for i in 0 .. ref_node.len() {
507+
for i in 0..ref_node.len() {
429508
result_node[i] = ref_node[i];
430509
result_weight[i] = ref_weight[i];
431510
}
432511

433-
for i in ref_node.len() .. n {
434-
result_node[i] = -ref_node[n-i];
435-
result_weight[i] = ref_weight[n-i];
512+
for i in ref_node.len()..n {
513+
result_node[i] = -ref_node[n - i];
514+
result_weight[i] = ref_weight[n - i];
436515
}
437516
}
438-
_ => unreachable!()
517+
_ => unreachable!(),
439518
}
440519
(result_weight, result_node)
441520
}
@@ -922,7 +1001,7 @@ const LEGENDRE_WEIGHT_25: [f64; 13] = [
9221001
0.054904695975835192,
9231002
0.040939156701306313,
9241003
0.026354986615032137,
925-
0.011393798501026288,
1004+
0.011393798501026288,
9261005
];
9271006
const LEGENDRE_WEIGHT_26: [f64; 13] = [
9281007
0.118321415279262277,

0 commit comments

Comments
 (0)