Input and output
Define input and output functions for the ACFlow toolkit.
Contents
Index
ACFlow.read_cmplx_data
ACFlow.read_real_data
ACFlow.write_backward
ACFlow.write_barycentric
ACFlow.write_complete
ACFlow.write_goodness
ACFlow.write_hamiltonian
ACFlow.write_misfit
ACFlow.write_model
ACFlow.write_passed
ACFlow.write_pole
ACFlow.write_probability
ACFlow.write_prony
ACFlow.write_spectrum
ACFlow.write_statistics
Read Data
ACFlow.read_real_data
— Functionread_real_data(finput::AbstractString, ngrid::I64)
Read input data. This function is used for imaginary time data. The input file should contain three columns. The first column is the imaginary time grid, the second column is the value, the third column is the standard deviation σ. Here, ngrid
specifies the number of grid points.
Arguments
- finput -> Filename for the input data.
- ngrid -> Number of grid points.
Returns
- rd -> A RawData struct.
See also: read_cmplx_data
.
ACFlow.read_cmplx_data
— Functionread_cmplx_data(finput::AbstractString, ngrid::I64)
Read input data. This function is used for Matsubara frequency data. The input should contain four columns or five columns. The first column is the Matsubara freqency grid, the second and third columns are the values (real part and imaginary part), the four and fifth columns are the standard deviations σ for the real and imaginary parts, respectively. If there are only four columns, it means that the real and imaginary parts share the same standard deviations.
Arguments
- finput -> Filename for the input data.
- ngrid -> Number of grid points.
Returns
- rd -> A RawData struct.
See also: read_real_data
.
read_cmplx_data(
finput::AbstractString,
ngrid::I64,
only_real_part::Bool
)
Read input data. This function is used for Matsubara frequency data. The input file only contains three columns. The first column is the Matsubara frequency grid, the second column is the real part or imaginary part of the data (which is specified by the argument only_real_part
), and the third column is the standard deviation σ. This function is for bosonic correlation function.
Arguments
- finput -> Filename for the input data.
- ngrid -> Number of grid points.
- onlyrealpart -> See above explanations.
Returns
- rd -> A RawData struct.
See also: read_real_data
.
Write Data
ACFlow.write_spectrum
— Functionwrite_spectrum(am::AbstractMesh, Aout::Vector{F64})
Write spectrum A(ω) or A(ω) / ω to Aout.data
. The grid is defined in am
, and the spectral data are contained in Aout
.
Note that for the MaxEnt, StochAC, StochSK, and StochOM solvers, Aout
is actually A(ω) / ω, instead of A(ω). However, for the BarRat, NevanAC, and StochPX solvers, Aout
is just A(ω).
Arguments
- am -> Real frequency mesh.
- Aout -> Spectral function. See above explanations.
Returns
N/A
write_spectrum(am::AbstractMesh, αₗ::Vector{F64}, Aout::Array{F64,2})
Write α-resolved spectrum A(ω) to Aout.data.alpha
. The grid is defined in am
, the α-resolved spectrum is contained in Aout
, αₗ
is the list for the α parameters. This function is called by the StochAC
solver.
Note that Aout
is actually A(ω) / ω, instead of A(ω).
Arguments
- am -> Real frequency mesh.
- αₗ -> List for α parameters.
- Aout -> α-dependent spectral function.
Returns
N/A
ACFlow.write_backward
— Functionwrite_backward(ag::AbstractGrid, G::Vector{F64})
We can use the calculated spectrum in real axis to reproduce the input data in imaginary axis. This function will write the reproduced data to repr.data
, which can be compared with the original data. Here, G
is the reproduced data.
Arguments
- ag -> Grid for input data, τ or iωₙ.
- G -> Reconstructed Green's function, G(τ) or G(iωₙ).
Returns
N/A
See also: reprod
.
ACFlow.write_complete
— Functionwrite_complete(am::AbstractMesh, G::Vector{C64})
Write the full data at real axis to Gout.data
. am
denotes the real axis, G
is the calculated Green's function data. Note that its real part is obtained via the so-called Kramers-Kronig transformation.
Arguments
- am -> Real frequency mesh, ω.
- G -> Retarded Green's function, G(ω).
Returns
N/A
See also: kramers
.
ACFlow.write_misfit
— Functionwrite_misfit(α_vec::Vector{F64}, χ²_vec::Vector{F64})
Write log10(α)-log10(χ²)
data to chi2.data
, which could be used to judge whether the obtained optimal α parameter is reasonable. It is used by the MaxEnt
solver only.
Arguments
- α_vec -> List for α parameters.
- χ²_vec -> α-dependent goodness-of-fit functional.
Returns
N/A
See also: write_goodness
.
ACFlow.write_goodness
— Functionwrite_goodness(Θ_vec::Vector{F64}, χ²_vec::Vector{F64})
Write log10(Θ)-log10(χ²)
data to goodness.data
, which could be used to judge whether the obtained optimal Θ parameter is reasonable. This function is only useful for the StochSK
solver.
Arguments
- Θ_vec -> List for Θ parameters.
- χ²_vec -> Θ-dependent goodness-of-fit functional.
Returns
N/A
See also: write_misfit
.
ACFlow.write_model
— Functionwrite_model(am::AbstractMesh, D::Vector{F64})
Write the default model function to model.data
. This function is usually for the MaxEnt
solver.
Arguments
- am -> Real frequency mesh, ω.
- D -> Default model, m(ω).
Returns
N/A
ACFlow.write_prony
— Functionwrite_prony(𝑁ₚ::I64, Γₚ::Vector{C64}, Ωₚ::Vector{C64})
Write Prony approximation to the input correlator. This information can be used to reconstruct or interpolate the correlator. This function is only useful for the BarRat
solver.
Arguments
- 𝑁ₚ -> Number of nodes for Prony approximation.
- Γₚ -> Nodes for Prony approximation, $γ_i$.
- Ωₚ -> Weights for Prony approximation, $w_i$.
Returns
N/A
write_prony(ag::AbstractGrid, G::Vector{C64})
Write Matsubara Green's function approximated by the Prony approximation. This function is only useful for the BarRat
solver.
Arguments
- ag -> Grid for input data, iωₙ.
- G -> Approximated Green's function, G(iωₙ).
Returns
N/A
ACFlow.write_barycentric
— Functionwrite_barycentric(
nodes::Vector{C64},
values::Vector{C64},
weights::Vector{C64}
)
Write barycentric rational function approximation to the input correlator. This information can be used to reconstruct or interpolate the correlator. This function is only useful for the BarRat
solver.
Arguments
- nodes -> Nodes of the rational function, $z_i$.
- values -> Values of the rational function, $r(z_i)$.
- weights -> Weights of the rational function, $w_i$.
Returns
N/A
ACFlow.write_hamiltonian
— Functionwrite_hamiltonian(α_vec::Vector{F64}, Uα::Vector{F64})
Write α-U(α)
data to hamil.data
, which could be used to judge whether the obtained optimal α parameter is reasonable. This function is only useful for the StochAC
solver.
Arguments
- α_vec -> List for α parameters.
- Uα -> α-dependent Hamiltonian.
Returns
N/A
ACFlow.write_passed
— Functionwrite_passed(passed::Vector{I64}, med::F64, αgood::F64)
Write indices of selected solutions which should be used to calculate the averaged spectrum. Here, passed
means the indices, med
is the median value of χ², and αgood
is the factor that is used to filter the solutions. This function is only useful for the StochOM
and the StochPX
solvers.
Arguments
- passed -> Indices for selected solutions.
- med -> Median value of χ².
- αgood -> Predefined parameter used to filter the solutions (spectra).
Returns
N/A
ACFlow.write_pole
— Functionwrite_pole(
Pᵥ::Vector{Vector{I64}},
Aᵥ::Vector{Vector{F64}},
𝕊ᵥ::Vector{Vector{F64}},
χ²ᵥ::Vector{F64},
fmesh::AbstractMesh
)
Write positions, amplitudes, and signs of poles to pole.data
. This function is only useful for the StochPX
solver.
Note that the util/ppole.jl
script can read poles' information from the pole.data
file, and construct new spectral functions with different eta
parameters.
Arguments
- Pᵥ -> Positions of the poles.
- Aᵥ -> Amplitudes of the poles.
- 𝕊ᵥ -> Signs of the poles.
- χ²ᵥ -> Goodness-of-fit functionals for all the solutions.
- fmesh -> A dense mesh for the poles.
Returns
N/A
ACFlow.write_probability
— Functionwrite_probability(α_vec::Vector{F64}, p_vec::Vector{F64})
Write p(α)
data to prob.data
. This function is only useful for the MaxEnt
solver (bryan
algorithm).
Arguments
- α_vec -> List for α parameters.
- p_vec -> α-dependent probabilities.
Returns
N/A
ACFlow.write_statistics
— Function