About the Package

SmoQyDEAC.jl utilizes the Differential Evolution for Analytic Continuation algorithm developed by Nathan S. Nichols, Paul Sokol, and Adrian Del Maestro arXiv:2201.04155.

This package takes imaginary time or Matsubara frequency correlation functions from condensed matter Monte Carlo simulations and provides the associated spectral function on the real axis.

Installation

NOTE: This package is in the experimental phase of development and is not yet published to the Julia General registry.

Open a Julia REPL environment and run the following commands:

julia> ]
pkg> add SmoQyDEAC

Running SmoQyDEAC

SmoQyDEAC has a simple API interface with a two callable functions DEAC_Binned and DEAC_Std for binned data and data with the standard error, respectively.

API

SmoQyDEAC.DEAC_BinnedFunction
DEAC_Binned(correlation_function::AbstractMatrix,
     β::Float64,
     input_grid::Vector{Float64},
     out_ωs::Vector{Float64},
     kernel_type::String,
     num_bins::Int64,
     runs_per_bin::Int64,
     output_file::String,
     checkpoint_file::String;

     find_fitness_floor::Bool = false,
     population_size::Int64 = 8,
     base_seed::Integer = 8675309,
     number_of_generations::Int64 = 1000000,
     keep_bin_data = true,
     autoresume_from_checkpoint = false,
     verbose::Bool = false,
     fitness = 1.0,
     
     crossover_probability::Float64 = 0.9,
     self_adapting_crossover_probability::Float64 = 0.1,
     differential_weight::Float64 = 0.9,
     self_adapting_differential_weight_probability::Float64 = 0.1,
     self_adapting_differential_weight::Float64 = 0.9,
     user_mutation! = nothing,
     eigenvalue_ratio_min::Float64 = 1e-8
     )

Runs the DEAC algorithm on data passed in correlation_function using $\Chi^2$ fitting using the eigenvalues of the covariance matrix

Arguments

  • correlation_function::AbstractMatrix: Input data in τ/ωₙ space, shape [Bins,τ/ωₙ]
  • β::Float64: Inverse temperature
  • input_grid::Vector{Float64}: Evenly spaced values in τ from 0 to β, including end points
  • out_ωs::Vector{Float64}: Evenly spaced energies for AC output
  • kernel_type::String: See below for allowable kernels
  • num_bins::Int64: Bins to generate
  • runs_per_bin::Int64: Number of runs per bin for statistics
  • output_file::String: File to store output dictionary. jld2 format recommended
  • checkpoint_file::String: File to store checkpoint data. jld2 format recommended

Optional Arguments

  • find_fitness_floor::Bool = false: Find the floor for fitness and append value to target fitnesses
  • population_size::Int64 = 8: DEAC population size. Must be ≥ 6
  • base_seed::Int64 = 8675309: Seed
  • number_of_generations::Int64 = 100000: Maximum number of mutation loops
  • keep_bin_data::Bool = true: Save binned data or not
  • autoresume_from_checkpoint::Bool = false: Resume from checkpoint if possible
  • verbose::Bool = false: Print stats per run
  • fitness = 1.0: single value or an array of values for χ² fitting

Optional algorithm arguments

  • crossover_probability::Float64 = 0.9: Starting likelihood of crossover
  • self_adapting_crossover_probability::Float64 = 0.1: Chance of crossover probability changing
  • differential_weight::Float64 = 0.9: Weight for second and third mutable indices
  • self_adapting_differential_weight_probability::Float64 = 0.1: Likelihood of SAD changing
  • self_adapting_differential_weight::Float64 = 0.9: SAD
  • user_mutation! = nothing: User passed function to add additional mutation to each iteration. See below for more information
  • eigenvalue_ratio_min::Float64 = 1e-8: Cutoff to mask out eigenvectors with eigenvalues ≈ 0.0 that arise due to symmetries in input data

Each run will use its own seed. E.g. if you run 10 bins with 100 runs per bin, you will use seeds base_seed:base_seed+999. You may increment your base seed by 1000, use another output file name, and generate more statistics later.

Output

SmoQyDEAC will save a dictionary to the file name passed via the output_file argument. The same dictionary will be returned by the function.

Checkpointing

SmoQyDEAC will periodically save a file named according to the parameter passed incheckpoint_file. JLD2 format is recommended, e.g. a file with the .jld2 extension. After completing every bin this file will be saved. After the last bin is finished the file will be deleted. When autoresume_from_checkpoint=true SmoQyDEAC will attempt to resume from the checkpoint. If the arguments passed do not match those in the checkpoint the code will exit.

source
SmoQyDEAC.DEAC_StdFunction
DEAC_Std(correlation_function::AbstractVector,
     correlation_function_error::AbstractVector,
     β::Float64,
     input_grid::Vector{Float64},
     out_ωs::Vector{Float64},
     kernel_type::String,
     num_bins::Int64,
     runs_per_bin::Int64,
     output_file::String,
     checkpoint_file::String;

     find_fitness_floor::Bool=true
     population_size::Int64=8,
     base_seed::Integer=8675309,
     number_of_generations::Int64=1000000,
     keep_bin_data=true,
     autoresume_from_checkpoint=false,
     verbose::Bool=false,
     
     fitness = 1.0,
     crossover_probability::Float64 = 0.9,
     self_adapting_crossover_probability::Float64 = 0.1,
     differential_weight::Float64 = 0.9,
     self_adapting_differential_weight_probability::Float64 = 0.1,
     self_adapting_differential_weight::Float64 = 0.9,
     user_mutation! = nothing 
     )

Runs the DEAC algorithm on data passed in correlation_function using $\Chi^2$ fitting from the error passed in by correlation_function_error.

Arguments

  • correlation_function::AbstractVector: Input data in τ/ωₙ space
  • correlation_function_error::AbstractVector: Error associated with input data
  • β::Float64: Inverse temperature
  • input_grid::Vector{Float64}: Evenly spaced values in τ from 0 to β, including end points or ωₙ
  • out_ωs::Vector{Float64}: Evenly spaced energies for AC output
  • kernel_type::String: See below for allowable kernels
  • num_bins::Int64: Bins to generate
  • runs_per_bin::Int64: Number of runs per bin for statistics
  • output_file::String: File to store output dictionary. jld2 format recommended
  • checkpoint_file::String: File to store checkpoint data. jld2 format recommended

Optional Arguments

  • find_fitness_floor::Bool = false: Find a reasonable floor for fitness and append that value to target fitnesses.
  • population_size::Int64 = 8: DEAC population size. Must be ≥ 6
  • base_seed::Int64 = 8675309: Seed
  • number_of_generations::Int64 = 100000: Maximum number of mutation loops
  • keep_bin_data::Bool = true: Save binned data or not
  • autoresume_from_checkpoint::Bool = false: Resume from checkpoint if possible
  • verbose::Bool = false: Print stats per run
  • fitness = 1.0: single value or an array of values for χ² fitting

Optional algorithm arguments

  • crossover_probability::Float64 = 0.9: Starting likelihood of crossover
  • self_adapting_crossover_probability::Float64 = 0.1: Chance of crossover probability changing
  • differential_weight::Float64 = 0.9: Weight for second and third mutable indices
  • self_adapting_differential_weight_probability::Float64 = 0.1: Likelihood of SAD changing
  • self_adapting_differential_weight::Float64 = 0.9: SAD
  • user_mutation! = nothing: User passed function to add additional mutation to each iteration. See below for more information

Each run will use its own seed. E.g. if you run 10 bins with 100 runs per bin, you will use seeds base_seed:base_seed+999. You may increment your base seed by 1000, use another output file name, and generate more statistics later.

Output

SmoQyDEAC will save a dictionary to the file name passed via the output_file argument. The same dictionary will be returned by the function.

Checkpointing

SmoQyDEAC will periodically save a file named according to the parameter passed incheckpoint_file. JLD2 format is recommended, e.g. a file with the .jld2 extension. After completing every bin this file will be saved. After the last bin is finished the file will be deleted. When autoresume_from_checkpoint=true SmoQyDEAC will attempt to resume from the checkpoint. If the arguments passed do not match those in the checkpoint the code will exit.

source

To add additional mutations to the base DEAC algorithm you utilize the user_mutation! functionality. user_mutation! is a user defined function which takes three (3) parameters, (population_new, population_old, rng). The result in population_new is treated as a new trial population, and if the fitness improves for a population, the trial population is stored.

population_new and population_old are of the shape [1:nω, 1:population_size]. Note, there is no guarantee that population_new is a copy of population_old, so the user should update the entire array. See Example 3: User Mutation for an example implementation.

Output

Default Keys

Both API functions return a Dict{String,Any} object as well as save that dictionary to the location specified in the parameter output_file. The dictionary has the following keys

  • A: 2D array of shape (n$\omega$,nFitness), where nFitness is the number of fitness targets in the run.
  • fitness: 1D array with all fitnesses associated with the run in descending order
  • σ: The calculated standard error for the run. NOTE: DEAC is non-ergodic, so this does not correspond actual error bars!
  • ωs: The $\omega$ values corresponding to the first dimenson of A.
  • avg_generations: number of generations to reach lowest target fitness
  • runtime: Total run time in seconds
  • zeroth_moment: 1D array of the calculated zeroth moments for the reported range of $\omega$s for each fitness
  • zeroth_moment_σ: Standard error for zeroth_moment.
  • full_eigenvalues: Only pertinent for binned data. When finding the eigenbasis for the covariance matrix near-zero eigenvalues may arise due to linear correlations in your data. Those below the value eigenvalueratiomin*maximum(eigenvalues) will be ignored and their eigenvectors not used.

Bin data

If keep_bin_data==true then the following keys are also in the output dictionary

  • bin_data: 3D array of size (n$\omega$,nBins,nFitness) with data from each bin.
  • bin_zeroth_moment: Same as zeroth_moment but on a bin by bin basis.

Kernels

The following are the supported kernels

  • time_fermionic$=\frac{e^{-\tau\omega}}{1+e^{-\beta\omega}}$
  • time_bosonic$=\frac{e^{-\tau\omega}}{1-e^{-\beta\omega}}$
  • time_bosonic_symmetric$=\frac{1}{2}\frac{e^{-\tau\omega}+e^{-(\beta-\tau)\omega}}{1-e^{-\beta\omega}}$
  • time_bosonic_symmetric_w$=\frac{\omega}{2}\frac{e^{-\tau\omega}+e^{-(\beta-\tau)\omega}}{1-e^{-\beta\omega}}$
  • frequency_fermionic$=\frac{1}{i\omega_n-\omega}$
  • frequency_bosonic$=\frac{1}{i\omega_n-\omega}$
  • frequency_bosonic_symmetric$=\frac{2 \omega}{\omega_n^2+\omega^2}$

Multithreading

SmoQyDEAC utilizes Julia's Threads.@threads multithreading capability. To take advantage of this run your code using $ julia --threads=auto yourscript.jl. auto will automatically use any available cores or hyperthreads. You can set the value to a fixed number as you wish.

Tips, tricks and caveats

  • You will likely never need to adjust any of the Optional Algorithm Arguments from the default. Many are merely initial values and will be updated stochastically early and often in the code.
  • The DEAC algorithm can have edge effects where it places spectral weight on the first or last ω point. This occurs when there is spectral weight just outside of your range of ωs. The solution is simply expanding the range of your output energies.
  • For bosonic correlations SmoQyDEAC returns the spectral function, e.g. $\mathcal{A}(\omega)$ not $\frac{\mathcal{A}(\omega)}{\omega}$ as some MaxEnt codes do.
  • Different simulation codes may report correlation functions slightly differently. E.g. for SmoQyDQMC phonon_greens $=\langle X(\tau)X(0)\rangle$ not the actual phonon green's function of $-2\Omega_0\langle X(\tau)X(0)\rangle$. While the negative sign will cancel out by our choice of Kernel, you may need to postprocess the spectral function you recover. In this case $\mathcal{A}\rightarrow \dfrac{\mathcal{A}}{2\Omega_0}$.

Funding

This work was supported by the U.S. Department of Energy, Office of Science, Office of Basic Energy Sciences, under Award Number DE-SC0022311. N.S.N. was supported by the Argonne Leadership Computing Facility, which is a U.S. Department of Energy Office of Science User Facility operated under contract DE-AC02-06CH11357.

Citation

If you found this library to be useful in the course of academic work, please consider citing us:

@Article{10.21468/SciPostPhysCodeb.39,
	title={{SmoQyDEAC.jl: A differential evolution package for the analytic continuation of imaginary time correlation functions}},
	author={James Neuhaus and Nathan S. Nichols and Debshikha Banerjee and Benjamin Cohen-Stead and Thomas A. Maier and Adrian Del Maestro and Steven Johnston},
	journal={SciPost Phys. Codebases},
	pages={39},
	year={2024},
	publisher={SciPost},
	doi={10.21468/SciPostPhysCodeb.39},
	url={https://scipost.org/10.21468/SciPostPhysCodeb.39},
}

Contact Us

For question and comments regarding this package, please email either James Neuhaus at jneuhau1@utk.edu or Professor Steven Johnston at sjohn145@utk.edu.

Publication List

*