Public interfaces
Documentation for SpinMonteCarlo.jl
's public interface (exported).
Driver
SpinMonteCarlo.runMC
— FunctionrunMC(param::Parameter)
runMC(params::AbstractArray{Parameter}
;
parallel::Bool=false,
autoID::Bool=true)
Runs Monte Carlo simulation(s) and returns calculated observables.
Restart
If a checkpoint file named "$(param["Checkpoint Filename Prefix"])_$(param["ID"]).dat"
exists and param["Checkpoint Interval"] > 0.0
, runMC
loads this file and restarts the pending simulation. NOTE: Restart will fail if the version or the system image of julia change (see the doc of Serialization.serialize
).
Keyward arguments
autoID
: If true,"ID"
s will be set (overwritten) asparams[i]["ID"] = i
.parallel
: If true, runs simulations in parallel (usespmap
instead ofmap
).
Required keys in param
- "Model"
- "Lattice"
- "Update Method"
Optional keys in param
- "MCS": The number of Monte Carlo steps after thermalization
- Default:
8192
- Default:
- "Thermalization": The number of Monte Carlo steps for thermalization
- Default:
MCS>>3
- Default:
- "Binning Size": The size of binning
- Default:
0
- Default:
- "Number of Bins": The number of bins
- Default:
0
- If both "Binning Size" and "Number of Bins" are not given, "Binning Size" is set to
floor(sqrt(MCS))
.
- Default:
- "Seed": The initial seed of the random number generator,
MersenneTwister
- Default: determined randomly (see the doc of
Random.seed!
)
- Default: determined randomly (see the doc of
- "Checkpoint Filename Prefix": See the "Restart" section.
- Default:
"cp"
- Default:
- "ID": See the "Restart" section.
- Default:
0
- Default:
- "Checkpoint Interval": Time interval between writing checkpoint file in seconds.
- Default:
0.0
, this means that NO checkpoint file will be loaded and saved.
- Default:
Lattice
PDMats.dim
— Functiondim(lat::Lattice)
dim(model::Model)
Returns the dimension of lattice.
Base.size
— Functionsize(lat::Lattice, [dim::Integer])
size(model::Model, [dim::Integer])
Returns the size of lattice.
SpinMonteCarlo.sites
— Functionsites(lat::Lattice, [sitetype::Integer])
sites(model::Model, [sitetype::Integer])
Returns an iterator over sites with sitetype
(if omitted, over all sites)
SpinMonteCarlo.bonds
— Functionbonds(lat::Lattice, [bondtype::Integer])
bonds(model::Model, [bondtype::Integer])
Returns an iterator over bonds with bondtype
(if omitted, over all bonds)
SpinMonteCarlo.numsites
— Functionnumsites(lat::Lattice)
numsites(model::Model)
Returns the number of all sites.
numsites(lat::Lattice, sitetype::Integer)
numsites(model::Model, sitetype::Integer)
Returns the number of sitetype
sites.
SpinMonteCarlo.numbonds
— Functionnumbonds(lat::Lattice)
numbonds(model::Model)
Returns the number of all bonds.
numbonds(lat::Lattice, bondtype::Integer)
numbonds(model::Model, bondtype::Integer)
Returns the number of bondtype
bonds.
SpinMonteCarlo.neighbors
— Functionneighbors(site::Site)
neighbors(lat::Lattice, site::Integer)
neighbors(model::Model, site::Integer)
Returns the index pairs of neighbor sites and bonds of site
.
SpinMonteCarlo.neighborsites
— Functionneighborsites(site::Site)
neighborsites(lat::Lattice, site::Integer)
neighborsites(model::Model, site::Integer)
Returns the indecies of neighbor sites of site
.
SpinMonteCarlo.neighborbonds
— Functionneighborbonds(site::Site)
neighborbonds(lat::Lattice, site::Integer)
neighborbonds(model::Model, site::Integer)
Returns the indecies of neighbor bonds of site
.
SpinMonteCarlo.source
— Functionsource(bond::Bond)
source(lat::Lattice, bond::Integer)
source(model::Model, bond::Integer)
Returns the source site index of bond
.
SpinMonteCarlo.target
— Functiontarget(bond::Bond)
target(lat::Lattice, bond::Integer)
target(model::Model, bond::Integer)
Returns the target site of bond
.
SpinMonteCarlo.sitetype
— Functionsitetype(site::Site)
sitetype(lat::Lattice, site::Integer)
sitetype(model::Model, site::Integer)
Returns the type of site
SpinMonteCarlo.bondtype
— Functionbondtype(bond::Bond)
bondtype(lat::Lattice, bond::Integer)
bondtype(model::Model, bond::Integer)
Returns the type of bond
.
SpinMonteCarlo.sitecoordinate
— Functionsitecoordinate(site::Site)
sitecoordinate(lat::Lattice, site::Integer)
sitecoordinate(model::Model, site::Integer)
Returns the coordinate of the site
in the Cartesian system
SpinMonteCarlo.bonddirection
— Functionbonddirection(bond::Bond)
bonddirection(lat::Lattice, bond::Integer)
bonddirection(model::Model, bond::Integer)
Returns the unnormalized direction of the bond
as vector in the Cartesian system
SpinMonteCarlo.cellcoordinate
— Functioncellcoordinate(site::Site)
cellcoordinate(lat::Lattice, site::Integer)
cellcoordinate(model::Model, site::Integer)
Returns the coordinate of the cell
including site
in the Lattice system
Model
SpinMonteCarlo.Ising
— TypeIsing model with energy $E = -\sum_{ij} J_{ij} \sigma_i \sigma_j$, where $\sigma_i$ takes value of 1 (up spin) or -1 (down spin).
SpinMonteCarlo.Potts
— TypeQ
state Potts model with energy $E = -\sum_{i,j} \delta_{\sigma_i, \sigma_j}$, where $\sigma_i$ takes an integer value from $1$ to $Q$ and $\delta$ is a Kronecker's delta. Order parameter (total magnetization) is defined as \begin{equation} M = \frac{Q-1}{Q}N1 - \frac{1}{Q}(N-N1), \end{equation} where $N$ is the number of sites and $N_1$ is the number of $\sigma=1$ spins.
SpinMonteCarlo.Clock
— TypeQ
state clock model with energy $E = -\sum_{ij} J_{ij} \cos(\theta_i - \theta_j)$, where $\theta_i = 2\pi \sigma_i/Q$ and $\sigma_i$ takes an integer value from $1$ to $Q$.
SpinMonteCarlo.XY
— TypeXY model with energy $E = -\sum_{ij} J_{ij} \cos(\theta_i - \theta_j)$, where $\theta_i = 2\pi \sigma_i$ and $\sigma_i \in [0, 1)$.
SpinMonteCarlo.AshkinTeller
— TypeAshkinTeller model with energy $E = -\sum_{ij} J^\sigma_{ij} \sigma_i \sigma_j + J^\tau_{ij} \tau_i \tau_j + K_{ij} \sigma_i \sigma_j \tau_i \tau_j$, where $\sigma_i$ and $\tau_i$ takes value of 1 (up spin) or -1 (down spin).
SpinMonteCarlo.QuantumXXZ
— TypeSpin-$S$ XXZ model represented as the following Hamiltonian,
$\mathcal{H} = \sum_{i,j} \left[ J_{ij}^z S_i^z S_j^z + \frac{J_{ij}^{xy}}{2} (S_i^+ S_j^- + S_i^-S_j^+) \right] - \sum_i \Gamma_i S_i^x,$
where $S^x, S^y, S^z$ are $x, y$ and $z$ component of spin operator with length $S$, and $S^\pm \equiv S^x \pm iS^y$ are ladder operator. A state is represented by a product state (spins at $\tau=0$) of local $S^z$ diagonal basis and an operator string (perturbations). A local spin with length $S$ is represented by a symmetrical summation of $2S$ sub spins with length $1/2$.
Update method
An index of model parameter (e.g., Js
) is corresponding to sitetype
or bondtype
.
SpinMonteCarlo.local_update!
— Functionlocal_update!(model, param)
Updates spin configuration by local spin flip and Metropolice algorithm
SpinMonteCarlo.SW_update!
— FunctionSW_update!(model, param::Parameter)
Updates spin configuration by Swendsen-Wang algorithm
SpinMonteCarlo.Wolff_update!
— FunctionWolff_update!(model, param::Parameter)
Updates spin configuration by Wolff algorithm
SpinMonteCarlo.loop_update!
— Functionloop_update!(model, param::Parameter)
Updates spin configuration by loop algorithm
loop_update!(model, param::Parameter)
loop_update!(model, T::Real,
Jz::AbstractArray,
Jxy::AbstractArray,
Gamma:AbstractArray)
Updates spin configuration by loop algorithm under the temperature T = param["T"]
and coupling constants Jz, Jxy
and transverse field Gamma
Estimator
SpinMonteCarlo.simple_estimator
— Functionsimple_estimator(model::Ising, T::Real, Js::AbstractArray)
Returns the following observables as Dict{String, Any}
Observables
"Energy"
- energy density
"Energy^2"
- square of energy density
"Magnetization"
- magnetization density
"|Magnetization|"
- absolute value of magnetization density
"Magnetization^2"
- square of magnetization density
"Magnetization^4"
- quadruple of magnetization density
simple_estimator(model::Clock, T::Real, Js::AbstractArray)
simple_estimator(model::XY, T::Real, Js::AbstractArray)
Returns the following observables as Dict{String, Any}
Observables
"Energy"
- Energy per spin (site)
"Energy^2"
"|Magnetization|"
- Absolute value of total magnetization per spin (order paremeter)
"|Magnetization|^2"
"|Magnetization|^4"
"Magnetization x"
- x component of total magnetization per spin (order paremeter)
"|Magnetization x|"
"Magnetization x^2"
"Magnetization x^4"
"Magnetization y"
- y component of total magnetization per spin (order paremeter)
"|Magnetization y|"
"Magnetization y^2"
"Magnetization y^4"
"Helicity Modulus x"
"Helicity Modulus y"
simple_estimator(model::AshkinTeller, T::Real, Jsigma, Jtau, K)
Returns the following observables as Dict{String, Any}
Observables
"Energy"
- energy density
"Energy^2"
- square of energy density
"|Magnetization|"
- absolute value of magnetization density, $|m| = \sqrt{ (\sum_i \sigma_i )^2 + (\sum_i \tau_i)^2 } / N$
"|Magnetization|^2"
- square of magnetization density
"|Magnetization|^4"
- quadruple of magnetization density
"Magnetization sigma"
- magnetization density (sigma spin)
"|Magnetization sigma|"
"Magnetization sigma^2"
"Magnetization sigma^4"
"Magnetization tau"
- magnetization density (tau spin)
"|Magnetization tau|"
"Magnetization tau^2"
"Magnetization tau^4"
SpinMonteCarlo.improved_estimator
— Functionimproved_estimator(model::Ising, T::Real, Js::AbstractArray, sw::SWInfo)
Returns the following observables as Dict{String, Any}
using cluster information sw
Observables
"Energy"
- energy density
"Energy^2"
- square of energy density
"Magnetization"
- magnetization density
"|Magnetization|"
- absolute value of magnetization density
"Magnetization^2"
- square of magnetization density
"Magnetization^4"
- quadruple of magnetization density
"Clustersize^2"
- $\sum_c r_c^2$, where $r_c$ is the size density of $c$-th cluster
"Clustersize^4"
- $\sum_c r_c^4$
"Clustersize^2 Clustersize^2"
- $\sum_{c\ne c'} r_c^2 r_{c'}^2$
improved_estimator(model::Potts, T::Real, Js::AbstractArray, sw::SWInfo)
Returns the following observables as Dict{String, Any}
using cluster information sw
Observables
"Energy"
- energy density
"Energy^2"
- square of energy density
"Magnetization"
- magnetization density
"|Magnetization|"
- absolute value of magnetization density
"|Magnetization|^2"
- square of magnetization density
"|Magnetization|^4"
- quadruple of magnetization density
"Clustersize^2"
- $\sum_c r_c^2$, where $r_c$ is the size density of $c$-th cluster
"Clustersize^4"
- $\sum_c r_c^4$
"Clustersize^2 Clustersize^2"
- $\sum_{c\ne c'} r_c^2 r_{c'}^2$
improved_estimator(model::QuantumXXZ, T::Real, Js::AbstractArray, uf::UnionFind)
Returns the following observables as Dict{String, Any}
using loop information uf
Observables
"Sign"
- Sign of the weight function
"Sign * Energy"
- Energy per spin (site)
"Sign * Energy^2"
"Sign * Magnetization"
- Total magnetization (Sz) per spin (site)
"Sign * |Magnetization|"
"Sign * Magnetization^2"
"Sign * Magnetization^4"
Observables
SpinMonteCarlo.MCObservable
— TypeMCObservable
abstract type representing observable in Monte Carlo calculation.
SpinMonteCarlo.ScalarObservable
— TypeScalarObservable <: MCObservable
abstract type representing scalar-type observable in Monte Carlo calculation.
SpinMonteCarlo.VectorObservable
— TypeVectorObservable <: MCObservable
abstract type representing vector-type observable in Monte Carlo calculation.
SpinMonteCarlo.MCObservableSet
— TypeMCObservableSet{Obs}
Alias of `Dict{String, Obs}` where `Obs` is a subtype of `MCObservable`.
SpinMonteCarlo.makeMCObservable!
— FunctionmakeMCObservable!(oset::MCObservableSet{Obs}, name::String)
Create an observable with the name `name` in the observable set `oset`.
SpinMonteCarlo.SimpleObservable
— TypeSimpleObservable
A simple observable which stores all the data in memory.
SpinMonteCarlo.SimpleVectorObservable
— TypeSimpleVectorObservable
A simple observable which stores all the data in memory.
SpinMonteCarlo.SimpleObservableSet
— TypeSimpleObservableSet
Alias of `MCObservableSet{SimpleObservable}`.
SpinMonteCarlo.SimpleVectorObservableSet
— TypeSimpleObservableSet
Alias of `MCObservableSet{SimpleVectorObservable}`.
SpinMonteCarlo.binning
— Functionbinning(obs::SimpleObservable; binsize::Int = 0, numbins::Int = 0)
binning(obs::SimpleVectorObservable; binsize::Int = 0, numbins::Int = 0)
Binning the observable `obs`.
Either `binsize` or `numbins` can be given. If both are given, `ArgumentError` is thrown.
If both are not given, `binsize` is set to `floor(sqrt(count(obs)))`.
SpinMonteCarlo.Jackknife
— TypeJackknife <: ScalarObservable
Jackknife resampling observable.
SpinMonteCarlo.JackknifeVector
— TypeJackknifeVector <: VectorObservable
Jackknife resampling observable for vector-type observable.
SpinMonteCarlo.JackknifeSet
— TypeJackknifeSet
Alias of `MCObservableSet{Jackknife}`.
SpinMonteCarlo.JackknifeVectorSet
— TypeJackknifeVectorSet
Alias of `MCObservableSet{JackknifeVector}`.
SpinMonteCarlo.jackknife
— Functionjackknife(obs::ScalarObservable)
Construct a Jackknife observable from a scalar observable
jackknife(obsset::MCObservableSet)
Construct a JackknifeSet from a MCObservableSet
jackknife(obs::VectorObservable)
Construct a JackknifeVector observable from a vector observable
Snapshots
SpinMonteCarlo.gen_snapshot!
— Functiongen_snapshot!(model, T, [N=1])
generate and return N
snapshots (spin configuration).
Arguments
model::Model
: modelT::Real
: temperatureN::Integer=1
: the number of snapshots to be generatedMCS::Integer=8192
: the number of Monte Carlo steps followed by snapshot generation
SpinMonteCarlo.gensave_snapshot!
— Functiongensave_snapshot!(io, model, T, [N=1])
gensave_snapshot!(filename, model, T, [N=1])
generate and write N
snapshots into io
or filename
.
Arguments
io::IO
: output stream where snapshots will be writtenfilename::String
: the name of file where snapshots will be writtenmodel::Model
: model to be simulatedT::Real
: temperatureN::Integer=1
: the number of snapshots to be generatedMCS::Integer=1
: the number of Monte Carlo steps followed by snapshot generationsep::String=" "
: field separator of snapshot
SpinMonteCarlo.load_snapshot
— Functionload_snapshot(io, splitter)
load_snapshot(filename, splitter)
load snapshot from io
or filename
. splitter
is the set of field splitter and will be passed Base.split
as the second argument (see the doc of Base.split
).
Utility
SpinMonteCarlo.Parameter
— TypeInput parameter of simulation
SpinMonteCarlo.convert_parameter
— Functionconvert_parameter(model, param)
Generates arguments of updater and estimator.
Example
julia> model = Ising(chain_lattice(4));
julia> p = convert_parameter(model, Parameter("J"=>1.0))
(1.0, [1.0, 1.0]) # T and Js
julia> p = convert_parameter(model, Parameter("J"=>[1.5, 0.5]))
(1.0, [1.5, 0.5]) # J can take a vector whose size is `numbondtypes(model)`
julia> model.spins
4-element Array{Int64,1}:
1
1
1
-1
julia> local_update!(model, p...);
julia> model.spins
4-element Array{Int64,1}:
1
1
-1
-1
SpinMonteCarlo.convert_parameter
— Methodconvert_parameter(model::Ising, param::Parameter)
Keynames:
- "T": a scalar (default: 1.0).
- "J": a vector with
numbondtypes(model)
elements (default: 1.0).
SpinMonteCarlo.convert_parameter
— Methodconvert_parameter(model::Potts, param::Parameter)
Keynames:
- "T": a scalar (default: 1.0).
- "J": a vector with
numbondtypes(model)
elements (default: 1.0).
SpinMonteCarlo.convert_parameter
— Methodconvert_parameter(model::Clock, param::Parameter)
Keynames:
- "T": a scalar (default: 1.0).
- "J": a vector with
numbondtypes(model)
elements (default: 1.0).
SpinMonteCarlo.convert_parameter
— Methodconvert_parameter(model::XY, param::Parameter)
Keynames:
- "T": a scalar (default: 1.0).
- "J": a vector with
numbondtypes(model)
elements (default: 1.0).
SpinMonteCarlo.convert_parameter
— Methodconvert_parameter(model::AshkinTeller, param::Parameter)
Keynames:
- "T": a scalar (default: 1.0).
- "Jsigma": a vector with
numbondtypes(model)
elements (default: 1.0). - "Jtau": a vector with
numbondtypes(model)
elements (default: 1.0). - "K": a vector with
numbondtypes(model)
elements (default: 0.0).
SpinMonteCarlo.convert_parameter
— Methodconvert_parameter(model::QuantumXXZ, param::Parameter)
Keynames:
- "T": a scalar (default: 1.0).
- "Jz": a vector with
numbondtypes(model)
elements (default: 1.0). - "Jxy": a vector with
numbondtypes(model)
elements (default: 1.0). - "Gamma": a vector with
numsitetypes(model)
elements (default: 0.0).