Spirit - Spin Simulation Framework¶
SPIRIT¶
SPIN SIMULATION FRAMEWORK
Core Library:
Service | System | Compiler | Status |
---|---|---|---|
Travis- CI | Ubuntu 14.04 macOS | GCC 6 Clang | master: |image0 |deve lop: |image1 | |
AppVeyo r | Windows | MSVC14 MSVC14.1 | master: |image2 |deve lop: |image3 | |
`Python package <https://pypi.org/project/spirit/>`__:
Branch | Python Package Coverage | Core Library Coverage |
---|---|---|
master: | ||
develop: |
The code is released under MIT License. If you intend to present and/or publish scientific results or visualisations for which you used Spirit, please read the REFERENCE.md
This is an open project and contributions and collaborations are always welcome!! See Contributing on how to contribute or write an email to g.mueller@fz-juelich.de For contributions and affiliations, see CONTRIBUTORS.md
Please note that a version of the Spirit Web interface is hosted by the Research Centre Jülich at http://juspin.de

Contents¶
Introduction¶
A modern framework for magnetism science on clusters, desktops & laptops and even your Phone¶
Spirit is a platform-independent framework for spin dynamics, written in C++11. It combines the traditional cluster work, using using the command-line, with modern visualisation capabilites in order to maximize scientists’ productivity.
“It is unworthy of excellent men to lose hours like slaves in the labour of calculation which could safely be relegated to anyone else if machines were used.” - Gottfried Wilhelm Leibniz
Our goal is to build such machines. The core library of the Spirit framework provides an easy to use API, which can be used from almost any programming language, and includes ready-to-use python bindings. A powerful desktop user interface is available, providing real-time visualisation and control of parameters.
Physics Features¶
- Atomistic Spin Lattice Heisenberg Model including also DMI and dipole-dipole
- Spin Dynamics simulations obeying the Landau-Lifschitz-Gilbert equation
- Direct Energy minimisation with different solvers
- Minimum Energy Path calculations for transitions between different spin configurations, using the GNEB method
Highlights of the Framework¶
- Cross-platform: everything can be built and run on Linux, OSX and Windows
- Standalone core library with C API which can be used from almost any programming language
- Python package making complex simulation workflows easy
- Desktop UI with powerful, live 3D visualisations and direct control of most system parameters
- Modular backends including parallelisation on GPU (CUDA) and CPU (OpenMP)
Documentation¶
More details may be found at spirit-docs.readthedocs.io or in the Reference section including * Framework build instructions * Core build instructions * Core API Reference * Python API Reference * Input File Reference
There is also a Wiki, hosted by the Research Centre Jülich.
Getting started with the Desktop Interface¶
See BUILD.md on how to install the desktop user interface.

The user interface provides a powerful OpenGL visualisation window using the VFRendering library. It provides functionality to - Control Calculations - Locally insert Configurations (homogeneous, skyrmions, spin spiral, … ) - Generate homogeneous Transition Paths - Change parameters of the Hamiltonian - Change parameters of the Method and Solver - Configure the Visualization (arrows, isosurfaces, lighting, …)
See the UI-QT Reference for the key bindings of the various features.
Unfortunately, distribution of binaries for the Desktop UI is not possible due to the restrictive license on QT-Charts.
Getting started with the Python Package¶
To install the Spirit python package, either build and install from source or simply use
pip install spirit
With this package you have access to powerful Python APIs to run and control dynamics simulations or optimizations. This is especially useful for work on clusters, where you can now script your workflow, never having to re-compile when testing, debugging or adding features.
The most simple example of a spin dynamics simulation would be
from spirit import state, simulation
with state.State("input/input.cfg") as p_state:
simulation.PlayPause(p_state, "LLG", "SIB")
Where "SIB"
denotes the semi-implicit method B and the starting
configuration will be random.
To add some meaningful content, we can change the initial configuration by inserting a Skyrmion into a homogeneous background:
def skyrmion_on_homogeneous(p_state):
from spirit import configuration
configuration.PlusZ(p_state)
configuration.Skyrmion(p_state, 5.0, phase=-90.0)
If we want to calculate a minimum energy path for a transition, we need to generate a sensible initial guess for the path and use the GNEB method. Let us consider the collapse of a skyrmion to the homogeneous state:
from spirit import state, chain, configuration, transition, simulation
### Copy the system a few times
chain.Image_to_Clipboard(p_state)
for number in range(1,7):
chain.Insert_Image_After(p_state)
noi = chain.Get_NOI(p_state)
### First image is homogeneous with a Skyrmion in the center
configuration.PlusZ(p_state, idx_image=0)
configuration.Skyrmion(p_state, 5.0, phase=-90.0, idx_image=0)
simulation.PlayPause(p_state, "LLG", "VP", idx_image=0)
### Last image is homogeneous
configuration.PlusZ(p_state, idx_image=noi-1)
simulation.PlayPause(p_state, "LLG", "VP", idx_image=noi-1)
### Create transition of images between first and last
transition.Homogeneous(p_state, 0, noi-1)
### GNEB calculation
simulation.PlayPause(p_state, "GNEB", "VP")
where "VP"
denotes a direct minimization with the velocity
projection algorithm.
You may also use Spirit order to extract quantitative data, such as the energy.
def evaluate(p_state):
from spirit import system, quantities
M = quantities.Get_Magnetization(p_state)
E = system.Get_Energy(p_state)
return M, E
Obviously you may easily create significantly more complex workflows and use Python to e.g. pre- or post-process data or to distribute your work on a cluster and much more!
Contributing¶
Contributions are always welcome!
- Fork this repository
- Check out the develop branch:
git checkout develop
- Create your feature branch:
git checkout -b feature-something
- Commit your changes:
git commit -am 'Add some feature'
- Push to the branch:
git push origin feature-something
- Submit a pull request
Please keep your pull requests feature-specific and limit yourself to one feature per feature branch. Remember to pull updates from this repository before opening a new feature branch.
If you are unsure where to add you feature into the code, please do not hesitate to contact us.
There is no strict coding guideline, but please try to match your code style to the code you edited or to the style in the respective module.
We aim to adhere to the “git flow” branching model: http://nvie.com/posts/a-successful-git-branching-model/
Release versions (master
branch) are taggedmajor.minor.patch
, starting at1.0.0
Download the latest stable version from https://github.com/spirit-code/spirit/releases
The develop branch contains the latest updates, but is generally less consistently tested than the releases.
Spirit Desktop UI¶
The cross-platform QT desktop user interface provides a productive tool for Spin simulations, providing powerful real-time visualisations and access to simulation parameters, as well as other very useful features.
See the framework build instructions for information on how to build the user interface on your machine.
Physics features¶
Insert Configurations: - White noise - (Anti-) Skyrmions - Domains - Spin Spirals
You may manipulate the Hamiltonian as well as simulation parameters and your output file configuration:
You may start and stop simulation and directly interact with a running simulation. - LLG Simulation: Dynamics and Minimization - GNEB: create transitions and calculate minimum energy paths
By copying and inserting spin systems and manipulating them you may create arbitrary transitions between different spin states to use them in GNEB calculations. Furthermore you can choose different images to be climbing or falling images during your calculation.
Real-time visualisation¶
This feature is most powerful for 3D systems but shows great use for the analysis of dynamical processes and understanding what is happening in your system during a simulation instead of post-processing your data.
- Arrows, Surface (2D/3D), Isosurface
- Spins or Eff. Field
- Every n’th arrow
- Spin Sphere
- Directional & Position filters
- Colormaps
You can also create quite complicate visualisations by combining these different features in order to visualise complex states in 3D systems:

Note that a data plot is available to visualise your chain of spin systems. It can also show interpolated energies if you run a GNEB calculation.

Additional features¶
- Drag mode: drag, copy, insert, change radius
- Screenshot
- Read configuration or chain
- Save configuration or chain
Key bindings¶
Note that some of the keybindings may only work correctly on US keyboard layout.
UI Controls¶
Effect | Keystroke |
---|---|
Show this | F1 |
Toggle Settings | F2 |
Toggle Plots | F3 |
Toggle Debug | F4 |
Toggle “Dragging” mode | F5 |
Toggle large visualization | F10 / Ctrl+F |
Toggle full-screen window | F11 / Ctrl+Shift+F |
Screenshot of Visualization region | F12 / Home |
Toggle OpenGL Visualization | Ctrl+Shift+V |
Try to return focus to main UI | Esc |
Camera Controls¶
Effect | Keystroke |
---|---|
Rotate the camera around | Left mouse / W A S D ( Shift to go slow) |
Move the camera around | Left mouse / T F G H ( Shift to go slow) |
Zoom in on focus point | Scroll mouse ( Shift to go slow) |
Set the camera in X, Y or Z direction | X Y Z ( shift to invert) |
Control Simulations¶
Effect | Keystroke |
---|---|
Play/Pause | Space |
Cycle Method | Ctrl+M |
Cycle Solver | Ctrl+S |
Manipulate the current images¶
Effect | Keystroke |
---|---|
Random configuration | Ctrl+R |
Add tempered noise | Ctrl+N |
Insert last used configuration | Enter |
Visualisation¶
Effect | Keystroke |
---|---|
Use more/less data points of the vector field | +/- |
Regular Visualisation Mode | 1 |
Isosurface Visualisation Mode | 2 |
Slab (X,Y,Z) Visualisation Mode | 3 4 5 |
Cycle Visualisation Mode | / |
Move Slab | , / . ( Shift to go faster) |
Building Spirit’s Framework Components¶
The Spirit framework is designed to run across different platforms
and so the build process is set up with CMake
, which will generate
the appropriate build scripts for each platform.
Please be aware that our CMake scripts are written for our use cases and you may need to adapt some paths and options in the Root CMakeLists.txt.
General Build Process¶
The following assumes you are in the Spirit root directory.
Options¶
There are some important Options you may need to consider. You can find them under ### Build Flags ### in the Root `CMakeLists.txt <../CMakeLists.txt>`__. Otherwise, the developers’ defaults will be used.
Some Paths you can set under ### User Paths ### (just uncomment the corresponding line) are:
CMake Variable | Use |
---|---|
USER_COMPILER_CUSER_COMPILER_CXX | Name of the compiler you wish to use |
USER_PATH_COMPILER | Directory your compiler is located in |
USER_PATHS_IFF | use the default IFF (FZJ) cluster paths |
Clean¶
Clear the build directory using
./clean.sh
or
rm -rf build && mkdir build
Further helper scripts for clean-up are clean_log.sh
,
clean_output.sh
,
Generate Build Files¶
./cmake.sh
lets cmake generate makefiles for your system inside a
‘build’ folder. Simply call
./cmake.sh
or
cd build && cmake .. && cd ..
Passing -debug
to the script will cause it to create a debug
configuration, meaning that you will be able to properly debug the
entire application.
On Windows (no MSys) you can simply use the git bash to do this or use the CMake GUI. When using MSys etc., CMake will create corresponding MSys makefiles.
Building the Projects¶
To execute the build and linking of the executable, simply call
./make.sh -jN
or
cd build && make -jN && cd ..
where -jN
is optional, with N
the number of parallel build
processes you want to create.
On Windows (no MSys), CMake will by default have generated a Visual Studio Solution. Open the generated Solution in the Visual Studio IDE and build it there.
Running the Unit Tests¶
We use CMake
s CTest
for unit testing. You can run
ctest.sh
or
cd build && ctest --output-on-failure && cd ..
or execute any of the test executables manually. To execute the tests
from the Visual Studio IDE, simply rebuild the RUN_TESTS
project.
Installing Components¶
This is not yet supported! however, you can already run
./install.sh
or
cd build && make install && cd ..
Which on OSX should build a .app bundle.
Core Library¶
For detailed build instructions concerning the standalone core library or how to include it in your own project, see core/docs/BUILD.md. * Shared and static library * Python bindings * Julia bindings * Transpiling to JavaScript * Unit Tests
The Root `CMakeLists.txt <../CMakeLists.txt>`__ has a few options you can set:
CMake Options | Use |
---|---|
SPIRIT_USE_CUDA | Use CUDA to spee d up nume rica lly inte nsiv e part s of the core |
SPIRIT_USE_OPENMP | Use Open MP to spee d up nume rica lly inte nsiv e part s of the core |
SPIRIT_SCALAR_TYPE | Shou ld be e.g. ``do uble `` or ``fl oat` `. Sets the C++ type for scal ar vari able s, arra ys etc. |
SPIRIT_BUILD_TEST | Buil d unit test s for the core libr ary |
SPIRIT_BUILD_FOR_CXX | Buil d the stat ic libr ary for C++ appl icat ions |
SPIRIT_BUILD_FOR_JULI A | Buil d the shar ed libr ary for Juli a |
SPIRIT_BUILD_FOR_PYTH ON | Buil d the shar ed libr ary for Pyth on |
SPIRIT_BUILD_FOR_JS | Buil d the Java Scri pt libr ary (use s a diff eren t tool chai n!) |
Desktop User Interface¶
Dependencies | Versions |
---|---|
OpenGL Drivers | >= 3.3 |
CMake | >= 3.1 |
QT | >= 5.7 including QT-Charts |
Note that in order to build with QT as a dependency on Windows, you
may need to add path/to/qt/qtbase/bin
to your PATH variable.
Necessary OpenGL drivers should be available through the regular drivers for any remotely modern graphics card.
CMake Options | Use |
---|---|
SPIRIT_BUILD_FOR_C XX | Buil d the C++ inte rfac es (con sole or QT) inst ead of othe rs |
UI_CXX_USE_QT | Buil d qt user inte rfac e inst ead of cons ole vers ion |
USER_PATH_QT | The path to your CMak e inst alla tion |
BUNDLE_APP | On OSX, crea te .app bund le (not yet full y func tion al) |
SPIRIT INPUT FILES¶
The following sections will list and explain the input file keywords.
- General Settings and Log
- Geometry
- Hamiltonian
- Method Output
- Method Parameters
- Pinning
- Disorder and Defects
General Settings and Log¶
### Add a tag to output files (for timestamp use "<time>")
output_file_tag some_tag
### Save input parameters on creation of State
log_input_save_initial 0
### Save input parameters on deletion of State
log_input_save_final 0
### Print log messages to the console
log_to_console 1
### Print messages up to (including) log_console_level
log_console_level 5
### Save the log as a file
log_to_file 1
### Save messages up to (including) log_file_level
log_file_level 5
Except for SEVERE
and ERROR
, only log messages up to
log_console_level
will be printed and only messages up to
log_file_level
will be saved. If log_to_file
, however is set to
zero, no file is written at all.
Log Levels | Integer | Description |
---|---|---|
ALL | 0 | Everything |
SEVERE | 1 | Only severe errors |
ERROR | 2 | Also non-fatal errors |
WARNING | 3 | Also warnings |
PARAMETER | 4 | Also input parameters |
INFO | 5 | Also info-messages |
DEBUG | 6 | Also deeper debug-info |
Geometry¶
The Geometry of a spin system is specified in form of a bravais lattice and a basis cell of atoms. The number of basis cells along each principal direction of the basis can be specified. Note: the default basis is a single atom at (0,0,0).
3D simple cubic example:
### The bravais lattice type
bravais_lattice sc
### Number of basis cells along principal
### directions (a b c)
n_basis_cells 100 100 10
2D honeycomb example:
### The bravais lattice type
bravais_lattice hex2d
### n No of spins in the basis cell
### 1.x 1.y 1.z position of spins within basis
### 2.x 2.y 2.z cell in terms of bravais vectors
basis
2
0 0 0
0.86602540378443864676 0.5 0
### Number of basis cells along principal
### directions (a b c)
n_basis_cells 100 100 1
The bravais lattice can be one of the following:
Bravais Lattice Type | Keyword | Comment |
---|---|---|
Simple cubic | sc | |
Body-centered cubic | bcc | |
Face-centered cubic | fcc | |
Hexagonal (2D) | hex2d | 60deg angle |
Hexagonal (2D) | hex2d60 | 60deg angle |
Hexagonal (2D) | hex2d120 | 120deg angle |
Hexagonal closely packed | hcp | 120deg, not yet implemented |
Hexagonal densely packed | hdp | 60deg, not yet implemented |
Rhombohedral | rho | not yet implemented |
Simple-tetragonal | stet | not yet implemented |
Simple-orthorhombic | so | not yet implemented |
Simple-monoclinic | sm | not yet implemented |
Simple triclinic | stri | not yet implemented |
Alternatively it can be input manually, either through vectors or as the bravais matrix:
### bravais_vectors or bravais_matrix
### a.x a.y a.z a.x b.x c.x
### b.x b.y b.z a.y b.y c.y
### c.x c.y c.z a.z b.z c.z
bravais_vectors
1.0 0.0 0.0
0.0 1.0 0.0
0.0 0.0 1.0
A lattice constant can be used for scaling:
### Scaling constant
lattice_constant 1.0
Hamiltonian¶
Note that you select the Hamiltonian you use with the hamiltonian
keyword.
Isotropic Heisenberg Hamiltonian:
Interactions are handled in terms of neighbours. You may specify shell-wise interaction parameters:
### Hamiltonian Type (heisenberg_neighbours, heisenberg_pairs, gaussian)
hamiltonian heisenberg_neighbours
### boundary_conditions (in a b c) = 0(open), 1(periodical)
boundary_conditions 1 1 0
### external magnetic field vector[T]
external_field_magnitude 25.0
external_field_normal 0.0 0.0 1.0
### µSpin
mu_s 2.0
### Uniaxial anisotropy constant [meV]
anisotropy_magnitude 0.0
anisotropy_normal 0.0 0.0 1.0
### Exchange constants [meV] for the respective shells
### Jij should appear after the >Number_of_neighbour_shells<
n_neigh_shells_exchange 2
jij 10.0 1.0
### Chirality of DM vectors (+/-1=bloch, +/-2=neel)
dm_chirality 2
### DM constant [meV]
n_neigh_shells_dmi 1
dij 6.0
### Dipole-Dipole radius
dd_radius 0.0
If you have a nontrivial basis cell, note that you should specify
mu_s
for all atoms in your basis cell.
Anisotropy: By specifying a number of anisotropy axes via
n_anisotropy
, one or more anisotropy axes can be set for the atoms
in the basis cell. Specify columns via headers: an index i
and an
axis Kx Ky Kz
or Ka Kb Kc
, as well as optionally a magnitude
K
.
Pair-wise Heisenberg Hamiltonian:
Interactions are specified pair-wise. Single-threaded applications can thus calculate interactions twice as fast as for the neighbour-wise case. You may specify shell-wise interaction parameters.
### Hamiltonian Type (heisenberg_neighbours, heisenberg_pairs, gaussian)
hamiltonian heisenberg_pairs
### Boundary_conditions (in a b c) = 0(open), 1(periodical)
boundary_conditions 1 1 0
### External magnetic field vector[T]
external_field_magnitude 25.0
external_field_normal 0.0 0.0 1.0
### µSpin
mu_s 2.0
### Uniaxial anisotropy constant [meV]
anisotropy_magnitude 0.0
anisotropy_normal 0.0 0.0 1.0
### Dipole-Dipole radius
dd_radius 0.0
### Pairs
n_interaction_pairs 3
i j da db dc Jij Dij Dijx Dijy Dijz
0 0 1 0 0 10.0 6.0 1.0 0.0 0.0
0 0 0 1 0 10.0 6.0 0.0 1.0 0.0
0 0 0 0 1 10.0 6.0 0.0 0.0 1.0
### Quadruplets
n_interaction_quadruplets 1
i j da_j db_j dc_j k da_k db_k dc_k l da_l db_l dc_l Q
0 0 1 0 0 0 0 1 0 0 0 0 1 3.0
If you have a nontrivial basis cell, note that you should specify
mu_s
for all atoms in your basis cell.
Anisotropy: By specifying a number of anisotropy axes via
n_anisotropy
, one or more anisotropy axes can be set for the atoms
in the basis cell. Specify columns via headers: an index i
and an
axis Kx Ky Kz
or Ka Kb Kc
, as well as optionally a magnitude
K
.
Pairs: Leaving out either exchange or DMI in the pairs is allowed and
columns can be placed in arbitrary order. Note that instead of
specifying the DM-vector as Dijx Dijy Dijz
, you may specify it as
Dija Dijb Dijc
if you prefer. You may also specify the magnitude
separately as a column Dij
, but note that if you do, the vector
(e.g. Dijx Dijy Dijz
) will be normalized.
Quadruplets: Columns for these may also be placed in arbitrary order.
Separate files: The anisotropy, pairs and quadruplets can be placed
into separate files, you can use anisotropy_from_file
,
pairs_from_file
and quadruplets_from_file
.
If the headers for anisotropies, pairs or quadruplets are at the top of
the respective file, it is not necessary to specify n_anisotropy
,
n_interaction_pairs
or n_interaction_quadruplets
respectively.
### Pairs
interaction_pairs_file input/pairs.txt
### Quadruplets
interaction_quadruplets_file input/quadruplets.txt
Gaussian Hamiltonian:
This is a testing Hamiltonian consisting of the superposition of gaussian potentials. It does not contain interactions.
hamiltonian gaussian
### Number of Gaussians
n_gaussians 2
### Gaussians
### a is the amplitude, s is the width, c the center
### the directions c you enter will be normalized
### a1 s1 c1.x c1.y c1.z
### ...
gaussians
1 0.2 -1 0 0
0.5 0.4 0 0 -1
Method Output¶
For llg
and equivalently mc
and gneb
, you can specify which
output you want your simulations to create. They share a few common
output types, for example:
llg_output_any 1 # Write any output at all
llg_output_initial 1 # Save before the first iteration
llg_output_final 1 # Save after the last iteration
Note in the following that step
means after each N
iterations
and denotes a separate file for each step, whereas archive
denotes
that results are appended to an archive file at each step.
LLG:
llg_output_energy_step 0 # Save system energy at each step
llg_output_energy_archive 1 # Archive system energy at each step
llg_output_energy_spin_resolved 0 # Also save energies for each spin
llg_output_energy_divide_by_nspins 1 # Normalize energies with number of spins
llg_output_configuration_step 1 # Save spin configuration at each step
llg_output_configuration_archive 0 # Archive spin configuration at each step
MC:
mc_output_energy_step 0
mc_output_energy_archive 1
mc_output_energy_spin_resolved 0
mc_output_energy_divide_by_nspins 1
mc_output_configuration_step 1
mc_output_configuration_archive 0
GNEB:
gneb_output_energies_step 0 # Save energies of images in chain
gneb_output_energies_interpolated 1 # Also save interpolated energies
gneb_output_energies_divide_by_nspins 1 # Normalize energies with number of spins
gneb_output_chain_step 0 # Save the whole chain at each step
Method Parameters¶
Again, the different Methods share a few common parameters. On the example of the LLG Method:
### Maximum wall time for single simulation
### hh:mm:ss, where 0:0:0 is infinity
llg_max_walltime 0:0:0
### Force convergence parameter
llg_force_convergence 10e-9
### Number of iterations
llg_n_iterations 2000000
### Number of iterations after which to save
llg_n_iterations_log 2000
LLG:
### Seed for Random Number Generator
llg_seed 20006
### Damping [none]
llg_damping 0.3E+0
### Time step dt
llg_dt 1.0E-3
### Temperature [K]
llg_temperature 0
llg_temperature_gradient_direction 1 0 0
llg_temperature_gradient_inclination 0.0
### Spin transfer torque parameter proportional to injected current density
llg_stt_magnitude 0.0
### Spin current polarisation normal vector
llg_stt_polarisation_normal 1.0 0.0 0.0
MC:
### Seed for Random Number Generator
mc_seed 20006
### Temperature [K]
mc_temperature 0
### Acceptance ratio
mc_acceptance_ratio 0.5
GNEB:
### Constant for the spring force
gneb_spring_constant 1.0
### Number of energy interpolations between images
gneb_n_energy_interpolations 10
Pinning¶
Note that for this feature you need to build with
SPIRIT_ENABLE_PINNING
set to ON
in cmake.
For each lattice direction a
b
and c
, you have two choices
for pinning. For example to pin n
cells in the a
direction, you
can set both pin_na_left
and pin_na_right
to different values or
set pin_na
to set both to the same value. To set the direction of
the pinned cells, you need to give the pinning_cell
keyword and one
vector for each basis atom.
You can for example do the following to create a U-shaped pinning in x-direction:
# Pin left side of the sample (2 rows)
pin_na_left 2
# Pin top and bottom sides (2 rows each)
pin_nb 2
# Pin the atoms to x-direction
pinning_cell
1 0 0
To specify individual pinned sites (overriding the above pinning settings), insert a list into your input. For example:
### Specify the number of pinned sites and then the directions
### ispin S_x S_y S_z
n_pinned 3
0 1.0 0.0 0.0
1 0.0 1.0 0.0
2 0.0 0.0 1.0
You may also place it into a separate file with the keyword
pinned_from_file
, e.g.
### Read pinned sites from a separate file
pinned_from_file input/pinned.txt
The file should either contain only the pinned sites or you need to
specify n_pinned
inside the file.
Disorder and Defects¶
Note that for this feature you need to build with
SPIRIT_ENABLE_DEFECTS
set to ON
in cmake.
Disorder is not yet implemented.
To specify defects, be it vacancies or impurities, you may fix atom types for sites of the whole lattice by inserting a list into your input. For example:
### Atom types: type index 0..n or or vacancy (type < 0)
### Specify the number of defects and then the defects
### ispin itype
n_defects 3
0 -1
1 -1
2 -1
You may also place it into a separate file with the keyword
defects_from_file
, e.g.
### Read defects from a separate file
defects_from_file input/defects.txt
The file should either contain only the defects or you need to specify
n_defects
inside the file.
BUILD THE SPIRIT LIBRARY¶
C/C++ Library¶
Dependencies | Versions |
---|---|
gcc | >= 5.1 (C++11 stdlib) or any modern compiler |
CMake | >= 3.1 |
This is pretty much a standalone library and should be easy to implement
into existing projects in which CMake is already used. It can be built
as shared or static and API headers are located in
path/to/Spirit/core/include/Spirit
CMake Options | Use |
---|---|
SPIRIT_USE_CUDA | Use CUDA to spee d up nume rica lly inte nsiv e part s of the core |
SPIRIT_USE_OPENMP | Use Open MP to spee d up nume rica lly inte nsiv e part s of the core |
SPIRIT_SCALAR_TYPE | Shou ld be e.g. ``do uble `` or ``fl oat` `. Sets the C++ type for scal ar vari able s, arra ys etc. |
SPIRIT_BUILD_TEST | Buil d unit test s for the core libr ary |
SPIRIT_BUILD_FOR_CXX | Buil d the stat ic libr ary for C++ appl icat ions |
SPIRIT_BUILD_FOR_JULI A | Buil d the shar ed libr ary for Juli a |
SPIRIT_BUILD_FOR_PYTH ON | Buil d the shar ed libr ary for Pyth on |
SPIRIT_BUILD_FOR_JS | Buil d the Java Scri pt libr ary (use s a diff eren t tool chai n!) |
Note that the CMake Options for the core library are also set in the root CMakeLists.txt.
Integrating into other CMake Projects¶
This is currently untested, but you should be able to use Spirit
either via ExternalProject
or AddSubdirectory
. You should thus
be able to simply copy the core
directory as Spirit
into the
thirdparty folder of your Project and use it as is done commonly.
Note that setting qhull_LIBS
and qhull_INCLUDE_DIRS
if qhull is
present on your disk removes the need for Spirit to clone and build it
separately.
A set of CMake variables is pushed to the parent scope by the core CMakeLists. These include: - SPIRIT_LIBRARIES - SPIRIT_LIBRARIES_STATIC - SPIRIT_INCLUDE_DIRS
Unit Tests¶
Dependencies | Versions |
---|---|
gcc | >= 4.8.1 (C++11 stdlib) or any modern compiler |
CMake | >= 2.8.12 |
Python | >= 2.7 or >= 3.5 (for the python tests) |
The core library includes a set of unit tests in order to make sure
certain conditions are fulfilled by the librarys functions. We use
CMake
s CTest
for unit testing. You can run
cd build && ctest --output-on-failure && cd ..
or execute any of the test executables manually. To execute the tests
from the Visual Studio IDE, simply rebuild the RUN_TESTS
project.
CMake Options | Use |
---|---|
CORE_BUILD_TEST | Build unit tests for the core library |
BUILD_UI_PYTHON | Python tests are only built if this is ON |
Python Package¶
Dependencies | Versions |
---|---|
Python | >= 2.7 or >= 3.5 |
CMake Options | Use |
---|---|
BUILD_UI_PYTHON | Buil d the shar ed libr ary for pyth on |
CORE_BUILD_TEST | Will buil d pyth on unit test s if the shar ed lib for pyth on is buil t |
You may use the python package simply by adding a
path/to/Spirit/core/python
to your PYTHONPATH.
To build a proper python package including binary wheel, you need to build the shared library for Python and then execute
cd path/to/Spirit/core/python
python setup.py sdist bdist_wheel
The resulting package can then be installed with e.g.
pip install -e /path/to/package/dist/spirit-<tags>.whl
Julia Package¶
These bindings are currently completely experimental. However, it seems it would be easy to create Julia bindings analogous to the Python Bindings.
Dependencies | Versions |
---|---|
None |
CMake Options | Use |
---|---|
BUILD_UI_JULIA | Build the shared library for julia |
SPIRIT API¶
This will list the available API functions of the Spirit library.
The API is exposed as a C interface and may thus also be used from other
languages, such as Python, Julia or even JavaScript (see ui-web). The
API revolves around a simulation State
which contains all the
necessary data and keeps track of running Solvers.
The API exposes functions for: * Control of simulations * Manipulation of parameters * Extracting information * Generating spin configurations and transitions * Logging messages * Reading spin configurations * Saving datafiles
State Managment¶
To create a new state with one chain containing a single image, initialized by an input file, and run the most simple example of a spin dynamics simulation:
#import "Spirit/State.h"
#import "Spirit/Simulation.h"
const char * cfgfile = "input/input.cfg"; // Input file
State * p_state = State_Setup(cfgfile); // State setup
Simulation_PlayPause(p_state, "LLG", "SIB"); // Start a LLG simulation using the SIB solver
State_Delete(p_state) // State cleanup
A new state can be created with State_Setup()
, where you can pass a
config file specifying your initial system parameters. If
you do not pass a config file, the implemented defaults are used. Note
that you currently cannot change the geometry of the systems in your
state once they are initialized.
The State struct is passed around in an application to make the simulation’s state available.
State manipulation function | Retur n | Effe ct |
---|---|---|
State_Setup( const char * config_file ) |
``Sta te *` ` | Crea te new stat e by pass ing a conf ig file |
State_Update( State * ) |
voi
d |
Upda te the stat e to hold curr ent valu es |
State_Delete( State * ) |
voi
d |
Dele te a stat e |
State_To_Config( State *, const char * config_file, const ch
ar * original_config_file) |
voi
d |
Writ e a conf ig file whic h will resu lt in the same stat e if used in ``St ate_ Setu p()` ` |
State_DateTime( State * ) |
``con st ch ar *` ` | Get date time tag of the crea tion of the stat e |
System¶
System Information | Return | Effect |
---|---|---|
System_Get_Index( State *) |
int |
Returns System’s Index |
System_Get_NOS( State *, int idx_image
, int idx_chain ) |
int |
Return System’s number of spins |
System Data | Return | Effe ct |
---|---|---|
System_Get_Spin_Directions( State *, int idx_image, int id
x_chain ) |
scal
ar * |
Get Syst em’s spin dire ctio n |
System_Get_Spin_Effective_Field( State *, int idx_image, i
nt idx_chain ) |
scal
ar * |
Get Syst em’s spin effe ctiv e fiel d |
System_Get_Rx( State *, int idx_image, int idx_chain) |
floa
t |
Get a Syst em’s reac tion coor dina te in it’s chai n |
``System_Get_Energy( State *, int idx_image, int idx_chain ) `` | floa
t |
Get Syst em’s ener gy |
System_Get_Energy_Array( State * energies, float * energie
s, int idx_image, int idx_chain ) |
``void `` | Ener gy Arra y (Sho uld NOT be used ) |
System Output | Effect |
---|---|
System_Print_Energy_Array( State *, int idx_image, int idx_chai
n) |
Print on the consol e State’ s energy array |
System Update | Effect |
---|---|
System_Update_Data( State *, int idx_image, int idx_chain) |
Update State’ s data. Used mainly for plotti ng |
Simulation¶
With Simulation_*
functions one can control and get information from
the State’s Solvers and their Methods.
Simulation Basics | Effect |
---|---|
Simulation_SingleShot( State *, const char * c_method_type, con
st char * c_solver_type, int n_iterations, int n_iterations_log,
int idx_image, int idx_chain ) |
Execut es a single Optimi zation iterat ion with a given method (CAUTI ON! does not check for alread y runnin g simula tions) |
Simulation_PlayPause( State *, const char * c_method_type, cons
t char * c_solver_type, int n_iterations, int n_iterations_log, i
nt idx_image, int idx_chain ) |
Play/P ause functi onalit y |
Simulation_Stop_All( State * ) |
Stop all State’ s simula tions |
Simulation Data | Return | Effe ct |
---|---|---|
Simulation_Get_MaxTorqueComponent( State *, int idx_i
mage, int idx_chain ) |
float |
Get Simu lati on’s maxi mum torq ue comp onen t |
Simulation_Get_IterationsPerSecond( State *, int idx_
image, int idx_chain ) |
float |
Get Simu lati on’s iter atio ns per seco nd |
Simulation_Get_Solver_Name( State *, int idx_image, i
nt idx_chain ) |
const cha
r * |
Get Solv er’s name |
Simulation_Get_Method_Name( State *, int idx_image, i
nt idx_chain ) |
const cha
r * |
Get Meth od’s name |
Simulation Running Checking | Return |
---|---|
Simulation_Running_Any_Anywhere( State * ) |
bool |
Simulation_Running_LLG_Anywhere( State * ) |
bool |
Simulation_Running_GNEB_Anywhere( State * ) |
bool |
Simulation_Running_LLG_Chain( State *state, int idx_chain ) |
bool |
Simulation_Running_Any( State *, int idx_image, int idx_chain ) |
bool |
Simulation_Running_LLG( State *, int idx_image, int idx_chain) |
bool |
Simulation_Running_GNEB( State *, int idx_chain ) |
bool |
Simulation_Running_MMF( State * ) |
bool |
Geometry¶
Geometry Data | Return |
---|---|
Geometry_Get_Spin_Positions( State *, int idx_image, int idx_c
hain ) |
scala
r * |
Geometry_Get_Bounds( State *, float min[3], float max[3], int
idx_image, int idx_chain ) |
``void` ` |
Geometry_Get_Center( State *, float center[3], int idx_image,
int idx_chain ) |
``void` ` |
``void` ` | |
Geometry_Get_Basis_Vectors( State *, float a[3], float b[3], f
loat c[3], int idx_image, int idx_chain ) |
int |
Geometry_Get_N_Basis_Atoms( State *, int idx_image, int idx_ch
ain ) |
``void` ` |
Geometry_Get_N_Cells( State *, int n_cells[3], int idx_image,
int idx_chain ) |
``void` ` |
Geometry_Get_Translation_Vectors( State *, float ta[3], float
tb[3], float tc[3], int idx_image, int idx_chain ) |
``void` ` |
Geometry_Get_Dimensionality( State *, int idx_image, int idx_c
hain ) |
int |
Geometry_Get_Triangulation( State *, const int **indices_ptr,
int n_cell_step, int idx_image, int idx_chain ) |
int |
Transitions¶
Transitions | Return | Effect |
---|---|---|
void |
||
Transition_Add_Noise_Temperature( State *, float tempe
rature, int idx_1, int idx_2, int idx_chain) |
void |
Quantities¶
Quantities | Return | Effect |
---|---|---|
Quantity_Get_Magnetization( State *, float m[3], int i
dx_image, int idx_chain ) |
void |
|
Quantity_Get_Topological_Charge( State *, int idx_imag
e, int idx_chain ) |
``float` ` |
Parameters¶
LLG Parameters Set | Retu rn | Effect |
---|---|---|
Parameters_Set_LLG_Time_Step( State *, float dt, int idx_i
mage, int idx_chain ) |
vo
id |
|
Parameters_Set_LLG_Damping( State *, float damping, int id
x_image, int idx_chain ) |
vo
id |
|
Parameters_Set_LLG_N_Iterations( State *, int n_iterations
, int idx_image, int idx_chain ) |
vo
id |
GNEB Parameters Set | Retu rn | Effect |
---|---|---|
Parameters_Set_GNEB_Spring_Constant( State *, float spring
_constant, int idx_image, int idx_chain ) |
vo
id |
|
Parameters_Set_GNEB_Climbing_Falling( State *, int image_t
ype, int idx_image, int idx_chain ) |
vo
id |
|
Parameters_Set_GNEB_N_Iterations( State *, int n_iteration
s, int idx_chain ) |
vo
id |
LLG Parameters Get | Retu rn | Effect |
---|---|---|
Parameters_Get_LLG_Time_Step( State *, float * dt, int idx
_image, int idx_chain) |
vo
id |
|
Parameters_Get_LLG_Damping( State *, float * damping, int
idx_image, int idx_chain) |
vo
id |
|
Parameters_Get_LLG_N_Iterations( State *, int idx_image, i
nt idx_chain) |
vo
id |
GNEB Parameters Get | Retu rn | Effect |
---|---|---|
Parameters_Get_GNEB_Spring_Constant( State *, float * spri
ng_constant, int idx_image, int idx_chain ) |
vo
id |
|
Parameters_Get_GNEB_Climbing_Falling( State *, int * image
_type, int idx_image, int idx_chain ) |
vo
id |
|
``Parameters_Get_GNEB_N_Iterations( State *, int idx_chain ) `` | in
t |
|
Parameters_Get_GNEB_N_Energy_Interpolations( State *, int
idx_chain ) |
in
t |
Chain¶
Get Chain’s information
Chain Info | Return | Effect |
---|---|---|
Chain_Get_Index( State * ) |
int |
|
Chain_Get_NOI( State *, int idx_chain ) |
int |
Move between images (change active_image
)
Chain moving | Return |
---|---|
Chain_prev_Image( State *, int idx_chain ) |
bool |
Chain_next_Image( State *, int idx_chain ) |
bool |
Chain_Jump_To_Image( State *, int idx_image, int idx_chain ) |
bool |
Insertion/deletion and replacement of images are done by
Chain control | Return |
---|---|
Chain_Image_to_Clipboard( State *, int idx_image, int idx_chain ) |
void |
Chain_Replace_Image( State *, int idx_image, int idx_chain ) |
void |
Chain_Insert_Image_Before( State *, int idx_image, int idx_chain ) |
void |
Chain_Insert_Image_After( State *, int idx_image, int idx_chain ) |
void |
Chain_Push_Back( State *, int idx_chain ) |
void |
Chain_Delete_Image( State *, int idx_image, int idx_chain ) |
bool |
Chain_Pop_Back( State *, int idx_chain ) |
bool |
Chain data | Return |
---|---|
Chain_Get_Rx( State *, float* Rx, int idx_chain ) |
``void `` |
Chain_Get_Rx_Interpolated( State *, float * Rx_interpolated, in
t idx_chain ) |
``void `` |
Chain_Get_Energy( State *, float* energy, int idx_chain ) |
``void `` |
Chain_Get_Energy_Interpolated( State *, float* E_interpolated,
int idx_chain ) |
``void `` |
Chain data update | Return |
---|---|
Chain_Update_Data( State *, int idx_chain ) |
void |
Chain_Setup_Data( State *, int idx_chain ) |
void |
Hamiltonian¶
Set Hamiltonian’s parameters
Hamiltonian Set | Retur n |
---|---|
Hamiltonian_Set_Boundary_Conditions( State *, const bool* period
ical, int idx_image, int idx_chain ) |
voi
d |
Hamiltonian_Set_mu_s( State *, float mu_s, int idx_image, int id
x_chain ) |
voi
d |
Hamiltonian_Set_Field( State *, float magnitude, const float* no
rmal, int idx_image, int idx_chain ) |
voi
d |
Hamiltonian_Set_Exchange( State *, int n_shells, const float* ji
j, int idx_image, int idx_chain ) |
voi
d |
Hamiltonian_Set_DMI( State *, float dij, int idx_image, int idx_
chain ) |
voi
d |
Hamiltonian_Set_BQE( State *, float dij, int idx_image, int idx_
chain ) |
voi
d |
Hamiltonian_Set_FSC( State *, float dij, int idx_image, int idx_
chain ) |
voi
d |
Hamiltonian_Set_Anisotropy( State *, float magnitude, const floa
t* normal, int idx_image, int idx_chain ) |
voi
d |
Hamiltonian_Set_STT( State *, float magnitude, const float * nor
mal, int idx_image, int idx_chain ) |
voi
d |
Hamiltonian_Set_Temperature( State *, float T, int idx_image, in
t idx_chain ) |
voi
d |
Get Hamiltonian’s parameters
Hamiltonian Get | Return |
---|---|
Hamiltonian_Get_Name( State *, int idx_image, int idx_chain
) |
const c
har * |
Hamiltonian_Get_Boundary_Conditions( State *, bool* periodic
al, int idx_image, int idx_chain ) |
void |
Hamiltonian_Get_mu_s( State *, float* mu_s, int idx_image, i
nt idx_chain ) |
void |
Hamiltonian_Get_Field( State *, float* magnitude, float* nor
mal, int idx_image, int idx_chain ) |
void |
Hamiltonian_Get_Exchange( State *, int* n_shells, float* jij
, int idx_image, int idx_chain ) |
void |
Hamiltonian_Get_Anisotropy( State *, float* magnitude, float
* normal, int idx_image, int idx_chain ) |
void |
|
void |
|
void |
|
void |
Hamiltonian_Get_STT( State *, float* magnitude, float* norma
l, int idx_image, int idx_chain ) |
void |
Hamiltonian_Get_Temperature( State *, float* T, int idx_imag
e, int idx_chain ) |
void |
Constants¶
Constants | Return | Description |
---|---|---|
Constants_mu_B() |
scalar |
The Bohr Magneton [meV / T] |
Constants_k_B() |
scalar |
The Boltzmann constant [meV / K] |
Log¶
Log Utilities | Retur n | Effect |
---|---|---|
Log_Send( State *, int level, int sender, const char * m
essage, int idx_image, int idx_chain ) |
voi
d |
Send a Log message |
Log_Get_Entries( State * ) |
std
::vec
tor<U
tilit
y::Lo
gEntr
y> |
Get the entries from the Log and write new number of entries into given int |
Log_Append( State * ) |
voi
d |
Append the Log to it’s file |
Log_Dump( State * ) |
voi
d |
Dump the Log into it’s file |
Log_Get_N_Entries( State * ) |
``int `` | Get the number of Log entries |
Log_Get_N_Errors( State * ) |
``int `` | Get the number of errors in the Log |
Log_Get_N_Warnings( State * ) |
``int `` | Get the number of warning s in the Log |
Log macro variables for Levels
Log Levels | value |
---|---|
Log_Level_All |
0 |
Log_Level_Severe |
1 |
Log_Level_Error |
2 |
Log_Level_Warning |
3 |
Log_Level_Parameter |
4 |
Log_Level_Info |
5 |
Log_Level_Debug |
6 |
Log macro variables for Senders
Log Senders | value |
---|---|
Log_Sender_All |
0 |
Log_Sender_IO |
1 |
Log_Sender_GNEB |
2 |
Log_Sender_LLG |
3 |
Log_Sender_MC |
4 |
Log_Sender_MMF |
5 |
Log_Sender_API |
6 |
Log_Sender_UI |
7 |
IO¶
Macros of File Formats for Vector Fields | values | Description |
---|---|---|
IO_Fileformat_Regular |
0 | sx sy sz (separated by whitespace) |
IO_Fileformat_Regular_Pos |
1 | px py pz sx sy sz (separated by whitespace) |
IO_Fileformat_CSV |
2 | sx, sy, sz (separated by commas) |
IO_Fileformat_CSV_Pos |
3 | px, py, pz, sx, sy, (sz separated by commas) |
IO_Fileformat_OVF_bin8 |
4 | OOMMF vector field (OVF) v2.0 file format |
IO_Fileformat_OVF_text |
6 |
Read and Write functions
Images IO | Retur n |
---|---|
``int `` | |
IO_Image_Write( State *state, const char *file, int format, cons
t char *comment, int idx_image, int idx_chain ) |
voi
d |
IO_Image_Append( State *state, const char *file, int format, con
st char *comment, int idx_image, int idx_chain ) |
voi
d |
IO_Image_Read( State *state, const char *file, int format, int i
dx_image_infile, int idx_image_inchain, int idx_chain ) |
voi
d |
Chains IO | Retur n |
---|---|
IO_Chain_Write( State *state, const char *file, int format, con
st char *comment, int idx_chain ) |
voi
d |
IO_Chain_Append( State *state, const char *file, int format, con
st char *comment, int idx_chain ) |
voi
d |
IO_Chain_Read( State *state, const char *file, int format, int s
tarting_image, int ending_image, int insert_idx, int idx_chain ) |
voi
d |
Energies from System
and Chain
System Energies | Return |
---|---|
|
``void` ` |
IO_Image_Write_Energy( State *, const char* file, int idx_imag
e, int idx_chain ) |
``void` ` |
Chain Energies | Return |
---|---|
IO_Chain_Write_Energies( State *, const char* file, int idx_ch
ain ) |
``void` ` |
IO_Chain_Write_Energies_Interpolated( State *, const char* fil
e, int idx_chain ) |
``void` ` |
SPIRIT Python API¶
State¶
To create a new state with one chain containing a single image, initialized by an input file, and run the most simple example of a spin dynamics simulation:
from spirit import state
from spirit import simulation
cfgfile = "input/input.cfg" # Input File
with state.State(cfgfile) as p_state: # State setup
simulation.PlayPause(p_state, "LLG", "SIB") # Start a LLG simulation using the SIB solver
or call setup and delete manually:
from spirit import state
from spirit import simulation
cfgfile = "input/input.cfg" # Input File
p_state = state.setup(cfgfile) # State setup
simulation.PlayPause(p_state, "LLG", "SIB") # Start a LLG simulation using the SIB solver
state.delete(p_state) # State cleanup
You can pass a config file specifying your initial system parameters. If you do not pass a config file, the implemented defaults are used. Note that you currently cannot change the geometry of the systems in your state once they are initialized.
State manipulation | Returns |
---|---|
setup( configfile="", quiet=False ) |
None |
delete(p_state ) |
None |
Chain¶
For having more images one can copy the active image in the Clipboard and then insert in a specified position of the chain.
chain.Image_to_Clipboard(p_state ) # Copy p_state to Clipboard
chain.Insert_Image_After(p_state ) # Insert the image from Clipboard right after the currently active image
For getting the total number of images in the chain
number_of_images = chain.Get_NOI(p_state )
Get Info | Returns | Description |
---|---|---|
Get_Index(p_state ) |
int |
Get Chain index |
Get_NOI(p_state, idx_chain=-1) |
int |
Get Chain number of images |
Get_Rx(p_state, idx_chain=-1) |
``Array `` | Get Rx |
Get_Rx_Interpolated(p_state, idx_c
hain=-1) |
``Array (float) `` | Get Rx interpolated |
``Get_Energy(p_state, idx_chain=-1)` ` | ``Array (float) `` | Get Energy of every System in Chain |
Get_Energy_Interpolated(p_state, i
dx_chain=-1) |
``Array (float) `` | Get interpolated Energy of every System in Chain |
Image Manipulation | Returns | Description |
---|---|---|
Next_Image(p_state, idx_chain=-
1) |
``None` ` | Switch active to next image of chain (one with largest index). If the current active is the last there is no effect. |
Prev_Image(p_state, idx_chain=-
1) |
``None` ` | Switch active to previous image of chain (one with smaller index). If the current active is the first one there is no effect |
Jump_To_Image(p_state, idx_imag
e=-1, idx_chain=-1) |
``None` ` | Switch active to specific image of chain. If this image does not exist there is no effect. |
Image_to_Clipboard(p_state, idx
_image=-1, idx_chain=-1) |
``None` ` | Copy active image to clipboard |
Replace_Image(p_state, idx_imag
e=-1, idx_chain=-1) |
``None` ` | Replace active image in chain. If the image does not exist there is no effect. |
Insert_Image_Before(p_state, id
x_image=-1, idx_chain=-1) |
``None` ` | Inserts clipboard image before the current active image. Active image index is increment by one. |
Insert_Image_After(p_state, idx
_image=-1, idx_chain=-1) |
``None` ` | Insert clipboard image after the current active image. Active image has the same index. |
Push_Back(p_state, idx_chain=-1
) |
``None` ` | Insert clipboard image at end of chain (after the image with the largest index). |
Delete_Image(p_state, idx_image
=-1, idx_chain=-1) |
``None` ` | Delete active image. If index is specified delete the corresponding image. If the image does not exist there is no effect. |
``Pop_Back(p_state, idx_chain=-1) `` | ``None` ` | Delete image at end of chain. |
Data | Returns | Description |
---|---|---|
Update_Data(p_state, idx
_chain=-1) |
None |
Update the chain’s data (interpolated energies etc.) |
Setup_Data(p_state, idx_
chain=-1) |
None |
Setup the chain’s data arrays |
System¶
System | Ret urn s | Description |
---|---|---|
Get_Index(p_state) |
``i nt` ` | Returns the index of the currently active image |
Get_NOS(p_state, idx_image=-
1, idx_chain=-1) |
``i nt` ` | Returns the number of spins |
|
``[ 3*N OS] `` | Returns an numpy.Array of size
3*NOS with the components of each
spin’s vector |
Get_Energy(p_state, idx_imag
e=-1, idx_chain=-1) |
f
loa
t |
Returns the energy of the system |
Update_Data(p_state, idx_ima
ge=-1, idx_chain=-1) |
``N one `` | Update the data of the state |
Print_Energy_Array(p_state,
idx_image=-1, idx_chain=-1) |
``N one `` | Print the energy array of the state |
Constants¶
Physical Constants | Returns | Description |
---|---|---|
mu_B() |
float |
The Bohr Magneton [meV / T] |
k_B() |
float |
The Boltzmann constant [meV / K] |
hbar() |
float |
Planck’s constant over 2pi [meV*ps / rad] |
mRy() |
float |
Millirydberg [mRy / meV] |
gamma() |
float |
The Gyromagnetic ratio of electron [rad / (ps*T)] |
g_e() |
float |
The Electron g-factor [unitless] |
Geometry¶
Get Geometry parameters | Returns | Description |
---|---|---|
Get_Bounds(p_state, idx_image=-1
, idx_chain=-1) |
``[3], [3]` ` | Get bounds (minimum and maximum arrays) |
Get_Center(p_state, idx_image=-1
, idx_chain=-1) |
``float, fl oat, float` ` | Get center |
Get_Basis_Vectors(p_state, idx_i
mage=-1, idx_chain=-1) |
[3],[3],[
3] |
Get basis vectors |
Get_N_Cells(p_state, idx_image=-
1, idx_chain=-1) |
|
Get number of cells in each dimension |
|
[3],[3],[
3] |
Get translation vectors |
Get_Dimensionality(p_state, idx_
image=-1, idx_chain=-1) |
int |
Get dimensionality of the system |
Get_Spin_Positions(p_state, idx_
image=-1, idx_chain=-1) |
[3*NOS] |
Get Spin positions |
Get_Atom_Types(p_state, idx_imag
e=-1, idx_chain=-1) |
[NOS] |
Get atom types |
Hamiltonian¶
Set Parameters | Returns | Description |
---|---|---|
Set_Field(p_state, magnitude, direc
tion, idx_image=-1, idx_chain=-1) |
None |
Set external magnetic field |
Set_Anisotropy(p_state, magnitude,
direction, idx_image=-1, idx_chain=-1
) |
None |
Set anisotropy |
Log¶
Log manipulation | Return s | Description |
---|---|---|
Send(p_state, level, sender, message, idx_ima
ge=-1, idx_chain=-1) |
``None `` | Send a Log message |
Append(p_state) |
``None `` | Append Log to file |
Parameters¶
LLG¶
Set LLG Parameters | Returns |
---|---|
setIterations(p_state, n_iterations, n_iterations_log, idx_i
mage=-1, idx_chain=-1) |
None |
setDirectMinimization(p_state, use_minimization, idx_image=-
1, idx_chain=-1) |
None |
setConvergence(p_state, convergence, idx_image=-1, idx_chain
=-1) |
None |
setTimeStep(p_state, dt, idx_image=-1, idx_chain=-1) |
None |
setDamping(p_state, damping, idx_image=-1, idx_chain=-1) |
None |
setSTT(p_state, use_gradient, magnitude, direction, idx_imag
e=-1, idx_chain=-1) |
None |
setTemperature(p_state, temperature, idx_image=-1, idx_chain
=-1) |
None |
Get LLG Parameters | Returns |
---|---|
getIterations(p_state, idx_image=-1, idx_chain=-1) |
int, int |
getDirect_Minimization(p_state, idx_image=-1, idx_chain=-1) |
int |
getConvergence(p_state, idx_image=-1, idx_chain=-1) |
float |
getTimeStep(p_state, idx_image=-1, idx_chain=-1) |
float |
getDamping(p_state, idx_image=-1, idx_chain=-1) |
float |
getSTT(p_state, idx_image=-1, idx_chain=-1) |
float, [3], bool |
getTemperature(p_state, idx_image=-1, idx_chain=-1) |
float |
GNEB¶
Set GNEB Parameters | Returns |
---|---|
setIterations(p_state, n_iterations, n_iterations_log, idx_im
age=-1, idx_chain=-1) |
None |
setConvergence(p_state, convergence, idx_image=-1, idx_chain=
-1) |
None |
setSpringConstant(p_state, c_spring, idx_image=-1, idx_chain=
-1) |
None |
setClimbingFalling(p_state, image_type, idx_image=-1, idx_cha
in=-1) |
None |
setImageTypeAutomatically(p_state, idx_chain=-1) |
None |
Get GNEB Parameters | Returns |
---|---|
getIterations(p_state, idx_chain=-1) |
int, int |
getConvergence(p_state, idx_image=-1, idx_chain=-1) |
float |
getSpringConstant(p_state, idx_image=-1, idx_chain=-1) |
float |
getClimbingFalling(p_state, idx_image=-1, idx_chain=-1) |
int |
getEnergyInterpolations(p_state, idx_chain=-1) |
int |
Quantities¶
Get Physical Quantities | Returns |
---|---|
Get_Magnetization(p_state, idx_image=-1, idx_chain=-1) |
[3*float] |
Simulation¶
The available method_type
s are:
Method | Argument |
---|---|
Landau-Lifshitz-Gilbert | "LLG" |
Geodesic Nudged Elastic Band | "GNEB" |
Monte-Carlo | "MC" |
The available solver_type
s are:
Solver | Argument |
---|---|
Semi-Implicit Method B | "SIB" |
Heun Method | "Heun" |
Depondt Method | "Depondt" |
Velocity Projection | "VP" |
Nonlinear Conjugate Gradient | "NCG" |
Note that the VP and NCG Solvers are only meant for direct minimization and not for dynamics.
Simulation state | Retur ns |
---|---|
SingleShot(p_state, method_type, solver_type, n_iterations=-1, n
_iterations_log=-1, idx_image=-1, idx_chain=-1) |
Non
e |
PlayPause(p_state, method_type, solver_type, n_iterations=-1, n_
iterations_log=-1, idx_image=-1, idx_chain=-1) |
Non
e |
Stop_All(p_state) |
Non
e |
Running_Image(p_state, idx_image=-1, idx_chain=-1) |
``Boo lean` ` |
Running_Chain(p_state, idx_chain=-1) |
``Boo lean` ` |
Running_Collection(p_state) |
``Boo lean` ` |
Running_Anywhere_Chain(p_state, idx_chain=-1) |
``Boo lean` ` |
Running_Anywhere_Collection(p_state) |
``Boo lean` ` |
Transition¶
Transition options | Ret urn s | Description |
---|---|---|
Homogeneous(p_state, idx_1, id
x_2, idx_chain=-1) |
``N one `` | Generate homogeneous transition between two images of a chain |
chain=-1)`` |
``N one `` | Add some temperature-scaled noise to a transition between two images of a chain |
Input/Output¶
Macros of File Formats for Vector Fields | values | Description |
---|---|---|
IO_Fileformat_Regular |
0 | sx sy sz (separated by whitespace) |
IO_Fileformat_Regular_Pos |
1 | px py pz sx sy sz (separated by whitespace) |
IO_Fileformat_CSV |
2 | sx, sy, sz (separated by commas) |
IO_Fileformat_CSV_Pos |
3 | px, py, pz, sx, sy, (sz separated by commas) |
IO_Fileformat_OVF_bin8 |
4 | OOMMF vector field (OVF) v2.0 file format |
IO_Fileformat_OVF_text |
6 |
For Image | Description |
---|---|
Image_Read(p_state, filename, fileformat=0, idx_i
mage=-1, idx_chain=-1) |
Read an image from disk |
Image_Write(p_state, filename, fileformat=0, comm
ent=" ", idx_image=-1, idx_chain=-1) |
Write an image to disk |
Image_Append(p_state, filename, fileformat=0, com
ment=" ", idx_image=-1, idx_chain=-1) |
Append an image to an existing file |
For Chain | Description |
---|---|
Chain_Read(p_state, filename, fileformat=0, id
x_chain=-1) |
Read a chain of images from disk |
Chain_Write(p_state, filename, fileformat=0, c
omment=" ", idx_chain=-1) |
Write a chain of images to disk |
Chain_Append(p_state, filename, fileformat=0,
comment=" ", idx_chain=-1) |
Append a chain of images to disk |
Contributors¶
Gideon Müller¶
- RWTH Aachen
- University of Iceland
- PGI-1/IAS-1 at Forschungszentrum Jülich
General code design and project setup (including CMake). Implementation of the core library and UIs, most notably: - GNEB and MMF methods - Velocity projection solver - CUDA and OpenMP parallelizations of backend - C API and Pyhton bindings - C++ QT GUI and initial OpenGL code
(Oct. 2014 - ongoing)
Daniel Schürhoff¶
- RWTH Aachen
- PGI-1/IAS-1 at Forschungszentrum Jülich
Implementation of the initial core library, notably translating from Fortran90 to C++ and addition of STT to the SIB solver. Work on QT GUI and Python bindings.
(Oct. 2015 - Sept. 2016)
Nikolai Kiselev¶
- PGI-1/IAS-1 at Forschungszentrum Jülich
Scientific advice, general help and feedback, initial (Fortran90) implementations of: - isotropic Heisenberg Hamiltonian - Neighbour calculations - SIB solver - Monte Carlo methods
(2007 - ongoing)
Florian Rhiem¶
- Scientific IT-Systems, PGI/JCNS at Forschungszentrum Jülich
Implementation of C++ OpenGL code (VFRendering library), as well as JavaScript Web UI and WebGL code. Code design improvements, including the C API and CMake.
(Jan. 2016 - ongoing)
Stefanos Mavros¶
- RWTH Aachen
Work on unit testing and documentation and implementation of the Depondt solver.
(April 2017 - ongoing)
Constantin Disselkamp¶
- RWTH Aachen
Implementation and testing of gradient approximation of spin transfer torque.
(April 2017 - Juli 2017)
Pavel Bessarab¶
- …
Help with the initial GNEB implementation. Initial implementation of the HTST method.
(April 2015 - ongoing)
Ingo Heimbach¶
- Scientific IT-Systems, PGI/JCNS at Forschungszentrum Jülich
Implementation of the initial OpenGL code. Code design suggestions and other general help.
(Jan. 2016 - ongoing)
Mathias Redies, Maximilian Merte, Rene Suckert¶
- RWTH Aachen
Initial CUDA implementation and tests. Code optimizations, suggestions and feedback.
(Sept. 2016 - Dec. 2016)
David Bauer¶
- RWTH Aachen
- PGI-1/IAS-1 at Forschungszentrum Jülich
Initial (Fortran90) implementations of the isotropic Heisenberg Hamiltonian, Neighbour calculations and the SIB solver.
(Oct. 2007 - Sept. 2008)
Reference¶
The Spirit framework is a scientific project. If you present and/or publish scientific results or visualisations that used Spirit, you should add a reference.
The Framework¶
If you used any components of this framework please add a reference to our GitHub page. You may use e.g. the following TeX code:
\bibitem{spirit}
{Spirit spin simulation framework} (see spirit-code.github.io)
Specific Methods¶
The following need only be cited if used.
Depondt Solver
This Heun-like method for solving the LLG equation including the stochastic term has been published by Depondt et al.: http://iopscience.iop.org/0953-8984/21/33/336005 You may use e.g. the following TeX code:
\bibitem{Depondt}
Ph. Depondt et al. \textit{J. Phys. Condens. Matter} \textbf{21}, 336005 (2009).
SIB Solver
This stable method for solving the LLG equation efficiently and including the stochastic term has been published by Mentink et al.: http://iopscience.iop.org/0953-8984/22/17/176001 You may use e.g. the following TeX code:
\bibitem{SIB}
J. H. Mentink et al. \textit{J. Phys. Condens. Matter} \textbf{22}, 176001 (2010).
VP Solver
This intuitive direct minimization routine has been published as supplementary material by Bessarab et al.: http://www.sciencedirect.com/science/article/pii/S0010465515002696 You may use e.g. the following TeX code:
\bibitem{VP}
P. F. Bessarab et al. \textit{Comp. Phys. Comm.} \textbf{196}, 335 (2015).
GNEB Method
This specialized nudged elastic band method for calculating transition paths of spin systems has been published by Bessarab et al.: http://www.sciencedirect.com/science/article/pii/S0010465515002696 You may use e.g. the following TeX code:
\bibitem{GNEB}
P. F. Bessarab et al. \textit{Comp. Phys. Comm.} \textbf{196}, 335 (2015).
Included Dependencies¶
Vector Field Rendering¶
libvfrendering is a C++ library for rendering vectorfields using OpenGL. Originally developed for spirit and based on WegGLSpins.js, it has an extendable architecture and currently offers renderer implementations for: - glyph-based vector field representations as arrows - colormapped surface and isosurface rendering - mapping vector directions onto a sphere
The library is still very much a work-in-progress, so its API is not yet stable and there are still several features missing that will be added in later releases. If you miss a feature or have another idea on how to improve libvfrendering, please open an issue or pull request!

Getting Started¶
To use libvfrendering, you need to perform the following steps:
- Include <VFRendering/View.hxx>
- Create a VFRendering::Geometry
- Read or calculate the vector directions
- Pass geometry and directions to a VFRendering::View
- Draw the view in an existing OpenGL context
1. Include <VFRendering/View.hxx>¶
When using libvfrendering, you will mostly interact with View
objects, so it should be enough to #include <VFRendering/View.hxx>
.
2. Create a VFRendering::Geometry¶
The geometry describes the positions on which you evaluated the vector field and how they might form a grid (optional, e.g. for isosurface and surface rendering). You can pass the positions directly to the constructor or call one of the class’ static methods.
As an example, this is how you could create a simple, cartesian 30x30x30 geometry, with coordinates between -1 and 1:
auto geometry = VFRendering::Geometry::cartesianGeometry(
{30, 30, 30},
{-1.0, -1.0, -1.0},
{1.0, 1.0, 1.0}
);
3. Read or calculate the vector directions¶
This step highly depends on your use case. The directions are stored as a ``std::vector<glm::vec3>``, so they can be created in a simple loop:
std::vector<glm::vec3> directions;
for (int iz = 0; iz < 10; iz++) {
for (int iy = 0; iy < 10; iy++) {
for (int ix = 0; ix < 10; ix++) {
// calculate direction for ix, iy, iz
directions.push_back(glm::normalize({ix-4.5, iy-4.5, iz-4.5}));
}
}
}
As shown here, the directions should be in C order when using the
VFRendering::Geometry
static methods. If you do not know
glm, think of a glm::vec3
as a struct
containing three floats x, y and z.
4. Create a VFRendering::VectorField¶
This class simply contains geometry and directions.
VFRendering::VectorField vf(geometry, directions);
To update the VectorField data, use VectorField::update
. If the
directions changed but the geometry is the same, you can use the
VectorField::updateVectors
method or VectorField::updateGeometry
vice versa.
5. Create a VFRendering::View and a Renderer¶
The view object is what you will interact most with. It provides an interface to the various renderers and includes functions for handling mouse input.
You can create a new view and then initialize the renderer(s)
(as an example, we use the VFRendering::ArrowRenderer
):
VFRendering::View view;
auto arrow_renderer_ptr = std::make_shared<VFRendering::ArrowRenderer>(view, vf);
view.renderers( {{ arrow_renderer_ptr, {0, 0, 1, 1} }} );
5. Draw the view in an existing OpenGL context¶
To actually see something, you need to create an OpenGL context using a toolkit of your choice, e.g. Qt or GLFW. After creating the context, pass the framebuffer size to the setFramebufferSize method. You can then call the draw method of the view to render the vector field, either in a loop or only when you update the data.
view.draw();
For a complete example, including an interactive camera, see demo.cxx.
Python Package¶
The Python package has bindings which correspond directly to the C++ class and function names. To use pyVFRendering, you need to perform the following steps:
import pyVFRendering as vfr
- Create a
vfr.Geometry
- Read or calculate the vector directions
- Pass geometry and directions to a
vfr.View
- Draw the view in an existing OpenGL context
1. import¶
In order to import pyVFRendering as vfr
, you can either
pip install pyVFRendering
or download and build it yourself.
You can build with python3 setup.py build
, which will generate a
library somewhere in your build
subfolder, which you can import
in python. Note that you may need to add the folder to your
PYTHONPATH
.
2. Create a pyVFRendering.Geometry¶
As above:
geometry = vfr.Geometry.cartesianGeometry(
(30, 30, 30), # number of lattice points
(-1.0, -1.0, -1.0), # lower bound
(1.0, 1.0, 1.0) ) # upper bound
3. Read or calculate the vector directions¶
This step highly depends on your use case. Example:
directions = []
for iz in range(n_cells[2]):
for iy in range(n_cells[1]):
for ix in range(n_cells[0]):
# calculate direction for ix, iy, iz
directions.append( [ix-4.5, iy-4.5, iz-4.5] )
4. Pass geometry and directions to a pyVFRendering.View¶
You can create a new view and then pass the geometry and directions by calling the update method:
view = vfr.View()
view.update(geometry, directions)
If the directions changed but the geometry is the same, you can use the updateVectors method.
5. Draw the view in an existing OpenGL context¶
To actually see something, you need to create an OpenGL context using a framework of your choice, e.g. Qt or GLFW. After creating the context, pass the framebuffer size to the setFramebufferSize method. You can then call the draw method of the view to render the vector field, either in a loop or only when you update the data.
view.setFramebufferSize(width*self.window().devicePixelRatio(), height*self.window().devicePixelRatio())
view.draw()
For a complete example, including an interactive camera, see demo.py.
Renderers¶
libvfrendering offers several types of renderers, which all inherit
from VFRendering::RendererBase
. The most relevant are the
VectorFieldRenderer
s:
- VFRendering::ArrowRenderer, which renders the vectors as colored arrows
- VFRendering::SphereRenderer, which renders the vectors as colored spheres
- VFRendering::SurfaceRenderer, which renders the surface of the geometry using a colormap
- VFRendering::IsosurfaceRenderer, which renders an isosurface of the vectorfield using a colormap
- VFRendering::VectorSphereRenderer, which renders the vectors as dots on a sphere, with the position of each dot representing the direction of the vector
In addition to these, there also the following renderers which do not
require a VectorField
: - VFRendering::CombinedRenderer, which can be
used to create a combination of several renderers, like an isosurface
rendering with arrows - VFRendering::BoundingBoxRenderer, which is used
for rendering bounding boxes around the geometry rendered by an
VFRendering::ArrorRenderer, VFRendering::SurfaceRenderer or
VFRendering::IsosurfaceRenderer - VFRendering::CoordinateSystemRenderer,
which is used for rendering a coordinate system, with the axes colored
by using the colormap
To control what renderers are used, you can use
VFRendering::View::renderers
, where you can pass it a
std::vector
s of std::pair
s of renderers as
std::shared_ptr<VFRendering::RendererBase>
(i.e. shared pointers)
and viewports as glm::vec4
.
Options¶
To modify the way the vector field is rendered, libvfrendering offers a variety of options. To set these, you can create an VFRendering::Options object.
As an example, to adjust the vertical field of view, you would do the following:
VFRendering::Options options;
options.set<VFRendering::View::Option::VERTICAL_FIELD_OF_VIEW>(30);
view.updateOptions(options);
If you want to set only one option, you can also use View::setOption:
view.setOption<VFRendering::View::Option::VERTICAL_FIELD_OF_VIEW>(30);
If you want to set an option for an individual Renderer, you can use the methods RendererBase::updateOptions and RendererBase::setOption in the same way.
Whether this way of setting options should be replaced by getters/setters will be evaluated as the API becomes more stable.
Currently, the following options are available:
Index | Type | Default value | Header file | Documentation |
---|---|---|---|---|
View::Op tion::BO UNDING_ BOX_MIN | glm::vec 3 | {-1, -1, -1} | View.hxx | VFRendering::Uti lities::Options: :Option< View::Option::BO UNDING_BOX_MIN > |
View::Op tion::BO UNDING_ BOX_MAX | glm::vec 3 | {1, 1, 1} | View.hxx | VFRendering::Uti lities::Options: :Option< View::Option::BO UNDING_BOX_MAX > |
View::Op tion::SY STEM_CE NTER | glm::vec 3 | {0, 0, 0} | View.hxx | VFRendering::Uti lities::Options: :Option< View::Option::SY STEM_CENTER > |
View::Op tion::VE RTICAL_ FIELD_O F_VIEW | float | 45.0 | View.hxx | VFRendering::Uti lities::Options: :Option< View::Option::VE RTICAL_FIELD_O F_VIEW > |
View::Op tion::BA CKGROUND _COLOR | glm::vec 3 | {0, 0, 0} | View.hxx | VFRendering::Uti lities::Options: :Option< View::Option::BA CKGROUND_COLOR > |
View::Op tion::CO LORMAP_ IMPLEMEN TATION | std::str ing | VFRendering::Uti lities::getColor mapImplementatio n(VFRendering::U tilities::Colorm ap::DEFAULT) | View.hxx | VFRendering::Uti lities::Options: :Option< View::Option::CO LORMAP_IMPLEMEN TATION > |
View::Op tion::IS _VISIBL E_IMPLE MENTATIO N | std::str ing | bool is_visible(vec3 position, vec3 direction) { return true; } | View.hxx | VFRendering::Uti lities::Options: :Option< View::Option::IS _VISIBLE_IMPLE MENTATION > |
View::Op tion::CA MERA_PO SITION | glm::vec 3 | {14.5, 14.5, 30} | View.hxx | VFRendering::Uti lities::Options: :Option< View::Option::CA MERA_POSITION > |
View::Op tion::CE NTER_PO SITION | glm::vec 3 | {14.5, 14.5, 0} | View.hxx | VFRendering::Uti lities::Options: :Option< View::Option::CE NTER_POSITION > |
View::Op tion::UP _VECTOR | glm::vec 3 | {0, 1, 0} | View.hxx | VFRendering::Uti lities::Options: :Option< View::Option::UP _VECTOR > |
ArrowRen derer::O ption::C ONE_RAD IUS | float | 0.25 | ArrowRenderer. hxx | VFRendering::Uti lities::Options: :Option< ArrowRenderer::O ption::CONE_RAD IUS > |
ArrowRen derer::O ption::C ONE_HEI GHT | float | 0.6 | ArrowRenderer. hxx | VFRendering::Uti lities::Options: :Option< ArrowRenderer::O ption::CONE_HEI GHT > |
ArrowRen derer::O ption::C YLINDER_RADIUS | float | 0.125 | ArrowRenderer. hxx | VFRendering::Uti lities::Options: :Option< ArrowRenderer::O ption::CYLINDER_RADIUS > |
ArrowRen derer::O ption::C YLINDER_HEIGHT | float | 0.7 | ArrowRenderer. hxx | VFRendering::Uti lities::Options: :Option< ArrowRenderer::O ption::CYLINDER_HEIGHT > |
ArrowRen derer::O ption::L EVEL_OF _DETAIL | unsigned int | 20 | ArrowRenderer. hxx | VFRendering::Uti lities::Options: :Option< ArrowRenderer::O ption::LEVEL_OF _DETAIL > |
Bounding BoxRende rer::Opt ion::COL OR | glm::vec 3 | {1.0, 1.0, 1.0} | BoundingBoxRen derer.hxx | VFRendering::Uti lities::Options: :Option< BoundingBoxRende rer::Option::COL OR > |
Coordina teSystem Renderer ::Option ::AXIS_ LENGTH | glm::vec 3 | {0.5, 0.5, 0.5} | CoordinateSyst emRenderer.hxx | VFRendering::Uti lities::Options: :Option< CoordinateSystem Renderer::Option ::AXIS_LENGTH > |
Coordina teSystem Renderer ::Option ::ORIGIN | glm::vec 3 | {0.0, 0.0, 0.0} | CoordinateSyst emRenderer.hxx | VFRendering::Uti lities::Options: :Option< CoordinateSystem Renderer::Option ::ORIGIN > |
Coordina teSystem Renderer ::Option ::NORMAL IZE | bool | false | CoordinateSyst emRenderer.hxx | VFRendering::Uti lities::Options: :Option< CoordinateSystem Renderer::Option ::NORMALIZE > |
Coordina teSystem Renderer ::Option ::LEVEL_OF_DET AIL | unsigned int | 100 | CoordinateSyst emRenderer.hxx | VFRendering::Uti lities::Options: :Option< CoordinateSystem Renderer::Option ::LEVEL_OF_DET AIL > |
Coordina teSystem Renderer ::Option ::CONE_ HEIGHT | float | 0.3 | CoordinateSyst emRenderer.hxx | VFRendering::Uti lities::Options: :Option< CoordinateSystem Renderer::Option ::CONE_HEIGHT > |
Coordina teSystem Renderer ::Option ::CONE_ RADIUS | float | 0.1 | CoordinateSyst emRenderer.hxx | VFRendering::Uti lities::Options: :Option< CoordinateSystem Renderer::Option ::CONE_RADIUS > |
Coordina teSystem Renderer ::Option ::CYLIND ER_HEIG HT | float | 0.7 | CoordinateSyst emRenderer.hxx | VFRendering::Uti lities::Options: :Option< CoordinateSystem Renderer::Option ::CYLINDER_HEIG HT > |
Coordina teSystem Renderer ::Option ::CYLIND ER_RADI US | float | 0.07 | CoordinateSyst emRenderer.hxx | VFRendering::Uti lities::Options: :Option< CoordinateSystem Renderer::Option ::CYLINDER_RADI US > |
Isosurfa ceRender er::Opti on::ISOV ALUE | float | 0.0 | IsosurfaceRend erer.hxx | VFRendering::Uti lities::Options: :Option< IsosurfaceRender er::Option::ISOV ALUE > |
Isosurfa ceRender er::Opti on::LIGH TING_IM PLEMENTA TION | std::str ing | float lighting(vec3 position, vec3 normal) { return 1.0; } | IsosurfaceRend erer.hxx | VFRendering::Uti lities::Options: :Option< IsosurfaceRender er::Option::LIGH TING_IMPLEMENTA TION > |
Isosurfa ceRender er::Opti on::VALU E_FUNCT ION | std::fun ction | [] (const glm::vec3& position, const glm::vec3& direction) { return direction.z; } | IsosurfaceRend erer.hxx | VFRendering::Uti lities::Options: :Option< IsosurfaceRender er::Option::VALU E_FUNCTION > |
VectorSp hereRend erer::Op tion::PO INT_SIZ E_RANGE | glm::vec 2 | {1.0, 4.0} | VectorSphereRe nderer.hxx | VFRendering::Uti lities::Options: :Option< VectorSphereRend erer::Option::PO INT_SIZE_RANGE > |
VectorSp hereRend erer::Op tion::IN NER_SPH ERE_RAD IUS | float | 0.95 | VectorSphereRe nderer.hxx | VFRendering::Uti lities::Options: :Option< VectorSphereRend erer::Option::IN NER_SPHERE_RAD IUS > |
VectorSp hereRend erer::Op tion::US E_SPHER E_FAKE_PERSPEC TIVE | bool | true | VectorSphereRe nderer.hxx | VFRendering::Uti lities::Options: :Option< VectorSphereRend erer::Option::US E_SPHERE_FAKE_PERSPECTIVE > |
ToDo¶
- A EGS plugin for combining libvfrendering with existing EGS plugins.
- Methods for reading geometry and directions from data files
- Documentation
See the issues for further information and adding your own requests.
Eigen¶
Eigen is a C++ template library for linear algebra: matrices, vectors, numerical solvers, and related algorithms.
For more information go to http://eigen.tuxfamily.org/.
Spectra¶
**Spectra** stands for Sparse Eigenvalue Computation Toolkit as a Redesigned ARPACK. It is a C++ library for large scale eigenvalue problems, built on top of Eigen, an open source linear algebra library.
Spectra is implemented as a header-only C++ library, whose only dependency, Eigen, is also header-only. Hence Spectra can be easily embedded in C++ projects that require calculating eigenvalues of large matrices.
Relation to ARPACK¶
ARPACK is a software written in FORTRAN for solving large scale eigenvalue problems. The development of Spectra is much inspired by ARPACK, and as the whole name indicates, Spectra is a redesign of the ARPACK library using C++ language.
In fact, Spectra is based on the algorithms described in the ARPACK Users’ Guide, but it does not use the ARPACK code, and it is NOT a clone of ARPACK for C++. In short, Spectra implements the major algorithms in ARPACK, but Spectra provides a completely different interface, and it does not depend on ARPACK.
Common Usage¶
Spectra is designed to calculate a specified number (k
) of
eigenvalues of a large square matrix (A
). Usually k
is much less
than the size of matrix (n
), so that only a few eigenvalues and
eigenvectors are computed, which in general is more efficient than
calculating the whole spectral decomposition. Users can choose
eigenvalue selection rules to pick up the eigenvalues of interest, such
as the largest k
eigenvalues, or eigenvalues with largest real
parts, etc.
To use the eigen solvers in this library, the user does not need to
directly provide the whole matrix, but instead, the algorithm only
requires certain operations defined on A
, and in the basic setting,
it is simply the matrix-vector multiplication. Therefore, if the
matrix-vector product A * x
can be computed efficiently, which is
the case when A
is sparse, Spectra will be very powerful for
large scale eigenvalue problems.
There are two major steps to use the Spectra library:
- Define a class that implements a certain matrix operation, for
example the matrix-vector multiplication
y = A * x
or the shift-solve operationy = inv(A - σ * I) * x
. Spectra has defined a number of helper classes to quickly create such operations from a matrix object. See the documentation of DenseGenMatProd, DenseSymShiftSolve, etc. - Create an object of one of the eigen solver classes, for example SymEigsSolver for symmetric matrices, and GenEigsSolver for general matrices. Member functions of this object can then be called to conduct the computation and retrieve the eigenvalues and/or eigenvectors.
Below is a list of the available eigen solvers in Spectra: - SymEigsSolver: For real symmetric matrices - GenEigsSolver: For general real matrices - SymEigsShiftSolver: For real symmetric matrices using the shift-and-invert mode - GenEigsRealShiftSolver: For general real matrices using the shift-and-invert mode, with a real-valued shift - GenEigsComplexShiftSolver: For general real matrices using the shift-and-invert mode, with a complex-valued shift - SymGEigsSolver: For generalized eigen solver for real symmetric matrices
Examples¶
Below is an example that demonstrates the use of the eigen solver for symmetric matrices.
#include <Eigen/Core>
#include <SymEigsSolver.h> // Also includes <MatOp/DenseSymMatProd.h>
#include <iostream>
using namespace Spectra;
int main()
{
// We are going to calculate the eigenvalues of M
Eigen::MatrixXd A = Eigen::MatrixXd::Random(10, 10);
Eigen::MatrixXd M = A + A.transpose();
// Construct matrix operation object using the wrapper class DenseGenMatProd
DenseSymMatProd<double> op(M);
// Construct eigen solver object, requesting the largest three eigenvalues
SymEigsSolver< double, LARGEST_ALGE, DenseSymMatProd<double> > eigs(&op, 3, 6);
// Initialize and compute
eigs.init();
int nconv = eigs.compute();
// Retrieve results
Eigen::VectorXd evalues;
if(eigs.info() == SUCCESSFUL)
evalues = eigs.eigenvalues();
std::cout << "Eigenvalues found:\n" << evalues << std::endl;
return 0;
}
Sparse matrix is supported via the SparseGenMatProd
class.
#include <Eigen/Core>
#include <Eigen/SparseCore>
#include <GenEigsSolver.h>
#include <MatOp/SparseGenMatProd.h>
#include <iostream>
using namespace Spectra;
int main()
{
// A band matrix with 1 on the main diagonal, 2 on the below-main subdiagonal,
// and 3 on the above-main subdiagonal
const int n = 10;
Eigen::SparseMatrix<double> M(n, n);
M.reserve(Eigen::VectorXi::Constant(n, 3));
for(int i = 0; i < n; i++)
{
M.insert(i, i) = 1.0;
if(i > 0)
M.insert(i - 1, i) = 3.0;
if(i < n - 1)
M.insert(i + 1, i) = 2.0;
}
// Construct matrix operation object using the wrapper class SparseGenMatProd
SparseGenMatProd<double> op(M);
// Construct eigen solver object, requesting the largest three eigenvalues
GenEigsSolver< double, LARGEST_MAGN, SparseGenMatProd<double> > eigs(&op, 3, 6);
// Initialize and compute
eigs.init();
int nconv = eigs.compute();
// Retrieve results
Eigen::VectorXcd evalues;
if(eigs.info() == SUCCESSFUL)
evalues = eigs.eigenvalues();
std::cout << "Eigenvalues found:\n" << evalues << std::endl;
return 0;
}
And here is an example for user-supplied matrix operation class.
#include <Eigen/Core>
#include <SymEigsSolver.h>
#include <iostream>
using namespace Spectra;
// M = diag(1, 2, ..., 10)
class MyDiagonalTen
{
public:
int rows() { return 10; }
int cols() { return 10; }
// y_out = M * x_in
void perform_op(double *x_in, double *y_out)
{
for(int i = 0; i < rows(); i++)
{
y_out[i] = x_in[i] * (i + 1);
}
}
};
int main()
{
MyDiagonalTen op;
SymEigsSolver<double, LARGEST_ALGE, MyDiagonalTen> eigs(&op, 3, 6);
eigs.init();
eigs.compute();
if(eigs.info() == SUCCESSFUL)
{
Eigen::VectorXd evalues = eigs.eigenvalues();
std::cout << "Eigenvalues found:\n" << evalues << std::endl;
}
return 0;
}
Shift-and-invert Mode¶
When we want to find eigenvalues that are closest to a number σ
, for
example to find the smallest eigenvalues of a positive definite matrix
(in which case σ = 0
), it is advised to use the shift-and-invert
mode of eigen solvers.
In the shift-and-invert mode, selection rules are applied to
1/(λ - σ)
rather than λ
, where λ
are eigenvalues of A
.
To use this mode, users need to define the shift-solve matrix operation.
See the documentation of
SymEigsShiftSolver
for details.
Documentation¶
The API reference page contains the documentation of Spectra generated by Doxygen, including all the background knowledge, example code and class APIs.
More information can be found in the project page http://yixuan.cos.name/spectra.
{fmt}¶

fmt is an open-source formatting library for C++. It can be used as a safe alternative to printf or as a fast alternative to IOStreams.
Features¶
- Two APIs: faster concatenation-based write API and slower, but still very fast, replacement-based format API with positional arguments for localization.
- Write API similar to the one used by IOStreams but stateless allowing faster implementation.
- Format API with format string syntax similar to the one used by str.format in Python.
- Safe printf implementation including the POSIX extension for positional arguments.
- Support for user-defined types.
- High speed: performance of the format API is close to that of glibc’s printf and better than the performance of IOStreams. See Speed tests and Fast integer to string conversion in C++.
- Small code size both in terms of source code (the core library consists of a single header file and a single source file) and compiled code. See Compile time and code bloat.
- Reliability: the library has an extensive set of unit tests.
- Safety: the library is fully type safe, errors in format strings are reported using exceptions, automatic memory management prevents buffer overflow errors.
- Ease of use: small self-contained code base, no external dependencies, permissive BSD license
- Portability with consistent output across platforms and support for older compilers.
- Clean warning-free codebase even on high warning levels (-Wall -Wextra -pedantic).
- Support for wide strings.
- Optional header-only configuration enabled with the
FMT_HEADER_ONLY
macro.
See the documentation for more details.
Examples¶
This prints Hello, world!
to stdout:
fmt::print("Hello, {}!", "world"); // uses Python-like format string syntax
fmt::printf("Hello, %s!", "world"); // uses printf format string syntax
Arguments can be accessed by position and arguments’ indices can be repeated:
std::string s = fmt::format("{0}{1}{0}", "abra", "cad");
// s == "abracadabra"
fmt can be used as a safe portable replacement for itoa
:
fmt::MemoryWriter w;
w << 42; // replaces itoa(42, buffer, 10)
w << fmt::hex(42); // replaces itoa(42, buffer, 16)
// access the string using w.str() or w.c_str()
An object of any user-defined type for which there is an overloaded
std::ostream
insertion operator (operator<<
) can be formatted:
#include "fmt/ostream.h"
class Date {
int year_, month_, day_;
public:
Date(int year, int month, int day) : year_(year), month_(month), day_(day) {}
friend std::ostream &operator<<(std::ostream &os, const Date &d) {
return os << d.year_ << '-' << d.month_ << '-' << d.day_;
}
};
std::string s = fmt::format("The date is {}", Date(2012, 12, 9));
// s == "The date is 2012-12-9"
You can use the FMT_VARIADIC macro to create your own functions similar to format and print which take arbitrary arguments:
// Prints formatted error message.
void report_error(const char *format, fmt::ArgList args) {
fmt::print("Error: ");
fmt::print(format, args);
}
FMT_VARIADIC(void, report_error, const char *)
report_error("file not found: {}", path);
Note that you only need to define one function that takes fmt::ArgList
argument. FMT_VARIADIC
automatically defines necessary wrappers that
accept variable number of arguments.
Projects using this library¶
- 0 A.D.: A free, open-source, cross-platform real-time strategy game
- AMPL/MP: An open-source library for mathematical programming
- CUAUV: Cornell University’s autonomous underwater vehicle
- Drake: A planning, control, and analysis toolbox for nonlinear dynamical systems (MIT)
- Envoy: C++ L7 proxy and communication bus (Lyft)
- FiveM: a modification framework for GTA V
- HarpyWar/pvpgn: Player vs Player Gaming Network with tweaks
- KBEngine: An open-source MMOG server engine
- Keypirinha: A semantic launcher for Windows
- Kodi (formerly xbmc): Home theater software
- Lifeline: A 2D game
- MongoDB Smasher: A small tool to generate randomized datasets
- OpenSpace: An open-source astrovisualization framework
- PenUltima Online (POL): An MMO server, compatible with most Ultima Online clients
- quasardb: A distributed, high-performance, associative database
- readpe: Read Portable Executable
- redis-cerberus: A Redis cluster proxy
- Saddy: Small crossplatform 2D graphic engine
- Salesforce Analytics Cloud: Business intelligence software
- Scylla: A Cassandra-compatible NoSQL data store that can handle 1 million transactions per second on a single server
- Seastar: An advanced, open-source C++ framework for high-performance server applications on modern hardware
- spdlog: Super fast C++ logging library
- Stellar: Financial platform
- Touch Surgery: Surgery simulator
- TrinityCore: Open-source MMORPG framework
If you are aware of other projects using this library, please let me know by email or by submitting an issue.
Motivation¶
So why yet another formatting library?
There are plenty of methods for doing this task, from standard ones like the printf family of function and IOStreams to Boost Format library and FastFormat. The reason for creating a new library is that every existing solution that I found either had serious issues or didn’t provide all the features I needed.
Printf¶
The good thing about printf is that it is pretty fast and readily available being a part of the C standard library. The main drawback is that it doesn’t support user-defined types. Printf also has safety issues although they are mostly solved with __attribute__ ((format (printf, …)) in GCC. There is a POSIX extension that adds positional arguments required for i18n to printf but it is not a part of C99 and may not be available on some platforms.
IOStreams¶
The main issue with IOStreams is best illustrated with an example:
std::cout << std::setprecision(2) << std::fixed << 1.23456 << "\n";
which is a lot of typing compared to printf:
printf("%.2f\n", 1.23456);
Matthew Wilson, the author of FastFormat, referred to this situation with IOStreams as “chevron hell”. IOStreams doesn’t support positional arguments by design.
The good part is that IOStreams supports user-defined types and is safe although error reporting is awkward.
Boost Format library¶
This is a very powerful library which supports both printf-like format strings and positional arguments. The main its drawback is performance. According to various benchmarks it is much slower than other methods considered here. Boost Format also has excessive build times and severe code bloat issues (see Benchmarks).
FastFormat¶
This is an interesting library which is fast, safe and has positional arguments. However it has significant limitations, citing its author:
Three features that have no hope of being accommodated within the current design are:
- Leading zeros (or any other non-space padding)
- Octal/hexadecimal encoding
- Runtime width/alignment specification
It is also quite big and has a heavy dependency, STLSoft, which might be too restrictive for using it in some projects.
Loki SafeFormat¶
SafeFormat is a formatting library which uses printf-like format strings
and is type safe. It doesn’t support user-defined types or positional
arguments. It makes unconventional use of operator()
for passing
format arguments.
Tinyformat¶
This library supports printf-like format strings and is very small and fast. Unfortunately it doesn’t support positional arguments and wrapping it in C++98 is somewhat difficult. Also its performance and code compactness are limited by IOStreams.
Boost Spirit.Karma¶
This is not really a formatting library but I decided to include it here
for completeness. As IOStreams it suffers from the problem of mixing
verbatim text with arguments. The library is pretty fast, but slower
on integer formatting than fmt::Writer
on Karma’s own benchmark,
see Fast integer to string conversion in C++.
Benchmarks¶
Speed tests¶
The following speed tests results were generated by building
tinyformat_test.cpp
on Ubuntu GNU/Linux 14.04.1 with
g++-4.8.2 -O3 -DSPEED_TEST -DHAVE_FORMAT
, and taking the best of three
runs. In the test, the format string "%0.10f:%04d:%+g:%s:%p:%c:%%\n"
or
equivalent is filled 2000000 times with output sent to /dev/null
; for
further details see the source.
Library | Method | Run Time, s |
---|---|---|
EGLIBC 2.19 | printf | 1.30 |
libstdc++ 4.8.2 | std::ostream | 1.85 |
fmt 1.0 | fmt::print | 1.42 |
tinyformat 2.0.1 | tfm::printf | 2.25 |
Boost Format 1.54 | boost::format | 9.94 |
As you can see boost::format
is much slower than the alternative methods; this
is confirmed by other tests.
Tinyformat is quite good coming close to IOStreams. Unfortunately tinyformat
cannot be faster than the IOStreams because it uses them internally.
Performance of fmt is close to that of printf, being faster than printf on integer
formatting,
but slower on floating-point formatting which dominates this benchmark.
Compile time and code bloat¶
The script bloat-test.py
from format-benchmark
tests compile time and code bloat for nontrivial projects.
It generates 100 translation units and uses printf()
or its alternative
five times in each to simulate a medium sized project. The resulting
executable size and compile time (g++-4.8.1, Ubuntu GNU/Linux 13.10,
best of three) is shown in the following tables.
Optimized build (-O3)
Method | Compile Time, s | Executable size, KiB | Stripped size, KiB |
---|---|---|---|
printf | 2.6 | 41 | 30 |
IOStreams | 19.4 | 92 | 70 |
fmt | 46.8 | 46 | 34 |
tinyformat | 64.6 | 418 | 386 |
Boost Format | 222.8 | 990 | 923 |
As you can see, fmt has two times less overhead in terms of resulting
code size compared to IOStreams and comes pretty close to printf
.
Boost Format has by far the largest overheads.
Non-optimized build
Method | Compile Time, s | Executable size, KiB | Stripped size, KiB |
---|---|---|---|
printf | 2.1 | 41 | 30 |
IOStreams | 19.7 | 86 | 62 |
fmt | 47.9 | 108 | 86 |
tinyformat | 27.7 | 234 | 190 |
Boost Format | 122.6 | 884 | 763 |
libc
, libstdc++
and libfmt
are all linked as shared
libraries to compare formatting function overhead only. Boost Format
and tinyformat are header-only libraries so they don’t provide any
linkage options.
Running the tests¶
Please refer to Building the library for the instructions on how to build the library and run the unit tests.
Benchmarks reside in a separate repository, format-benchmarks, so to run the benchmarks you first need to clone this repository and generate Makefiles with CMake:
$ git clone --recursive https://github.com/fmtlib/format-benchmark.git
$ cd format-benchmark
$ cmake .
Then you can run the speed test:
$ make speed-test
or the bloat test:
$ make bloat-test
License¶
fmt is distributed under the BSD license.
The Format String Syntax section in the documentation is based on the one from Python string module documentation adapted for the current library. For this reason the documentation is distributed under the Python Software Foundation license available in doc/python-license.txt. It only applies if you distribute the documentation of fmt.
Acknowledgments¶
The fmt library is maintained by Victor Zverovich (vitaut) and Jonathan Müller (foonathan) with contributions from many other people. See Contributors and Releases for some of the names. Let us know if your contribution is not listed or mentioned incorrectly and we’ll make it right.
The benchmark section of this readme file and the performance tests are taken from the excellent tinyformat library written by Chris Foster. Boost Format library is acknowledged transitively since it had some influence on tinyformat. Some ideas used in the implementation are borrowed from Loki SafeFormat and Diagnostic API in Clang. Format string syntax and the documentation are based on Python’s str.format. Thanks Doug Turnbull for his valuable comments and contribution to the design of the type-safe API and Gregory Czajkowski for implementing binary formatting. Thanks Ruslan Baratov for comprehensive comparison of integer formatting algorithms and useful comments regarding performance, Boris Kaul for C++ counting digits benchmark. Thanks to CarterLi for contributing various improvements to the code.