-
-
Couldn't load subscription status.
- Fork 56
feat: add PETScSNES
#482
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
Merged
Merged
feat: add PETScSNES
#482
Changes from all commits
Commits
Show all changes
6 commits
Select commit
Hold shift + click to select a range
8fe5ee6
feat: add `PETScSNES`
avik-pal 12baf76
feat: support automatic sparsity detection for PETSc
avik-pal 9e63995
test: add PETScSNES to the wrapper tests
avik-pal 370aa68
docs: add PETSc example
avik-pal 05bbc79
test: skip PETSc tests on windows
avik-pal c143b5c
docs: print the benchmark results
avik-pal 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 hidden or 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 hidden or 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 hidden or 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 hidden or 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,17 @@ | ||
| # PETSc.jl | ||
|
|
||
| This is a extension for importing solvers from PETSc.jl SNES into the SciML interface. Note | ||
| that these solvers do not come by default, and thus one needs to install the package before | ||
| using these solvers: | ||
|
|
||
| ```julia | ||
| using Pkg | ||
| Pkg.add("PETSc") | ||
| using PETSc, NonlinearSolve | ||
| ``` | ||
|
|
||
| ## Solver API | ||
|
|
||
| ```@docs | ||
| PETScSNES | ||
| ``` | ||
This file contains hidden or 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 hidden or 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,81 @@ | ||
| # [PETSc SNES Example 2](@id snes_ex2) | ||
|
|
||
| This implements `src/snes/examples/tutorials/ex2.c` from PETSc and `examples/SNES_ex2.jl` | ||
| from PETSc.jl using automatic sparsity detection and automatic differentiation using | ||
| `NonlinearSolve.jl`. | ||
|
|
||
| This solves the equations sequentially. Newton method to solve | ||
| `u'' + u^{2} = f`, sequentially. | ||
|
|
||
| ```@example snes_ex2 | ||
| using NonlinearSolve, PETSc, LinearAlgebra, SparseConnectivityTracer, BenchmarkTools | ||
|
|
||
| u0 = fill(0.5, 128) | ||
|
|
||
| function form_residual!(resid, x, _) | ||
| n = length(x) | ||
| xp = LinRange(0.0, 1.0, n) | ||
| F = 6xp .+ (xp .+ 1e-12) .^ 6 | ||
|
|
||
| dx = 1 / (n - 1) | ||
| resid[1] = x[1] | ||
| for i in 2:(n - 1) | ||
| resid[i] = (x[i - 1] - 2x[i] + x[i + 1]) / dx^2 + x[i] * x[i] - F[i] | ||
| end | ||
| resid[n] = x[n] - 1 | ||
|
|
||
| return | ||
| end | ||
| ``` | ||
|
|
||
| To use automatic sparsity detection, we need to specify `sparsity` keyword argument to | ||
| `NonlinearFunction`. See [Automatic Sparsity Detection](@ref sparsity-detection) for more | ||
| details. | ||
|
|
||
| ```@example snes_ex2 | ||
| nlfunc_dense = NonlinearFunction(form_residual!) | ||
| nlfunc_sparse = NonlinearFunction(form_residual!; sparsity = TracerSparsityDetector()) | ||
|
|
||
| nlprob_dense = NonlinearProblem(nlfunc_dense, u0) | ||
| nlprob_sparse = NonlinearProblem(nlfunc_sparse, u0) | ||
| ``` | ||
|
|
||
| Now we can solve the problem using `PETScSNES` or with one of the native `NonlinearSolve.jl` | ||
| solvers. | ||
|
|
||
| ```@example snes_ex2 | ||
| sol_dense_nr = solve(nlprob_dense, NewtonRaphson(); abstol = 1e-8) | ||
| sol_dense_snes = solve(nlprob_dense, PETScSNES(); abstol = 1e-8) | ||
| sol_dense_nr .- sol_dense_snes | ||
| ``` | ||
|
|
||
| ```@example snes_ex2 | ||
| sol_sparse_nr = solve(nlprob_sparse, NewtonRaphson(); abstol = 1e-8) | ||
| sol_sparse_snes = solve(nlprob_sparse, PETScSNES(); abstol = 1e-8) | ||
| sol_sparse_nr .- sol_sparse_snes | ||
| ``` | ||
|
|
||
| As expected the solutions are the same (upto floating point error). Now let's compare the | ||
| runtimes. | ||
|
|
||
| ## Runtimes | ||
|
|
||
| ### Dense Jacobian | ||
|
|
||
| ```@example snes_ex2 | ||
| @benchmark solve($(nlprob_dense), $(NewtonRaphson()); abstol = 1e-8) | ||
| ``` | ||
|
|
||
| ```@example snes_ex2 | ||
| @benchmark solve($(nlprob_dense), $(PETScSNES()); abstol = 1e-8) | ||
| ``` | ||
|
|
||
| ### Sparse Jacobian | ||
|
|
||
| ```@example snes_ex2 | ||
| @benchmark solve($(nlprob_sparse), $(NewtonRaphson()); abstol = 1e-8) | ||
| ``` | ||
|
|
||
| ```@example snes_ex2 | ||
| @benchmark solve($(nlprob_sparse), $(PETScSNES()); abstol = 1e-8) | ||
| ``` |
This file contains hidden or 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,120 @@ | ||
| module NonlinearSolvePETScExt | ||
|
|
||
| using FastClosures: @closure | ||
| using MPI: MPI | ||
| using NonlinearSolveBase: NonlinearSolveBase, get_tolerance | ||
| using NonlinearSolve: NonlinearSolve, PETScSNES | ||
| using PETSc: PETSc | ||
| using SciMLBase: SciMLBase, NonlinearProblem, ReturnCode | ||
| using SparseArrays: AbstractSparseMatrix | ||
|
|
||
| function SciMLBase.__solve( | ||
| prob::NonlinearProblem, alg::PETScSNES, args...; abstol = nothing, reltol = nothing, | ||
| maxiters = 1000, alias_u0::Bool = false, termination_condition = nothing, | ||
| show_trace::Val{ShT} = Val(false), kwargs...) where {ShT} | ||
| # XXX: https://petsc.org/release/manualpages/SNES/SNESSetConvergenceTest/ | ||
| termination_condition === nothing || | ||
| error("`PETScSNES` does not support termination conditions!") | ||
|
|
||
| _f!, u0, resid = NonlinearSolve.__construct_extension_f(prob; alias_u0) | ||
| T = eltype(prob.u0) | ||
| @assert T ∈ PETSc.scalar_types | ||
|
|
||
| if alg.petsclib === missing | ||
| petsclibidx = findfirst(PETSc.petsclibs) do petsclib | ||
| petsclib isa PETSc.PetscLibType{T} | ||
| end | ||
|
|
||
| if petsclibidx === nothing | ||
| error("No compatible PETSc library found for element type $(T). Pass in a \ | ||
| custom `petsclib` via `PETScSNES(; petsclib = <petsclib>, ....)`.") | ||
| end | ||
| petsclib = PETSc.petsclibs[petsclibidx] | ||
| else | ||
| petsclib = alg.petsclib | ||
| end | ||
| PETSc.initialized(petsclib) || PETSc.initialize(petsclib) | ||
|
|
||
| abstol = get_tolerance(abstol, T) | ||
| reltol = get_tolerance(reltol, T) | ||
|
|
||
| nf = Ref{Int}(0) | ||
|
|
||
| f! = @closure (cfx, cx, user_ctx) -> begin | ||
| nf[] += 1 | ||
| fx = cfx isa Ptr{Nothing} ? PETSc.unsafe_localarray(T, cfx; read = false) : cfx | ||
| x = cx isa Ptr{Nothing} ? PETSc.unsafe_localarray(T, cx; write = false) : cx | ||
| _f!(fx, x) | ||
| Base.finalize(fx) | ||
| Base.finalize(x) | ||
| return | ||
| end | ||
|
|
||
| snes = PETSc.SNES{T}(petsclib, | ||
| alg.mpi_comm === missing ? MPI.COMM_SELF : alg.mpi_comm; | ||
| alg.snes_options..., snes_monitor = ShT, snes_rtol = reltol, | ||
| snes_atol = abstol, snes_max_it = maxiters) | ||
|
|
||
| PETSc.setfunction!(snes, f!, PETSc.VecSeq(zero(u0))) | ||
|
|
||
| if alg.autodiff === missing && prob.f.jac === nothing | ||
| _jac! = nothing | ||
| njac = Ref{Int}(-1) | ||
| else | ||
| autodiff = alg.autodiff === missing ? nothing : alg.autodiff | ||
| if prob.u0 isa Number | ||
| _jac! = NonlinearSolve.__construct_extension_jac( | ||
| prob, alg, prob.u0, prob.u0; autodiff) | ||
| J_init = zeros(T, 1, 1) | ||
| else | ||
| _jac!, J_init = NonlinearSolve.__construct_extension_jac( | ||
| prob, alg, u0, resid; autodiff, initial_jacobian = Val(true)) | ||
| end | ||
|
|
||
| njac = Ref{Int}(0) | ||
|
|
||
| if J_init isa AbstractSparseMatrix | ||
| PJ = PETSc.MatSeqAIJ(J_init) | ||
| jac! = @closure (cx, J, _, user_ctx) -> begin | ||
| njac[] += 1 | ||
| x = cx isa Ptr{Nothing} ? PETSc.unsafe_localarray(T, cx; write = false) : cx | ||
| if J isa PETSc.AbstractMat | ||
| _jac!(user_ctx.jacobian, x) | ||
| copyto!(J, user_ctx.jacobian) | ||
| PETSc.assemble(J) | ||
| else | ||
| _jac!(J, x) | ||
| end | ||
| Base.finalize(x) | ||
| return | ||
| end | ||
| PETSc.setjacobian!(snes, jac!, PJ, PJ) | ||
| snes.user_ctx = (; jacobian = J_init) | ||
| else | ||
| PJ = PETSc.MatSeqDense(J_init) | ||
| jac! = @closure (cx, J, _, user_ctx) -> begin | ||
| njac[] += 1 | ||
| x = cx isa Ptr{Nothing} ? PETSc.unsafe_localarray(T, cx; write = false) : cx | ||
| _jac!(J, x) | ||
| Base.finalize(x) | ||
| J isa PETSc.AbstractMat && PETSc.assemble(J) | ||
| return | ||
| end | ||
| PETSc.setjacobian!(snes, jac!, PJ, PJ) | ||
| end | ||
| end | ||
|
|
||
| res = PETSc.solve!(u0, snes) | ||
|
|
||
| _f!(resid, res) | ||
| u_ = prob.u0 isa Number ? res[1] : res | ||
| resid_ = prob.u0 isa Number ? resid[1] : resid | ||
|
|
||
| objective = maximum(abs, resid) | ||
| # XXX: Return Code from PETSc | ||
| retcode = ifelse(objective ≤ abstol, ReturnCode.Success, ReturnCode.Failure) | ||
| return SciMLBase.build_solution(prob, alg, u_, resid_; retcode, original = snes, | ||
| stats = SciMLBase.NLStats(nf[], njac[], -1, -1, -1)) | ||
| end | ||
|
|
||
| end |
This file contains hidden or 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
Oops, something went wrong.
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.
Uh oh!
There was an error while loading. Please reload this page.