Core

Provide basic user's interfaces for the ACFlow toolkit.

Contents

Index

Solvers

ACFlow.solveMethod
solve(grid::Vector{F64}, Gval::Vector{T}, Gerr::Vector{T})

Solve the analytic continuation problem. The arguments grid, Gval, and Gerr are the grid, value, and error bar, respectively.

Arguments

  • grid -> Imaginary axis grid for correlators, τ or ωₙ.
  • Gval -> Function values for correlators, G(τ) or G(iωₙ).
  • Gerr -> Standard deviations for correlators.

Returns

  • mesh -> Real frequency mesh, ω, Vector{F64}.
  • Aout -> Spectral function, A(ω), Vector{F64}.
  • Gout -> Retarded Green's function, G(ω), Vector{C64}.
ACFlow.solveMethod
solve(grid::Vector{F64}, Gval::Vector{T}, err::T)

Solve the analytic continuation problem. The arguments grid, Gval, and err are the grid, value, and error bar, respectively.

Here, we just assume that the standard deviations for correlators are fixed to a constant value, err.

Arguments

  • grid -> Imaginary axis grid for correlators, τ or ωₙ.
  • Gval -> Function values for correlators, G(τ) or G(iωₙ).
  • err -> Standard deviations for correlators.

Returns

  • mesh -> Real frequency mesh, ω, Vector{F64}.
  • Aout -> Spectral function, A(ω), Vector{F64}.
  • Gout -> Retarded Green's function, G(ω), Vector{C64}.
ACFlow.solveMethod
solve(grid::Vector{F64}, Gval::Vector{T})

Solve the analytic continuation problem. The arguments grid and Gval are the grid and value, respectively. Furthermore, the error bar is set to a fixed value 1.0e-4.

Arguments

  • grid -> Imaginary axis grid for correlators, τ or ωₙ.
  • Gval -> Function values for correlators, G(τ) or G(iωₙ).

Returns

  • mesh -> Real frequency mesh, ω, Vector{F64}.
  • Aout -> Spectral function, A(ω), Vector{F64}.
  • Gout -> Retarded Green's function, G(ω), Vector{C64}.
ACFlow.solveMethod
solve(rd::RawData)

Solve the analytic continuation problem. The input data are encapsulated in a RawData struct. This function call is the actual interface to the desired analytic continuation solvers.

Arguments

  • rd -> A RawData struct that contains the grid, correator, and error bar.

Returns

  • mesh -> Real frequency mesh, ω, Vector{F64}.
  • Aout -> Spectral function, A(ω) or A(ω) / ω, Vector{F64}.
  • Gout -> Retarded Green's function, G(ω), Vector{C64}.

Examples

# Setup the configuration file
setup_args("ac.toml")

# Read the parameters
read_param()

# Call the solver.
mesh, Aout, Gout = solve(read_data())

See also: RawData.

Parameters

ACFlow.setup_paramFunction
setup_param(
    C::Dict{String,Any},
    S::Dict{String,Any},
    reset::Bool = true
)

Setup the configuration dictionaries via function call. Here C contains parameters for general setup, while S contains parameters for selected analytic continuation solver. If reset is true, then the configuration dictionaries will be reset to their default values at first. Later, C and S will be used to customized the dictionaries further.

Arguments

See above explanations.

Returns

N/A

See also: read_param.

ACFlow.read_paramFunction
read_param()

Setup the configuration dictionaries via an external file. The valid format of a configuration file is toml.

Arguments

N/A

Returns

N/A

See also: setup_param.

Data

ACFlow.read_dataFunction
read_data(only_real_part::Bool = true)

Read data in imaginary axis and return a RawData struct. The argument only_real_part is only useful for bosonic cases. If the kernel is bsymm (it means a symmetric bosonic kernel) and the grid is bfreq or bfrag (it means Matsubara frequency grid), the function values for input correators should be real in principle. In these cases, we should set only_real_part = true.

Arguments

See above explanations.

Returns

  • rd -> Raw input data that is encapsulated in a RawData struct.

See also: RawData.

ACFlow.make_dataFunction
make_data(rd::RawData; T::DataType = F64)

Convert RawData struct to GreenData struct. Note that RawData is provided by the users directly, while GreenData is more suitable for various analytic continuation solvers and algorithms. Note that the GreenData struct is accessed and manipulated by this code internally, while the RawData struct is exposed to the users.

Arguments

See above explanations.

Returns

  • gd -> A GreenData struct.

See also: RawData, GreenData.

Grids

ACFlow.make_gridFunction
make_grid(rd::RawData; T::DataType = F64)

Extract grid for input data from a RawData struct (RD). It will return a sub-type of the AbstractGrid struct.

Arguments

See above explanations.

Returns

  • grid -> Imaginary time or imaginary frequency grid.

See also: RawData, AbstractGrid.

Meshes

ACFlow.make_meshFunction
make_mesh(; T::DataType = F64)

Try to generate an uniform (linear) or non-uniform (non-linear) mesh for the spectral function in real axis. Notice that it supports arbitrary precision mesh. By default, the precision is F64. One can specify the precision by the argument T.

Arguments

See above explanations.

Returns

  • mesh -> Real frequency mesh. It should be a subtype of AbstractMesh.

See also: LinearMesh, TangentMesh, LorentzMesh.

Models

ACFlow.make_modelFunction
make_model(am::AbstractMesh)

Try to generate a default model function at given mesh am through various schemes.

Arguments

  • am -> Real frequency mesh.

Returns

  • model -> A normalized model function.

See also: AbstractMesh.

Kernels

ACFlow.make_kernelFunction
make_kernel(am::AbstractMesh, ag::AbstractGrid)

Try to generate various kernel functions.

Arguments

  • am -> Real frequency mesh.
  • ag -> Imaginary axis grid.

Returns

  • kernel -> Kernel function, a 2D array, (ntime,nmesh) or (nfreq,nmesh).

See also: AbstractMesh, AbstractGrid.

Postprocessing

ACFlow.reprodFunction
reprod(am::AbstractMesh, kernel::Matrix{F64}, A::Vector{F64})

Try to reproduce the input data, which can be compared with the raw data to see whether the analytic continuation is reasonable.

Arguments

  • am -> Real frequency mesh.
  • kernel -> The kernel function.
  • A -> The calculated spectral function, A(ω) or A(ω) / ω.

Returns

  • G -> Reconstructed correlators, G(τ) or G(iωₙ), Vector{F64}.

See also: AbstractMesh.

ACFlow.kramersFunction
kramers(am::AbstractMesh, A::Vector{F64})

Try to calculate the real part of the Green's function from its imaginary part via the Kramers-Kronig relations. The objective of this function is to get the full retarded Green's function.

Arguments

  • am -> Real frequency mesh.
  • A -> Spectral function at real frequency mesh, A(ω).

Returns

  • G -> Retarded Green's function, G(ω), Vector{C64}.