Public interfaces

Documentation for SpinMonteCarlo.jl's public interface (exported).

Driver

SpinMonteCarlo.runMCFunction
runMC(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) as params[i]["ID"] = i.
  • parallel: If true, runs simulations in parallel (uses pmap instead of map).

Required keys in param

  • "Model"
  • "Lattice"
  • "Update Method"

Optional keys in param

  • "MCS": The number of Monte Carlo steps after thermalization
    • Default: 8192
  • "Thermalization": The number of Monte Carlo steps for thermalization
    • Default: MCS>>3
  • "Seed": The initial seed of the random number generator, MersenneTwister
    • Default: determined randomly (see the doc of Random.seed!)
  • "Checkpoint Filename Prefix": See the "Restart" section.
    • Default: "cp"
  • "ID": See the "Restart" section.
    • Default: 0
  • "Checkpoint Interval": Time interval between writing checkpoint file in seconds.
    • Default: 0.0, this means that NO checkpoint file will be loaded and saved.
source

Lattice

PDMats.dimFunction
dim(lat::Lattice)
dim(model::Model)

Returns the dimension of lattice.

source
Base.sizeFunction
size(lat::Lattice, [dim::Integer])
size(model::Model, [dim::Integer])

Returns the size of lattice.

source
SpinMonteCarlo.sitesFunction
sites(lat::Lattice, [sitetype::Integer])
sites(model::Model, [sitetype::Integer])

Returns an iterator over sites with sitetype (if omitted, over all sites)

source
SpinMonteCarlo.bondsFunction
bonds(lat::Lattice, [bondtype::Integer])
bonds(model::Model, [bondtype::Integer])

Returns an iterator over bonds with bondtype (if omitted, over all bonds)

source
SpinMonteCarlo.numsitesFunction
numsites(lat::Lattice)
numsites(model::Model)

Returns the number of all sites.

source
numsites(lat::Lattice, sitetype::Integer)
numsites(model::Model, sitetype::Integer)

Returns the number of sitetype sites.

source
SpinMonteCarlo.numbondsFunction
numbonds(lat::Lattice)
numbonds(model::Model)

Returns the number of all bonds.

source
numbonds(lat::Lattice, bondtype::Integer)
numbonds(model::Model, bondtype::Integer)

Returns the number of bondtype bonds.

source
SpinMonteCarlo.neighborsFunction
neighbors(site::Site)
neighbors(lat::Lattice, site::Integer)
neighbors(model::Model, site::Integer)

Returns the index pairs of neighbor sites and bonds of site.

source
SpinMonteCarlo.neighborsitesFunction
neighborsites(site::Site)
neighborsites(lat::Lattice, site::Integer)
neighborsites(model::Model, site::Integer)

Returns the indecies of neighbor sites of site.

source
SpinMonteCarlo.neighborbondsFunction
neighborbonds(site::Site)
neighborbonds(lat::Lattice, site::Integer)
neighborbonds(model::Model, site::Integer)

Returns the indecies of neighbor bonds of site.

source
SpinMonteCarlo.sourceFunction
source(bond::Bond)
source(lat::Lattice, bond::Integer)
source(model::Model, bond::Integer)

Returns the source site index of bond.

source
SpinMonteCarlo.targetFunction
target(bond::Bond)
target(lat::Lattice, bond::Integer)
target(model::Model, bond::Integer)

Returns the target site of bond.

source
SpinMonteCarlo.sitetypeFunction
sitetype(site::Site)
sitetype(lat::Lattice, site::Integer)
sitetype(model::Model, site::Integer)

Returns the type of site

source
SpinMonteCarlo.bondtypeFunction
bondtype(bond::Bond)
bondtype(lat::Lattice, bond::Integer)
bondtype(model::Model, bond::Integer)

Returns the type of bond.

source
SpinMonteCarlo.sitecoordinateFunction
sitecoordinate(site::Site)
sitecoordinate(lat::Lattice, site::Integer)
sitecoordinate(model::Model, site::Integer)

Returns the coordinate of the site in the Cartesian system

source
SpinMonteCarlo.bonddirectionFunction
bonddirection(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

source
SpinMonteCarlo.cellcoordinateFunction
cellcoordinate(site::Site)
cellcoordinate(lat::Lattice, site::Integer)
cellcoordinate(model::Model, site::Integer)

Returns the coordinate of the cell including site in the Lattice system

source

Model

SpinMonteCarlo.IsingType

Ising 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).

source
SpinMonteCarlo.PottsType

Q 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.

source
SpinMonteCarlo.ClockType

Q 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$.

source
SpinMonteCarlo.XYType

XY 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)$.

source
SpinMonteCarlo.AshkinTellerType

AshkinTeller 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).

source
SpinMonteCarlo.QuantumXXZType

Spin-$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$.

source

Update method

An index of model parameter (e.g., Js) is corresponding to sitetype or bondtype.

SpinMonteCarlo.loop_update!Function
loop_update!(model, param::Parameter)

Updates spin configuration by loop algorithm

source
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

source

Estimator

SpinMonteCarlo.simple_estimatorFunction
simple_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
source
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"
source
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"
source
SpinMonteCarlo.improved_estimatorFunction
improved_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$
source
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$
source
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"
source

Snapshots

SpinMonteCarlo.gen_snapshot!Function
gen_snapshot!(model, T, [N=1])

generate and return N snapshots (spin configuration).

Arguments

  • model::Model : model
  • T::Real : temperature
  • N::Integer=1 : the number of snapshots to be generated
  • MCS::Integer=8192 : the number of Monte Carlo steps followed by snapshot generation
source
SpinMonteCarlo.gensave_snapshot!Function
gensave_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 written
  • filename::String : the name of file where snapshots will be written
  • model::Model : model to be simulated
  • T::Real : temperature
  • N::Integer=1 : the number of snapshots to be generated
  • MCS::Integer=1 : the number of Monte Carlo steps followed by snapshot generation
  • sep::String=" " : field separator of snapshot
source
SpinMonteCarlo.load_snapshotFunction
load_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).

source

Utility

SpinMonteCarlo.convert_parameterFunction
convert_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
source
SpinMonteCarlo.convert_parameterMethod
convert_parameter(model::Ising, param::Parameter)

Keynames:

  • "T": a scalar (default: 1.0).
  • "J": a vector with numbondtypes(model) elements (default: 1.0).
source
SpinMonteCarlo.convert_parameterMethod
convert_parameter(model::Potts, param::Parameter)

Keynames:

  • "T": a scalar (default: 1.0).
  • "J": a vector with numbondtypes(model) elements (default: 1.0).
source
SpinMonteCarlo.convert_parameterMethod
convert_parameter(model::Clock, param::Parameter)

Keynames:

  • "T": a scalar (default: 1.0).
  • "J": a vector with numbondtypes(model) elements (default: 1.0).
source
SpinMonteCarlo.convert_parameterMethod
convert_parameter(model::XY, param::Parameter)

Keynames:

  • "T": a scalar (default: 1.0).
  • "J": a vector with numbondtypes(model) elements (default: 1.0).
source
SpinMonteCarlo.convert_parameterMethod
convert_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).
source
SpinMonteCarlo.convert_parameterMethod
convert_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).
source