src/PET/AlgorithmProximal.jl

Fri, 03 May 2024 18:03:06 +0300

author
Neil Dizon <neil.dizon@helsinki.fi>
date
Fri, 03 May 2024 18:03:06 +0300
changeset 52
cb029cdb141a
parent 36
e4a8f662a1ac
permissions
-rw-r--r--

activation function for dual scscaling

####################################################################
# Predictive online PDPS for optical flow with known velocity field
####################################################################

__precompile__()

module AlgorithmProximal

identifier = "pdps_known_proximal"

using Printf

using AlgTools.Util
import AlgTools.Iterate
using ImageTools.Gradient
using ImageTools.Translate
  
using ImageTransformations
using Images, CoordinateTransformations, Rotations, OffsetArrays
using ImageCore, Interpolations

using ...Radon
using ...OpticalFlow: ImageSize,
                      Image,
                      petpdflow!

#########################
# Iterate initialisation
#########################

function init_rest(x::Image)
    imdim=size(x)
    y = zeros(2, imdim...)
    Δx = copy(x)
    Δy = copy(y)
    x̄ = copy(x)
    radonx = copy(x)
    return x, y, Δx, Δy, x̄, radonx
end

function init_iterates(xinit::Image)
    return init_rest(copy(xinit))
end

function init_iterates(dim::ImageSize)
    return init_rest(zeros(dim...))
end

#########################
# PETscan related
#########################
function petvalue(x, b, c)
    tmp = similar(b)
    radon!(tmp, x)
    return sum(@. tmp - b*log(tmp+c))
end

function petgrad!(res, x, b, c, S)
    tmp = similar(b)
    radon!(tmp, x)
    @. tmp = S .- b/(tmp+c)
    backproject!(res, S.*tmp)
end

function proj_nonneg!(y)
    @inbounds @simd for i=1:length(y)
        if y[i] < 0
            y[i] = 0
        end
    end
    return y
end

############
# Algorithm
############

function step_lengths(params, γ, R_K², L)
    ρ̃₀, τ₀, σ₀, σ̃₀ =  params.ρ̃₀, params.τ₀, params.σ₀, params.σ̃₀
    δ = params.δ
    ρ = isdefined(params, :phantom_ρ) ? params.phantom_ρ : params.ρ
    Λ = params.Λ
    Θ = params.dual_flow ? Λ : 1

    τ = τ₀/L
    @assert(1+γ*τ ≥ Λ)
    σ = σ₀*min((1-τ₀)/(τ*R_K²), 1/max(0, τ*R_K²/((1+γ*τ-Λ)*(1-δ))-ρ))
    q = δ*(1+σ*ρ)/Θ
    if 1 ≥ q
        σ̃ = σ̃₀*σ/q
        #ρ̃ = ρ̃₀*max(0, ((Θ*σ)/(2*δ*σ̃^2*(1+σ*ρ))+1/(2σ)-1/σ̃))
        ρ̃ = max(0, (1-q)/(2*σ))
    else
        σ̃ = σ̃₀*σ/(q*(1-√(1-1/q)))
        ρ̃ = 0
    end

    return τ, σ, σ̃, ρ̃
end

function solve( :: Type{DisplacementT};
               dim :: ImageSize,
               iterate = AlgTools.simple_iterate,
               params::NamedTuple) where DisplacementT

    ################################
    # Extract and set up parameters
    ################################
    α, ρ = params.α, params.ρ
    R_K² = ∇₂_norm₂₂_est²
    γ = 1
    L = params.L
    τ, σ, σ̃, ρ̃ = step_lengths(params, γ, R_K², L)
    println("Step length parameters: τ=$(τ), σ=$(σ), σ̃=$(σ̃), ρ̃=$(ρ̃)")

    λ = params.λ
    c = params.c*ones(params.radondims...)

    

    ######################
    # Initialise iterates
    ######################

    x, y, Δx, Δy, x̄, r∇ = init_iterates(dim)
    
    if params.L_experiment
        oldpetgradx = zeros(size(x)...)
        petgradx = zeros(size(x))
        oldx = ones(size(x))
    end

    ####################
    # Run the algorithm
    ####################
                        
    v = iterate(params) do verbose :: Function,
                           b :: Image,                   # noisy_sinogram
                           v_known :: DisplacementT,
                           theta_known :: DisplacementT,
                           b_true :: Image,
                           S :: Image

        ###################
        # Prediction steps
        ###################

        petpdflow!(x, Δx, y, Δy, v_known, theta_known, params.dual_flow)        # Usual flow
        
        if params.L_experiment
            @. oldx = x
        end

        ##############################
        # Proximal step of prediction
        ##############################
        ∇₂!(Δy, x)
        @. y = (y + σ̃*Δy)/(1 + σ̃*(ρ̃+ρ/α))
        #@. cc = y + 1000000*σ̃*Δy
        #@. y = (y + σ̃*Δy)/(1 + σ̃*(ρ̃+ρ/α)) + (1 - 1/(1 + ρ̃*σ̃))*cc
        proj_norm₂₁ball!(y, α)

        ############
        # PDPS step
        ############

        ∇₂ᵀ!(Δx, y)                    # primal step:
        @. x̄ = x                       # | save old x for over-relax
        petgrad!(r∇, x, b, c, S)       # | Calculate gradient of fidelity term

        @. x = x-(τ*λ)*r∇-τ*Δx         # |
        proj_nonneg!(x)                # | non-negativity constaint prox
        @. x̄ = 2*x - x̄                 # over-relax: x̄ = 2x-x_old
        ∇₂!(Δy, x̄)                     # dual step:
        @. y = y + σ*Δy                # |
        proj_norm₂₁ball!(y, α)         # |  prox

        #####################
        # L update if needed
        #####################
        if params.L_experiment
            petgrad!(petgradx, x, b, c, S)
            petgrad!(oldpetgradx, oldx, b, c, S)
            if norm₂(x-oldx)>1e-12
                L = max(0.9*norm₂(petgradx - oldpetgradx)/norm₂(x-oldx),L)
                println("Step length parameters: L=$(L)")
                τ = τ₀/L
                σ = σ₀*(1-τ₀)/(R_K²*τ)
            end
        end
       
        ################################
        # Give function value if needed
        ################################
        
        v = verbose() do
            ∇₂!(Δy, x)
            value = λ*petvalue(x, b, c) + params.α*norm₂₁(Δy)
            value, x, [NaN, NaN], nothing
        end
        
        v
    end

    return x, y, v
end

end # Module

mercurial