vibespin.models package
Submodules
models.clock_model module
2D q-state Clock Model simulation using Metropolis-Hastings and Wolff cluster algorithms.
- class models.clock_model.ClockSimulation(*, size: int, temp: float, J: float = 1.0, A: float = 1.0, q: int = 6, update: str = 'checkerboard', init_state: str = 'random', parallel: bool = False, seed: int | None = None)[source]
Bases:
MonteCarloSimulationSimulation of the 2D q-state clock model on a square lattice.
- class models.clock_model.DiscreteClockSimulation(*, size: int, temp: float, J: float = 1.0, q: int = 6, update: str = 'checkerboard', init_state: str = 'random', parallel: bool = False, seed: int | None = None)[source]
Bases:
MonteCarloSimulationSimulation of the 2D q-state clock model with discrete integer spins.
Each spin takes one of q evenly spaced orientations theta_k = 2*pi*k/q. The Hamiltonian is E = -J * sum_<i,j> cos(theta_i - theta_j).
Unlike
ClockSimulation, which uses continuous unit-vector spins plus an anisotropy term, this class represents spins as integers in {0, …, q-1} and enforces the discrete symmetry exactly.
- models.clock_model.clock_energy_numba(*, spins: ndarray, J: float, A: float, q: int, idx_next: ndarray) float[source]
Calculate the total energy of the Clock Model lattice.
- Parameters:
spins ((N, N, 2) array of unit vectors.)
J (Coupling constant.)
A (Anisotropy strength.)
q (Number of clock states.)
idx_next (Pre-calculated next-neighbor indices.)
- Returns:
energy
- Return type:
Total energy per site.
- models.clock_model.clock_step_numba(*, spins: ndarray, beta: float, J: float, A: float, q: int, idx_next: ndarray, idx_prev: ndarray) ndarray[source]
Perform one full Monte Carlo sweep of the Clock Model lattice. Uses a checkerboard update pattern for better Numba optimization.
- Parameters:
spins ((N, N, 2) array of unit vectors.)
beta (Inverse temperature 1/kT.)
J (Coupling constant.)
A (Anisotropy strength.)
q (Number of clock states.)
idx_next (Pre-calculated next-neighbor indices.)
idx_prev (Pre-calculated previous-neighbor indices.)
- Return type:
Updated spins array.
- models.clock_model.clock_step_parallel_numba(*, spins: ndarray, beta: float, J: float, A: float, q: int, idx_next: ndarray, idx_prev: ndarray) ndarray[source]
Parallel version of clock_step_numba.
- Parameters:
spins ((N, N, 2) array of unit vectors.)
beta (Inverse temperature 1/kT.)
J (Coupling constant.)
A (Anisotropy strength.)
q (Number of clock states.)
idx_next (Pre-calculated next-neighbor indices.)
idx_prev (Pre-calculated previous-neighbor indices.)
- Return type:
Updated spins array.
- models.clock_model.clock_step_random_numba(*, spins: ndarray, beta: float, J: float, A: float, q: int, idx_next: ndarray, idx_prev: ndarray) ndarray[source]
Perform one full Monte Carlo sweep of the Clock Model lattice using random sequential updates.
- Parameters:
spins ((N, N, 2) array of unit vectors.)
beta (Inverse temperature 1/kT.)
J (Coupling constant.)
A (Anisotropy strength.)
q (Number of clock states.)
idx_next (Pre-calculated next-neighbor indices.)
idx_prev (Pre-calculated previous-neighbor indices.)
- Return type:
Updated spins array.
- models.clock_model.clock_wolff_step_numba(*, spins: ndarray, beta: float, J: float, idx_next: ndarray, idx_prev: ndarray) ndarray[source]
Perform one Wolff cluster flip on the continuous Clock Model lattice.
Implements the Wolff-Evertz reflection scheme for O(2) spins using only the Heisenberg exchange coupling J. A random mirror-plane axis r̂ is sampled from S^1, and bonds are activated with probability
P_add = 1 - exp(-2 beta J sigma_i sigma_j)wheresigma_i = s_i · r̂. Cluster spins are reflected:s -> s - 2 (s · r̂) r̂.Limitation: the crystal-field anisotropy term
A cos(q phi)breaks the O(2) reflection symmetry required by the Fortuin-Kasteleyn bond construction. This kernel therefore satisfies detailed balance only for the exchange part of the Hamiltonian; use it at parameters where the anisotropy is weak relative to J, or set A=0 for pure XY-like studies.One call constitutes one cluster sweep.
parallel=Trueis silently ignored.- Parameters:
spins ((N, N, 2) array of unit vectors.)
beta (Inverse temperature 1/kT.)
J (Coupling constant.)
idx_next (Pre-calculated next-neighbor indices (PBC).)
idx_prev (Pre-calculated previous-neighbor indices (PBC).)
- Return type:
Updated spins array.
- models.clock_model.discrete_clock_energy_numba(*, spins: ndarray, J: float, q: int, cos_table: ndarray, idx_next: ndarray) float[source]
Total energy per site of the discrete clock model.
Sums -J * cos(theta_i - theta_j) over unique right/down neighbor pairs.
- Parameters:
spins ((N, N) array of integer spin states in {0, ..., q-1}.)
J (Coupling constant.)
q (Number of clock states.)
cos_table (Pre-computed cos(2*pi*d/q) for d in {0, ..., q-1}.)
idx_next (Pre-calculated next-neighbor indices.)
- Return type:
Energy per site.
- models.clock_model.discrete_clock_step_numba(*, spins: ndarray, beta: float, J: float, q: int, cos_table: ndarray, idx_next: ndarray, idx_prev: ndarray) ndarray[source]
One Metropolis sweep of the discrete clock model (checkerboard update).
- Parameters:
spins ((N, N) array of integer spin states in {0, ..., q-1}.)
beta (Inverse temperature 1/kT.)
J (Coupling constant.)
q (Number of clock states.)
cos_table (Pre-computed cos(2*pi*d/q) for d in {0, ..., q-1}.)
idx_next (Pre-calculated next-neighbor indices.)
idx_prev (Pre-calculated previous-neighbor indices.)
- Return type:
Updated spins array.
- models.clock_model.discrete_clock_step_parallel_numba(*, spins: ndarray, beta: float, J: float, q: int, cos_table: ndarray, idx_next: ndarray, idx_prev: ndarray) ndarray[source]
Parallel version of discrete_clock_step_numba.
- Parameters:
spins ((N, N) array of integer spin states in {0, ..., q-1}.)
beta (Inverse temperature 1/kT.)
J (Coupling constant.)
q (Number of clock states.)
cos_table (Pre-computed cos(2*pi*d/q) for d in {0, ..., q-1}.)
idx_next (Pre-calculated next-neighbor indices.)
idx_prev (Pre-calculated previous-neighbor indices.)
- Return type:
Updated spins array.
- models.clock_model.discrete_clock_step_random_numba(*, spins: ndarray, beta: float, J: float, q: int, cos_table: ndarray, idx_next: ndarray, idx_prev: ndarray) ndarray[source]
One Metropolis sweep of the discrete clock model (random sequential update).
N^2 single-spin flip attempts at uniformly random sites, giving physical stochastic dynamics suitable for kinetics studies.
- Parameters:
spins ((N, N) array of integer spin states in {0, ..., q-1}.)
beta (Inverse temperature 1/kT.)
J (Coupling constant.)
q (Number of clock states.)
cos_table (Pre-computed cos(2*pi*d/q) for d in {0, ..., q-1}.)
idx_next (Pre-calculated next-neighbor indices.)
idx_prev (Pre-calculated previous-neighbor indices.)
- Return type:
Updated spins array.
models.ising_model module
2D Ising Model simulation using Metropolis-Hastings and Wolff cluster algorithms.
- class models.ising_model.IsingSimulation(*, size: int, temp: float, J: float = 1.0, update: str = 'checkerboard', init_state: str = 'random', parallel: bool = False, seed: int | None = None)[source]
Bases:
MonteCarloSimulationSimulation of the 2D Ising model on a square lattice.
- run_with_cluster_sizes(*, n_steps: int) tuple[ndarray, ndarray, ndarray][source]
Run the simulation and additionally record the Wolff cluster size at every step.
For non-Wolff update schemes the cluster-size array is filled with zeros because local Metropolis flips exactly one spin at a time (cluster size = 1 by convention would also be valid, but zero makes it easy to detect misuse).
- Parameters:
n_steps – Number of MC steps to perform and record.
- Returns:
magnetization – Array of |M| per spin at each step.
energies – Array of energy per spin at each step.
cluster_sizes – Array of cluster sizes (number of spins flipped) at each step. Always zero for non-Wolff algorithms.
- models.ising_model.ising_energy_numba(*, spins: ndarray, J: float, idx_next: ndarray) float[source]
Calculate the total energy of the Ising lattice.
- Parameters:
spins ((N, N) array of spins.)
J (Coupling constant.)
idx_next (Pre-calculated next-neighbor indices.)
- Returns:
energy
- Return type:
Total energy per site.
- models.ising_model.ising_step_numba(*, spins: ndarray, beta: float, J: float, idx_next: ndarray, idx_prev: ndarray) ndarray[source]
Perform one full Monte Carlo sweep of the Ising lattice. Uses a checkerboard update pattern for better Numba optimization.
- Parameters:
spins ((N, N) array of spins (+1 or -1).)
beta (Inverse temperature 1/kT.)
J (Coupling constant.)
idx_next (Pre-calculated next-neighbor indices.)
idx_prev (Pre-calculated previous-neighbor indices.)
- Return type:
Updated spins array.
- models.ising_model.ising_step_parallel_numba(*, spins: ndarray, beta: float, J: float, idx_next: ndarray, idx_prev: ndarray) ndarray[source]
Parallel version of ising_step_numba using numba.prange.
- Parameters:
spins ((N, N) array of spins (+1 or -1).)
beta (Inverse temperature 1/kT.)
J (Coupling constant.)
idx_next (Pre-calculated next-neighbor indices.)
idx_prev (Pre-calculated previous-neighbor indices.)
- Return type:
Updated spins array.
- models.ising_model.ising_step_random_numba(*, spins: ndarray, beta: float, J: float, idx_next: ndarray, idx_prev: ndarray) ndarray[source]
Perform one full Monte Carlo sweep of the Ising lattice using random sequential updates (N^2 randomly chosen single-spin flip attempts).
Unlike the checkerboard sweep, sites are chosen uniformly at random with replacement, giving a more physical stochastic dynamics with no spatial bias. This is the standard Metropolis single-spin flip algorithm.
Boltzmann factors for the two possible positive energy changes (4J, 8J) are pre-calculated outside the inner loop to avoid repeated calls to
np.exp, matching the optimisation used inising_step_numba.- Parameters:
spins ((N, N) array of spins (+1 or -1).)
beta (Inverse temperature 1/kT.)
J (Coupling constant.)
idx_next (Pre-calculated next-neighbor indices.)
idx_prev (Pre-calculated previous-neighbor indices.)
- Return type:
Updated spins array.
- models.ising_model.ising_wolff_step_numba(*, spins: ndarray, beta: float, J: float, idx_next: ndarray, idx_prev: ndarray) tuple[source]
Perform one Wolff cluster flip on the Ising lattice.
A single cluster is grown by probabilistic DFS from a random seed site and then flipped in a single collective move. The bond acceptance probability
P_add = 1 - exp(-2 beta J)is derived from the Fortuin-Kasteleyn representation and guarantees detailed balance without rejections.One call constitutes one cluster sweep. The effective number of site updates per call scales as the mean cluster size
<C>and diverges at T_c, which is precisely where Metropolis suffers maximum critical slowing down.parallel=Trueis silently ignored: cluster growth is inherently sequential.- Parameters:
spins ((N, N) array of spins (+1 or -1).)
beta (Inverse temperature 1/kT.)
J (Coupling constant.)
idx_next (Pre-calculated next-neighbor indices (PBC).)
idx_prev (Pre-calculated previous-neighbor indices (PBC).)
- Returns:
spins (Updated spins array.)
cluster_size (Number of spins flipped in this step.)
models.simulation_base module
Base classes and shared Numba-accelerated kernels for Monte Carlo simulations.
- class models.simulation_base.MonteCarloSimulation(*, size: int, temp: float, init_state: str = 'random', seed: int | None = None)[source]
Bases:
ABCAbstract base class for 2D lattice Monte Carlo simulations. Provides infrastructure for equilibration, measurement runs, and statistical analysis.
- equilibrate(*, n_steps: int) None[source]
Perform equilibration steps without recording measurements.
- Parameters:
n_steps (Number of MC steps to perform.)
- models.simulation_base.calculate_vortex_density_numba(*, spins: ndarray, idx_next: ndarray) float[source]
Calculate vortex density n_v, i.e. the fraction of plaquettes with non-zero winding.
- Parameters:
spins ((N, N, 2) array of unit vectors.)
idx_next (Pre-calculated next-neighbor indices.)
- Return type:
Vortex density n_v in [0, 1].
- models.simulation_base.calculate_vorticity_numba(*, spins: ndarray, idx_next: ndarray) ndarray[source]
Calculate the vorticity (winding number) of each plaquette for 2D vector spins. Optimized by calculating angles first and then using a fast wrapping kernel.
- Parameters:
spins ((N, N, 2) array of unit vectors.)
idx_next (Pre-calculated next-neighbor indices.)
- Returns:
vorticity
- Return type:
(N, N) array containing winding numbers (+1, -1, or 0).
- models.simulation_base.get_helicity_data_numba(*, spins: ndarray, idx_next: ndarray) tuple[float, float][source]
Calculate sum of cos and sin of angle differences in x-direction for helicity modulus.
- Parameters:
spins ((N, N, 2) array of unit vectors.)
idx_next (Pre-calculated next-neighbor indices for PBCs.)
- Returns:
cos_sum (Sum of cosine of angle differences.)
sin_sum (Sum of sine of angle differences.)
models.xy_model module
2D XY Model simulation using Metropolis-Hastings and Wolff cluster algorithms.
- class models.xy_model.XYSimulation(*, size: int, temp: float, J: float = 1.0, update: str = 'checkerboard', init_state: str = 'random', parallel: bool = False, seed: int | None = None)[source]
Bases:
MonteCarloSimulationSimulation of the 2D XY model on a square lattice.
- models.xy_model.xy_energy_numba(*, spins: ndarray, J: float, idx_next: ndarray) float[source]
Calculate the total energy of the XY lattice.
- Parameters:
spins ((N, N, 2) array of unit vectors.)
J (Coupling constant.)
idx_next (Pre-calculated next-neighbor indices.)
- Returns:
energy
- Return type:
Total energy per site.
- models.xy_model.xy_step_numba(*, spins: ndarray, beta: float, J: float, idx_next: ndarray, idx_prev: ndarray) ndarray[source]
Perform one full Monte Carlo sweep of the XY lattice. Uses a checkerboard update pattern for better Numba optimization.
- Parameters:
spins ((N, N, 2) array of unit vectors.)
beta (Inverse temperature 1/kT.)
J (Coupling constant.)
idx_next (Pre-calculated next-neighbor indices.)
idx_prev (Pre-calculated previous-neighbor indices.)
- Return type:
Updated spins array.
- models.xy_model.xy_step_parallel_numba(*, spins: ndarray, beta: float, J: float, idx_next: ndarray, idx_prev: ndarray) ndarray[source]
Parallel version of xy_step_numba.
- Parameters:
spins ((N, N, 2) array of unit vectors.)
beta (Inverse temperature 1/kT.)
J (Coupling constant.)
idx_next (Pre-calculated next-neighbor indices.)
idx_prev (Pre-calculated previous-neighbor indices.)
- Return type:
Updated spins array.
- models.xy_model.xy_step_random_numba(*, spins: ndarray, beta: float, J: float, idx_next: ndarray, idx_prev: ndarray) ndarray[source]
Perform one full Monte Carlo sweep of the XY lattice using random sequential updates.
- Parameters:
spins ((N, N, 2) array of unit vectors.)
beta (Inverse temperature 1/kT.)
J (Coupling constant.)
idx_next (Pre-calculated next-neighbor indices.)
idx_prev (Pre-calculated previous-neighbor indices.)
- Return type:
Updated spins array.
- models.xy_model.xy_wolff_step_numba(*, spins: ndarray, beta: float, J: float, idx_next: ndarray, idx_prev: ndarray) ndarray[source]
Perform one Wolff cluster flip on the XY lattice (Wolff-Evertz reflection).
A random mirror-plane axis r̂ is sampled uniformly from S^1. Each spin projects onto r̂ as sigma_i = s_i · r̂. Bonds between neighbouring sites with aligned projections (sigma_i * sigma_j > 0) are activated with probability
P_add = 1 - exp(-2 beta J sigma_i sigma_j). All cluster spins are then reflected through the plane perpendicular to r̂:s -> s - 2 (s · r̂) r̂, which preserves unit length exactly and satisfies detailed balance for the pure exchange term.One call constitutes one cluster sweep.
parallel=Trueis silently ignored.- Parameters:
spins ((N, N, 2) array of unit vectors.)
beta (Inverse temperature 1/kT.)
J (Coupling constant.)
idx_next (Pre-calculated next-neighbor indices (PBC).)
idx_prev (Pre-calculated previous-neighbor indices (PBC).)
- Return type:
Updated spins array.