diff --git a/Project.toml b/Project.toml index 7af4146..54c3a99 100644 --- a/Project.toml +++ b/Project.toml @@ -1,8 +1,13 @@ name = "ModelOrderReduction" uuid = "207801d6-6cee-43a9-ad0c-f0c64933efa0" -authors = ["Bowen S. Zhu"] +authors = ["Bowen S. Zhu", "Rahul Manavalan"] version = "0.1.0" +[deps] +LinearAlgebra = "37e2e46d-f89d-539d-b4ee-838fcccc9c8e" +RandomizedLinAlg = "0448d7d9-159c-5637-8537-fd72090fea46" +TSVD = "9449cd9e-2762-5aa3-a617-5413e99d722e" + [compat] julia = "1.6" @@ -10,4 +15,4 @@ julia = "1.6" Test = "8dfed614-e22c-5e08-85e1-65c5234f0b40" [targets] -test = ["Test"] +test = ["Test"] \ No newline at end of file diff --git a/examples/POD/lorenz.jl b/examples/POD/lorenz.jl new file mode 100644 index 0000000..e32cec8 --- /dev/null +++ b/examples/POD/lorenz.jl @@ -0,0 +1,39 @@ +using ModelOrderReduction +using Plots +using OrdinaryDiffEq + +function lorenz_prob() + function lorenz!(du, u, p, t) + du[1] = p[1] * (u[2] - u[1]) + du[2] = u[1] * (p[2] - u[3]) - u[2] + du[3] = u[1] * u[2] - p[3] * u[3] + end + + u0 = [1, 0, 0] + p = [10, 28, 8 / 3] + tspan = (0, 100) + prob = ODEProblem(lorenz!, u0, tspan, p) + sol = solve(prob, Tsit5()) + sol +end + +sol = lorenz_prob() +solution = Array(sol) +plot(solution[1, :], solution[2, :], solution[3, :]) +# savefig("lorenz_attractor.png") + +## Two way POD +reducer = POD(solution, 2) +reduce!(reducer, SVD()) +bases = reducer.rbasis +plot(bases[:, 1], bases[:, 2], label = "POD2") +# savefig("pod2.png") + +## One way POD +reducer = POD(solution, 1) +reduce!(reducer, SVD()) +bases = reducer.rbasis +plot(bases[:, 1], label = "POD1") +# savefig("pod1.png") +plot(solution[:, 3], label = 'z') +# savefig("z.png") diff --git a/src/DataReduction/POD.jl b/src/DataReduction/POD.jl new file mode 100644 index 0000000..7c7dfce --- /dev/null +++ b/src/DataReduction/POD.jl @@ -0,0 +1,68 @@ +function matricize(VoV::Vector{Vector{FT}}) where {FT} + Matrix(reduce(hcat, VoV)) +end + +mutable struct POD{FT} <: AbstractDRProblem + snapshots::Union{Vector{Vector{FT}}, Matrix{FT}} + nmodes::Int + rbasis::Matrix{FT} + renergy::FT + + function POD(snaps::Vector{Vector{FT}}, nmodes::Int) where {FT} + errorhandle(matricize(snaps), nmodes) + new{eltype(snaps[1])}(snaps, nmodes, Array{FT, 2}(undef, size(snaps, 1), nmodes), + FT(0)) + end + + function POD(snaps::Matrix{FT}, nmodes::Int) where {FT} + errorhandle(snaps, nmodes) + new{eltype(snaps)}(snaps, nmodes, Array{FT, 2}(undef, size(snaps, 1), nmodes), + FT(0)) + end + + function POD(snaps::Adjoint{FT, Matrix{FT}}, nmodes::Int) where {FT} + POD(Matrix(snaps), nmodes, Array{FT, 2}(undef, size(snaps, 1), nmodes), FT(0)) + end +end + +function reduce!(pod::POD{FT}, ::SVD) where {FT} + op_matrix = pod.snapshots + if typeof(pod.snapshots) == Vector{Vector{FT}} + op_matrix = matricize(pod.snapshots) + end + u, s, v = svd(op_matrix) + pod.rbasis .= u[:, 1:(pod.nmodes)] + sr = s[1:(pod.nmodes)] + pod.renergy = sum(s) / (sum(s) + (size(op_matrix, 1) - pod.nmodes) * s[end]) + nothing +end + +function reduce!(pod::POD{FT}, ::TSVD) where {FT} + op_matrix = pod.snapshots + if typeof(pod.snapshots) == Vector{Vector{FT}} + op_matrix = matricize(pod.snapshots) + end + u, s, v = tsvd(op_matrix, pod.nmodes) + pod.renergy = sum(s) / (sum(s) + (size(op_matrix, 1) - pod.nmodes) * s[end]) + pod.rbasis .= u + nothing +end + +function reduce!(pod::POD{FT}, ::RSVD) where {FT} + op_matrix = pod.snapshots + if typeof(pod.snapshots) == Vector{Vector{FT}} + op_matrix = matricize(pod.snapshots) + end + u, s, v = rsvd(op_matrix, pod.nmodes) + pod.renergy = sum(s) / (sum(s) + (size(op_matrix, 1) - pod.nmodes) * s[end]) + pod.rbasis .= u + nothing +end + +function Base.show(io::IO, pod::POD) + print(io, "POD \n") + print(io, "Reduction Order = ", pod.nmodes, "\n") + print(io, "Snapshot size = (", size(pod.snapshots, 1), ",", size(pod.snapshots[1], 2), + ")\n") + print(io, "Relative Energy = ", pod.renergy, "\n") +end diff --git a/src/ErrorHandle.jl b/src/ErrorHandle.jl new file mode 100644 index 0000000..acba209 --- /dev/null +++ b/src/ErrorHandle.jl @@ -0,0 +1,5 @@ +function errorhandle(data::Matrix{FT}, modes::IT) where {FT, IT} + @assert size(data, 1)>1 "State vector is expected to be vector valued." + s = size(data, 2) + @assert (modes > 0)&(modes < s) "Number of modes should be in {1,2,...,$(s-1)}." +end diff --git a/src/ModelOrderReduction.jl b/src/ModelOrderReduction.jl index 5ae21ae..c94ac26 100644 --- a/src/ModelOrderReduction.jl +++ b/src/ModelOrderReduction.jl @@ -1,5 +1,17 @@ module ModelOrderReduction +#========================Data Reduction=========================================# +include("Types.jl") +include("ErrorHandle.jl") -# Write your package code here. +using LinearAlgebra +using TSVD +using RandomizedLinAlg +include("DataReduction/POD.jl") + +export SVD, TSVD, RSVD +export POD, reduce!, matricize +#========================Model Reduction========================================# + +#===============================================================================# end diff --git a/src/Types.jl b/src/Types.jl new file mode 100644 index 0000000..d877a07 --- /dev/null +++ b/src/Types.jl @@ -0,0 +1,8 @@ +abstract type AbstractReductionProblem end +abstract type AbstractMORProblem <: AbstractReductionProblem end +abstract type AbstractDRProblem <: AbstractReductionProblem end + +abstract type AbstractSVD end +struct SVD <: AbstractSVD end +struct TSVD <: AbstractSVD end +struct RSVD <: AbstractSVD end diff --git a/test/DataReduction.jl b/test/DataReduction.jl new file mode 100644 index 0000000..e177d89 --- /dev/null +++ b/test/DataReduction.jl @@ -0,0 +1,59 @@ +#--------- Data Reduction -----------------# + +function lorenz_prob() + function lorenz!(du, u, p, t) + du[1] = p[1] * (u[2] - u[1]) + du[2] = u[1] * (p[2] - u[3]) - u[2] + du[3] = u[1] * u[2] - p[3] * u[3] + end + + u0 = [1, 0, 0] + p = [10, 28, 8 / 3] + tspan = (0, 100) + prob = ODEProblem(lorenz!, u0, tspan, p) + sol = solve(prob, Tsit5()) + sol +end + +@testset "POD-Utils" begin + solution = lorenz_prob() + VoV = solution.u + M = matricize(VoV) + @test size(M, 1) == length(VoV[1]) # Parameters + @test size(M, 2) == length(VoV) # Time +end + +@testset "POD - Attractor Test" begin + sol = lorenz_prob() + solution = Array(sol) + + order = 2 + solver = SVD() + reducer = POD(solution, order) + reduce!(reducer, solver) + + reducer.renergy + + # Ad-hoc tests. To be checked with Chris. + @test size(reducer.rbasis, 2) == reducer.nmodes + @test size(reducer.rbasis, 1) == size(solution, 1) + @test reducer.renergy > 0.9 + + order = 2 + solver = TSVD() + reducer = POD(solution, order) + reduce!(reducer, solver) + + # Ad-hoc tests. To be checked with Chris. + @test size(reducer.rbasis, 2) == reducer.nmodes + @test size(reducer.rbasis, 1) == size(solution, 1) + + order = 2 + solver = RSVD() + reducer = POD(solution, order) + reduce!(reducer, solver) + + # Ad-hoc tests. To be checked with Chris. + @test size(reducer.rbasis, 2) == reducer.nmodes + @test size(reducer.rbasis, 1) == size(solution, 1) +end diff --git a/test/Project.toml b/test/Project.toml new file mode 100644 index 0000000..9aaa105 --- /dev/null +++ b/test/Project.toml @@ -0,0 +1,3 @@ +[deps] +OrdinaryDiffEq = "1dea7af3-3e70-54e6-95c3-0bf5283fa5ed" +Test = "8dfed614-e22c-5e08-85e1-65c5234f0b40" diff --git a/test/runtests.jl b/test/runtests.jl index 1f356e7..434234f 100644 --- a/test/runtests.jl +++ b/test/runtests.jl @@ -1,6 +1,6 @@ using ModelOrderReduction using Test +using OrdinaryDiffEq -@testset "ModelOrderReduction.jl" begin - # Write your tests here. -end +include("DataReduction.jl") +#---------- Model Reduction ----------------#