Skip to content
This repository was archived by the owner on May 15, 2025. It is now read-only.

Fix Issue #34: Add new method Alefeld, Potra, Shi (1995) #56

Merged
merged 19 commits into from
Mar 30, 2023

Conversation

huiyuxie
Copy link
Contributor

Summary

This pull request addresses issue #34 by adding a new method Alefeld to the SimpleNonlinearSolve.

Changes

  • Implemented the new method Alefeld(), a subtype of AbstractBracketingAlgorithm, in alefeld.jl
  • Included alefeld.jl in SimpleNonlinearSolve.jl and exported Alefeld
  • Added simple scalar function test like other bracketing methods (Falsi, Ridder, Brent) in basictests.jl

Remaining Bugs

When tried to run runtests.jl, the terminal reported ''UndefVarError: Alefeld not defined'' checked every detail about this repository but still could not fix it.
Also, tried @ChrisRackauckas 's suggestion ]dev SimpleNonlinearSolve and ]st the terminal showed ''Status ~/.julia/environments/JuliaDev/Project.toml
[6e4b80f9] BenchmarkTools v1.3.2
[2b5f629d] DiffEqBase v6.122.1
[f6369f11] ForwardDiff v0.10.35
[872c559c] NNlib v0.8.19
[1bc83da4] SafeTestsets v0.0.1
[727e6d20] SimpleNonlinearSolve v0.1.14 ~/.julia/dev/SimpleNonlinearSolve
[90137ffa] StaticArrays v1.5.19''

Next Steps

Add more complicated test to Alefeld() in basictests.jl and optimize the algorithm implementation.

Feedback

Any additional suggestions to help resolve the remaining bugs would be highly appreciated. Please let me know if there is anything I can do to improve the code or address the remaining bugs.

Project.toml Outdated
Comment on lines 8 to 19
BenchmarkTools = "6e4b80f9-dd63-53aa-95a3-0cdb28fa8baf"
DiffEqBase = "2b5f629d-d688-5b77-993f-72d75c75574e"
FiniteDiff = "6a86dc24-6348-571c-b903-95158fe2bd41"
ForwardDiff = "f6369f11-7733-5829-9624-2563aa707210"
LinearAlgebra = "37e2e46d-f89d-539d-b4ee-838fcccc9c8e"
NNlib = "872c559c-99b0-510c-b3b7-b6c96a88d5cd"
Reexport = "189a3867-3050-52da-a836-e630ba90ab69"
Requires = "ae029012-a4dd-5104-9daa-d747884805df"
SafeTestsets = "1bc83da4-3b8d-516f-aca4-4fe02f6d838f"
SciMLBase = "0bca4576-84f4-4d90-8ffe-ffa030f20462"
SnoopPrecompile = "66db9d55-30c0-4569-8b51-7e840670fc0c"
StaticArrays = "90137ffa-7385-5640-81b9-e52037218182"
Copy link
Member

Choose a reason for hiding this comment

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

These shouldn't be added as deps. They are only test deps.

huiyuxie and others added 5 commits March 29, 2023 09:14
Co-authored-by: Christopher Rackauckas <[email protected]>
Co-authored-by: Christopher Rackauckas <[email protected]>
Co-authored-by: Christopher Rackauckas <[email protected]>
@codecov
Copy link

codecov bot commented Mar 29, 2023

Codecov Report

Merging #56 (e44bab2) into main (9088070) will decrease coverage by 0.95%.
The diff coverage is 84.09%.

@@            Coverage Diff             @@
##             main      #56      +/-   ##
==========================================
- Coverage   92.46%   91.51%   -0.95%     
==========================================
  Files          15       16       +1     
  Lines         690      778      +88     
==========================================
+ Hits          638      712      +74     
- Misses         52       66      +14     
Impacted Files Coverage Δ
src/SimpleNonlinearSolve.jl 50.00% <ø> (ø)
src/alefeld.jl 84.09% <84.09%> (ø)

📣 We’re building smart automated test selection to slash your CI/CD build times. Learn more

@huiyuxie
Copy link
Contributor Author

huiyuxie commented Mar 30, 2023

I suspect that there might be something wrong with a small step in this paper. I debugged the failed cases and found that the errors came from Newton-Quadratic, one of the subroutines in this paper.

In Algorithm 4.2, Newton-Quadratic provides an interval [a, b] and an approximation c to bracket to approximate c again, which should ensure that c is within the interval (a, b). However, in some test cases, Newton-Quadratic returns c outside the interval (a, b). I further examined this subroutine function and found that on page 330, equation (4), the induction is incorrect.
Screenshot 2023-03-29 at 8 36 23 PM
I could not get Br_{i-1} in the second line from P(r_{i-1}) in the first line. According to the definition of P(x), P(r_{i-1}) should be f(a)+B(r_{i-1}-a) + A(r_{i-1}-a)(r_{i-1}-b). This equation even does not hold for base case where r_{0} equals a or b.

However, this mistake seems so straightforward and wired to me that I am unsure whether it is an actual error or if I have misunderstood something. I would appreciate it if someone @ranocha could help me double-check this equation. Thank you!

@huiyuxie
Copy link
Contributor Author

If the induction is indeed incorrect, I will directly use the formula in the first line of equation (4) rather than the second line.

@ChrisRackauckas
Copy link
Member

Co-authored-by: Christopher Rackauckas <[email protected]>
@ChrisRackauckas
Copy link
Member

I think you are correct. At first I was wondering if something factored out and canceled but no from the looks of the denominator that's not the case. Just looks like there was a typo.

@huiyuxie
Copy link
Contributor Author

huiyuxie commented Mar 30, 2023

Thank you so much for your help @ChrisRackauckas!

Here I have fixed two bugs in alefeld.jl:

  1. As I mentioned earlier, there seem to be some errors in the Newton-Quadratic function induction from the paper. I have now directly used the formula from the first line of equation (4) on page 330.
  2. The error "a\bar not defined" arose from the floating-point limit during the computation process, so I added a section to check whether c is equal to either a or b and constructed the solution with "retcode = ReturnCode.FloatingPointLimit".

Now the majority of test cases have passed, but some test cases still fail. I believe there are still some corner cases that the code cannot handle very well. I am in the process of further debugging.
Screenshot 2023-03-30 at 2 06 47 PM

@ChrisRackauckas
Copy link
Member

Tests seem to have passed. Looks good!

@ChrisRackauckas ChrisRackauckas merged commit 8eb8afc into SciML:main Mar 30, 2023
Sign up for free to subscribe to this conversation on GitHub. Already have an account? Sign in.
Labels
None yet
Projects
None yet
Development

Successfully merging this pull request may close these issues.

2 participants