RSCG.jl

The Reduced-Shifted Conjugate-Gradient Method method to calculate the element of the Green's function matrix. See, Y. Nagai, Y. Shinohara, Y. Futamura, and T. Sakurai,[arXiv:1607.03992v2 or DOI:10.7566/JPSJ.86.014708]. http://dx.doi.org/10.7566/JPSJ.86.014708
Author cometscome
Popularity
0 Stars
Updated Last
2 Years Ago
Started In
October 2018

Build Status

Coverage Status

RSCG.jl

This package can calculate the elements of the Green's function:

$$G_ij(σk) = ([σj I - A]^-1)_{ij},$$

with the use of the reduced-shifted conjugate gradient method (See, Y. Nagai, Y. Shinohara, Y. Futamura, and T. Sakurai,[arXiv:1607.03992v2 or DOI:10.7566/JPSJ.86.014708]). One can obtain G_{ij}(\sigma_k) with different frequencies \sigma_k, simultaneously.

The matrix should be symmetric or hermitian.

We can use Arrays, LinearMaps, and SparseArrays.

This software is written in Julia 1.0.

This software is released under the MIT License, see LICENSE.

Install

add RSCG

Example

Let us obtain the Green' functions G(z) on the complex plane.

M = 100
σ = zeros(ComplexF64,M)
σmin = -4.0*im-1.2
σmax = 4.0*im +4.3
for i=1:M
    σ[i] = (i-1)*(σmax-σmin)/(M-1) + σmin
end

We define the matrix:

using SparseArrays


function make_mat(n)
    A = spzeros(Float64,n,n)
    t = -1.0
    μ = -1.5
    for i=1:n
        dx = 1
        jp = i+dx
        jp += ifelse(jp > n,-n,0) #+1
        dx = -1
        jm = i+dx
        jm += ifelse(jm < 1,n,0) #-1
        A[i,jp] = t
        A[i,i] = -μ
        A[i,jm] = t
    end
    return A
end

n=1000
A1 = make_mat(n)

Or, we can also use LinearMaps.jl to define the matrix:

using LinearMaps

function set_diff(v)
    function calc_diff!(y::AbstractVector, x::AbstractVector)
        n = length(x)
        length(y) == n || throw(DimensionMismatch())
        μ = -1.5
        for i=1:n
            dx = 1
            jp = i+dx
            jp += ifelse(jp > n,-n,0) #+1
            dx = -1
            jm = i+dx
            jm += ifelse(jm < 1,n,0) #-1
            y[i] = v*(x[jp]+x[jm])-μ*x[i]
        end

        return y
    end
    (y,x) -> calc_diff!(y,x)
end

n=1000
Al = set_diff(-1.0)
A2 = LinearMap(Al, n; ismutating=true,issymmetric=true)

an element

If we want to obtain the element G_{ij}(σ_k),

i = 1
j = 1
Gij1 = greensfunctions(i,j,σ,A1) #SparseArrays
Gij2 = greensfunctions(i,j,σ,A2) #LinearMaps

elements

If we want to obtain the elements G_{ij}(σ_k) with different i,

vec_i = [1,4,8,43,98]
j = 1
vec_Gij1 = greensfunctions(vec_i,j,σ,A1) #SparseArrays
vec_Gij2 = greensfunctions(vec_i,j,σ,A2) #LinearMaps

Functions

greensfunctions(i::Integer,j::Integer,σ::Array{ComplexF64,1},A)

Inputs:

  • i :index of the Green's function

  • j :index of the Green's function

  • σ :frequencies

  • A :hermitian matrix. We can use Arrays,LinearMaps, SparseArrays

  • eps :residual (optional) Default:1e-12

  • maximumsteps : maximum number of steps (optional) Default:20000

Output:

  • Gij[1:M]: the matrix element Green's functions at M frequencies defined by \sigma_k.

greensfunctions(vec_left::Array{<:Integer,1},j::Integer,σ::Array{ComplexF64,1},A)

Inputs:

  • vec_left :i indices of the Green's function

  • j :index of the Green's function

  • σ :frequencies

  • A :hermitian matrix. We can use Arrays,LinearMaps, SparseArrays

  • eps :residual (optional) Default:1e-12

  • maximumsteps : maximum number of steps (optional) Default:20000

Output:

  • Gij[1:M,1:length(vec_left)]: the matrix element Green's functions at M frequencies defined by \sigma_k.

greensfunctions_col(j::Integer,σ::Array{ComplexF64,1},A)

Inputs:

  • σ :frequencies

  • A :hermitian matrix. We can use Arrays,LinearMaps, SparseArrays

  • b :input vector

  • eps :residual (optional) Default:1e-12

  • maximumsteps : maximum number of steps (optional) Default:20000

Output:

  • Gij[1:M,1:size(A)[1]]: the matrix elements Green's functions of j-th col at M frequencies defined by \sigma_k.

Used By Packages