diff --git a/src/ksp.jl b/src/ksp.jl index 86f1cc29..aa272152 100644 --- a/src/ksp.jl +++ b/src/ksp.jl @@ -27,7 +27,7 @@ LinearAlgebra.adjoint(ksp) = LinearAlgebra.Adjoint(ksp) @for_libpetsc begin function KSP{$PetscScalar}(comm::MPI.Comm; kwargs...) - initialize($PetscScalar) + @assert initialized($PetscScalar) opts = Options{$PetscScalar}(kwargs...) ksp = KSP{$PetscScalar}(C_NULL, comm, nothing, nothing, opts) @chk ccall((:KSPCreate, $libpetsc), PetscErrorCode, (MPI.MPI_Comm, Ptr{CKSP}), comm, ksp) diff --git a/src/mat.jl b/src/mat.jl index ce6ad622..4afbc9fb 100644 --- a/src/mat.jl +++ b/src/mat.jl @@ -40,7 +40,7 @@ end @for_libpetsc begin function MatSeqAIJ{$PetscScalar}(m::Integer, n::Integer, nnz::Vector{$PetscInt}) - initialize($PetscScalar) + @assert initialized($PetscScalar) comm = MPI.COMM_SELF mat = MatSeqAIJ{$PetscScalar}(C_NULL, comm) @chk ccall((:MatCreateSeqAIJ, $libpetsc), PetscErrorCode, @@ -50,7 +50,7 @@ end return mat end function MatSeqDense(A::Matrix{$PetscScalar}) - initialize($PetscScalar) + @assert initialized($PetscScalar) comm = MPI.COMM_SELF mat = MatSeqDense(C_NULL, comm, A) @chk ccall((:MatCreateSeqDense, $libpetsc), PetscErrorCode, diff --git a/src/options.jl b/src/options.jl index fbf77f9f..5dcf9399 100644 --- a/src/options.jl +++ b/src/options.jl @@ -25,7 +25,7 @@ scalartype(::Options{T}) where {T} = T @for_libpetsc begin function Options{$PetscScalar}() - initialize($PetscScalar) + @assert initialized($PetscScalar) opts = Options{$PetscScalar}(C_NULL) @chk ccall((:PetscOptionsCreate, $libpetsc), PetscErrorCode, (Ptr{CPetscOptions},), opts) finalizer(destroy, opts) diff --git a/src/snes.jl b/src/snes.jl index fde7d83a..014d1182 100644 --- a/src/snes.jl +++ b/src/snes.jl @@ -49,7 +49,7 @@ end @for_libpetsc begin function SNES{$PetscScalar}(comm::MPI.Comm; kwargs...) - initialize($PetscScalar) + @assert initialized($PetscScalar) opts = Options{$PetscScalar}(kwargs...) snes = SNES{$PetscScalar}(C_NULL, comm, opts, nothing, nothing, nothing, nothing, nothing) @chk ccall((:SNESCreate, $libpetsc), PetscErrorCode, (MPI.MPI_Comm, Ptr{CSNES}), comm, snes) @@ -157,4 +157,4 @@ end end -solve!(x::AbstractVector{T}, snes::SNES{T}) where {T} = parent(solve!(AbstractVec(x), snes)) \ No newline at end of file +solve!(x::AbstractVector{T}, snes::SNES{T}) where {T} = parent(solve!(AbstractVec(x), snes)) diff --git a/src/vec.jl b/src/vec.jl index 7d640877..257fbbeb 100644 --- a/src/vec.jl +++ b/src/vec.jl @@ -39,7 +39,7 @@ Base.parent(v::AbstractVec) = v.array @for_libpetsc begin function VecSeq(comm::MPI.Comm, X::Vector{$PetscScalar}; blocksize=1) - initialize($PetscScalar) + @assert initialized($PetscScalar) v = VecSeq(C_NULL, comm, X) @chk ccall((:VecCreateSeqWithArray, $libpetsc), PetscErrorCode, (MPI.MPI_Comm, $PetscInt, $PetscInt, Ptr{$PetscScalar}, Ptr{CVec}), diff --git a/test/runtests.jl b/test/runtests.jl index a22fd5f7..8fef3171 100644 --- a/test/runtests.jl +++ b/test/runtests.jl @@ -1,5 +1,6 @@ using Test using PETSc, MPI, LinearAlgebra, SparseArrays +PETSc.initialize() @testset "Tests" begin m,n = 20,20