Skip to content
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

Deprecate applying Gauss-Kronrod rules to surfaces #136

Merged
merged 14 commits into from
Nov 28, 2024
Merged

Deprecate applying Gauss-Kronrod rules to surfaces #136

merged 14 commits into from
Nov 28, 2024

Conversation

mikeingold
Copy link
Collaborator

@mikeingold mikeingold commented Nov 27, 2024

Changes

  • Deprecates support for integrating (non-specialization) surfaces using GaussKronrod rules, issuing a Base.depwarn deprecation warning.
  • Eliminates the _integral_gk_*d helper functions - consolidated into _integral(f, g, ::GaussKronrod) function instead.
  • Minor tweaks
    • Removed explicit FP conversions inside GaussLegendre method, relying instead on automatic promotion from rational to float, resulting in a >2x performance improvement.
    • Removed obsolete QuadGK direct dependency from tests.

Discussion

Surfaces are geometries where paramdim(geometry) == 2. Use of nested GaussKronrod rules for these geometries is strictly inferior to simply using an equivalent HAdaptiveCubature rule, where generally returns results faster and with more accuracy. Earlier in the development of this package, using GaussKronrod rules on surfaces provided a workaround since HAdaptiveCubature rules didn't yet support Unitful-valued integrand functions. Since that issue has been resolved, I think we should mark this usage as deprecated in v0.16 and plan to remove support in v0.17. Support for usage on several of the specialization geometries would be unaffected for now, where there remains some justification for their usage.

@mikeingold mikeingold changed the title Deprecatie applying Gauss-Kronrod rules to surfaces Deprecate applying Gauss-Kronrod rules to surfaces Nov 28, 2024
src/integral.jl Outdated Show resolved Hide resolved
@mikeingold
Copy link
Collaborator Author

@JoshuaLampert can you sanity check the _integral(f, geometry, ::GaussKronrod) function for me? For some reason the integrand symbol isn't being recognized, but I just don't see the issue. I tried commenting out the preceding depwarn in case there were any weird side effects, but saw no change.

Copy link
Member

@JoshuaLampert JoshuaLampert left a comment

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Looks mostly good to me. Please find some comments below.

src/integral.jl Outdated Show resolved Hide resolved
src/integral.jl Outdated Show resolved Hide resolved
src/integral.jl Outdated Show resolved Hide resolved
src/integral.jl Outdated Show resolved Hide resolved
docs/src/supportmatrix.md Outdated Show resolved Hide resolved
docs/src/supportmatrix.md Outdated Show resolved Hide resolved
mikeingold and others added 2 commits November 28, 2024 08:54
Co-authored-by: Joshua Lampert <[email protected]>
@mikeingold
Copy link
Collaborator Author

Thanks!

Copy link

codecov bot commented Nov 28, 2024

Codecov Report

All modified and coverable lines are covered by tests ✅

Project coverage is 95.86%. Comparing base (0cd335d) to head (f1e7d46).
Report is 1 commits behind head on main.

Additional details and impacted files
@@            Coverage Diff             @@
##             main     #136      +/-   ##
==========================================
- Coverage   95.90%   95.86%   -0.05%     
==========================================
  Files          17       17              
  Lines         293      290       -3     
==========================================
- Hits          281      278       -3     
  Misses         12       12              

☔ View full report in Codecov by Sentry.
📢 Have feedback on the report? Share it here.

Copy link
Contributor

Benchmark Results

main f1e7d46... main/f1e7d46fae4de3...
Differentials/Differential 0.239 ± 0.0021 μs 0.24 ± 0.001 μs 0.996
Differentials/Jacobian 0.205 ± 0.001 μs 0.205 ± 0.001 μs 1
Integrals/Meshes.BezierCurve/Scalar GaussKronrod 22.4 ± 0.04 ms 22.4 ± 0.04 ms 1
Integrals/Meshes.BezierCurve/Scalar GaussLegendre 21.4 ± 0.043 ms 21.4 ± 0.04 ms 1
Integrals/Meshes.BezierCurve/Scalar HAdaptiveCubature 22.6 ± 0.042 ms 22.6 ± 0.038 ms 1
Integrals/Meshes.BezierCurve/Vector GaussKronrod 22.4 ± 0.043 ms 22.4 ± 0.037 ms 1
Integrals/Meshes.BezierCurve/Vector GaussLegendre 21.5 ± 0.047 ms 21.4 ± 0.03 ms 1
Integrals/Meshes.BezierCurve/Vector HAdaptiveCubature 22.6 ± 0.047 ms 22.6 ± 0.053 ms 1
Integrals/Meshes.Segment/Scalar GaussKronrod 1.19 ± 0.0051 μs 1.2 ± 0.007 μs 0.992
Integrals/Meshes.Segment/Scalar GaussLegendre 0.0535 ± 0.00064 ms 8.5 ± 0.2 μs 6.3
Integrals/Meshes.Segment/Scalar HAdaptiveCubature 1.7 ± 0.045 μs 1.73 ± 0.045 μs 0.983
Integrals/Meshes.Segment/Vector GaussKronrod 3.38 ± 0.054 μs 3.41 ± 0.06 μs 0.992
Integrals/Meshes.Segment/Vector GaussLegendre 0.0638 ± 0.00086 ms 24.4 ± 0.61 μs 2.62
Integrals/Meshes.Segment/Vector HAdaptiveCubature 4.69 ± 0.09 μs 4.58 ± 0.077 μs 1.02
Integrals/Meshes.Sphere/Scalar GaussKronrod 0.0904 ± 0.00035 ms 0.0839 ± 0.00053 ms 1.08
Integrals/Meshes.Sphere/Scalar GaussLegendre 12.2 ± 0.095 ms 2.25 ± 0.0085 ms 5.42
Integrals/Meshes.Sphere/Scalar HAdaptiveCubature 0.0601 ± 0.00022 ms 0.0604 ± 0.00025 ms 0.996
Integrals/Meshes.Sphere/Vector GaussKronrod 0.117 ± 0.0012 ms 0.117 ± 0.00085 ms 0.999
Integrals/Meshes.Sphere/Vector GaussLegendre 13.5 ± 0.52 ms 3.61 ± 0.074 ms 3.75
Integrals/Meshes.Sphere/Vector HAdaptiveCubature 0.11 ± 0.0011 ms 0.11 ± 0.0011 ms 0.999
time_to_load 1.5 ± 0.011 s 1.53 ± 0.013 s 0.98

Benchmark Plots

A plot of the benchmark results have been uploaded as an artifact to the workflow run for this PR.
Go to "Actions"->"Benchmark a pull request"->[the most recent run]->"Artifacts" (at the bottom).

@mikeingold mikeingold marked this pull request as ready for review November 28, 2024 14:29
@mikeingold
Copy link
Collaborator Author

Good to go here. Looks like the tweaks to GaussLegendre methods achieved a surprisingly big performance improvement.

Copy link
Member

@JoshuaLampert JoshuaLampert left a comment

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Good to go here. Looks like the tweaks to GaussLegendre methods achieved a surprisingly big performance improvement.

That's nice! Thanks!

@mikeingold mikeingold merged commit cba39f0 into main Nov 28, 2024
12 checks passed
@mikeingold mikeingold deleted the tweaks branch November 28, 2024 16:35
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
Projects
None yet
Development

Successfully merging this pull request may close these issues.

2 participants