-
Jonas Greitemann authoredJonas Greitemann authored
Probing hidden spin order with interpretable machine learning
These are the source codes accompanying the papers:
- Jonas Greitemann, Ke Liu, and Lode Pollet: Probing hidden spin order with interpretable machine learning, Phys. Rev. B 99, 060404(R) (2019), open access via arXiv:1804.08557;
- Ke Liu, Jonas Greitemann, and Lode Pollet: Learning multiple order parameters with interpretable machines, Phys. Rev. B 99, 104410 (2019), open access via arXiv:1810.05538;
- Jonas Greitemann, Ke Liu, Ludovic D.C. Jaubert, Han Yan, Nic Shannon, and Lode Pollet: Identification of emergent constraints and hidden order in frustrated magnets using tensorial kernel methods of machine learning, Phys. Rev. B 100, 174408 (2019), open access via arXiv:1907.12322.
They aim to provide a reusable framework that may be used to apply the TK-SVM
method to third-party Monte Carlo simulations written using ALPSCore (→ Project
structure). The gauge theory code used for the
first two papers is included as part of this repository. The simulation of the
XXZ pyrochlore model for the third paper is not openly available, yet the
present version of the framework contains many of the supporting features. While
this version is still capable of producing the results presented in the first
two papers, you may want to check out the versions of the code that were
originally published along with the first two papers. The relevant commits are
tagged PRB-99-060404(R)
and PRB-99-104410
, respectively, i.e. run
$ git checkout PRB-99-060404(R)
$ git checkout PRB-99-104410
Additionally J. Greitemann's PhD thesis includes results for the simulation of the Kagome antiferromagnet which may also be reproduced using the present version of the code.
Contents
- Dependencies
- Project structure
- Building and installation
- Basic usage
- Runtime parameters
- Client code API
- License
A Changelog is provided separately.
Dependencies
These codes are based on the ALPSCore library. Refer to their website for installation instructions. ALPSCore requires the Boost and HDF5 libraries.
Further, we require a C++ compiler (tested with GCC 9 and Clang 8) with
C++14 support. For the solution of the SVM optimization problem, we rely on the
libsvm library which is included as part of our self-developed C++ wrapper
library which is kept in-tree as a git submodule. Argh! is used for
parsing command line options. Eigen 3 is used for linear algebra.
As of version 3.0, MPI is used for parallelization. A thread-safe MPI
implementation is hence required (MPI_THREAD_MULTIPLE
build option enabled;
tested with Open MPI 4.0 and MPICH 3.3).
Project structure
This repository contains three client codes which are for the most part regular Monte Carlo simulation codes based on the ALPSCore framework. In addition they implement the relevant API functions required by the TK-SVM framework.
The client codes are kept in their eponymous subdirectories:
-
ising
: This is a copy of the ALPSCore demo code for the 2D Ising model. It may be used to reproduce the results presented in P. Ponte and R.G. Melko, Phys. Rev. B 96, 205146 (2017). -
gauge
: Lattice gauge model proposed in K. Liu, J. Nissinen, R.-J. Slager, K. Wu, and J. Zaanen, PRX 6, 041025 (2016) which lets us realize arbitrary orientational orders and was used to produce the results in our first two papers. -
frustmag
: A generic code for classical Monte Carlo simulations of models of frustrated magnets. Several lattice geometries are available to choose from, additional interactions may be included with relative ease. Supports several update schemes, including parallel tempering. This implementation was used to generate data for the study of the Kagome antiferromagnet in Ch. 9 of J. Greitemann's PhD thesis.
Additional information on the runtime parameters pertaining to these client codes are documented in the README files accompanying these respective subdirectories.
This project is designed as a generic framework that works with existing Monte
Carlo codes with minimal adaptation. Any code outside the client-code
directories is designed to be agnostic with respect to the client codes. For
each of the client codes, several executables are generated, replacing the
main
function of a typical ALPSCore simulation and ending in the following
suffixes:
-
*-sample
: Run the Monte Carlo simulations at different points in the phase diagram, occasionally sampling the spin configuration and optionally applying the monomial mapping. The point in phase diagram space is recorded along with the raw (inLAZY
mode) or mapped (inEAGER
mode) configuration. The resulting optimization problem is saved in a checkpoint file ending in*.clone.h5
. -
*-learn
: Read one or more checkpoint files (*.clone.h5
), label the samples according to the chosen → Classifier, and perform the SVM optimization. The resultant SVM model is saved to disk in a file ending in*.out.h5
. -
*-test
: Run the Monte Carlo simulation along a set of parameter points in the phase diagram. If launched from an*.out.h5
file, the SVM model is evaluated and the decision function and its variance are measured as observables. Additionally, any other predefined reference observables are measured. If launched from an*.ini
file instead, only the latter observables are measured. -
*-coeffs
: Extract the coefficient matrix from the SVM model to infer the analytical order parameter. -
*-segregate-phases
: Carry out the graph analysis to find the phase diagram without prior knowledge of the correct phase classification.
Building and installation
Building and installing ALPSCore
Detailed instructions on how to build ALPSCore can be fournd in the project's wiki. The procedure revolves around the following:
$ cd alpscore
$ mkdir build.tmp && cd build.tmp
$ cmake ..
$ make -jN
$ make test
$ make install
Replace N
with the number of threads you want to use to build, e.g. -j8
.
You may want to specify additional flags to cmake
:
-
-DCMAKE_INSTALL_PREFIX=$HOME/.local
, or another custom install location. This is required if you don't have permission to write to the default install prefix path (/usr/local
). Mind that ALPSCore installs a CMake script that has to be picked up by CMake when building our codes. Thus, any non-standard install location needs to be matched by a-DCMAKE_PREFIX_PATH=<...>
flag when configuring the client codes. - If a local version of boost has been installed, you may need topoint CMake
to it by specifying
-DBOOST_ROOT=/path/to/boost/install
. Otherwise your local version may not be found, or be shadowed by an incompatible version.
Building client codes
The TK-SVM framework does not follow a traditional build model, in that it does
not provide a library for the client codes to link against. This is rooted in
the fact that it provides main
functions for the five executables mentioned in
the foregoing section. Instead, the root-level
CMakeLists.txt
exposes the relevant source files as CMake
variables ${TKSVM_SAMPLE_SRC}
(and analogous for the other executables), the
libraries as ${TKSVM_LIBRARIES}
and the include directories as
${TKSVM_INCLUDE_DIRS}
.
These may then be used from the CMakeLists.txt
file of each of the client
codes. It must also add the compile definitions exposed in
${TKSVM_DEFINITIONS}
and additionally define a variable TKSVM_SIMINCL
holding the name of the header of the main ALPSCore simulation class which
should define a type (alias) sim_base
to that class.
Refer to the build file of one of the client codes (e.g.
gauge/CMakeLists.txt
) for an example. Third-party
client codes may rather keep this repository as a submodule in their own tree,
in which case the line
add_subdirectory(.. tksvm)
simply becomes
add_subdirectory(svm-order-params)
Refer to the READMEs of the client codes for build instructions on them
specifically, e.g. Building the gauge
client
code.
Build configuration
Most of the behavior of the programs can be controlled through runtime parameters. However, some options have to be selected at compile time via the CMake variables listed below. To customize a variable, pass the appropriate flag on to CMake, e.g.
$ cmake -DCONFIG_MAPPING=LAZY ..
or use the interactive ccmake
configurator.
Variable name | Possible values |
---|---|
ALPSCore_DIR |
nonstandard/path/to/share/ALPSCore |
Eigen3_DIR |
nonstandard/path/to/share/eigen3/cmake |
CMAKE_BUILD_TYPE |
Release (default), Debug
|
CMAKE_INSTALL_PREFIX |
path to install directory (executables will be copied into bin/ ) |
SVM__ENABLE_TESTS |
OFF (default), ON
|
CONFIG_MAPPING |
EAGER (default), LAZY
|
Basic usage
This section demonstrates the basic usage of the framework, drawing from the
examples provided along with the gauge
client code as an example, though the
principles hold true for all applications.
Parameter files
Simulation parameters are stored in INI-style files.
We provide example parameter files for each client code, in their respective
params
directories (e.g. for the gauge code).
Each parameter file comes with a brief comment which explains its purpose and
includes the shell commands to run that particular simulation.
They also contain estimates for the run time of the individual steps.
In the following subsections, we will give a more detailed description of some aspects of the learning, testing, and analysis stages. For a (somewhat) comprehensive listing of supported parameters, refer to the section entitled → Runtime parameters.
Sampling: collecting Monte Carlo snapshots from across the parameter space
The *-sample
binaries run the Monte Carlo simulation at a number
of points in the phase diagram. The configurations are sampled periodically and
stored as the optimization problem in a file ending in .clone.h5
. For
example,
$ mpirun -n 4 gauge-sample Td-hyperplane.ini
will carry out Monte Carlo simulations of the gauge model with the tetrahedral Td symmetry group.
Note: Both the *-sample
and *-test
programs rely on MPI for
parallelization. Replace 4
with the number of MPI processes you want to run
in parallel.
Time limit: If the required amount of Monte Carlo steps cannot be carried
out within timelimit
seconds, the simulation terminates
prematurely and the incomplete data are written to Td-hyperplane.clone.h5
. One
can resume the simulation with
$ mpirun -n 4 gauge-sample Td-hyperplane.clone.h5
Sweep policy: The set of parameter points at which samples are collected is determined by a so-called sweep policy. Multiple, rather flexible, options exist which are documented and elaborated in the section entitled → Sweep through the phase diagram.
EAGER
vs. LAZY
mode
Dependent upon the build-time option CONFIG_MAPPING
, the *-sample
program
may or may not already carry out the tensorial feature mapping on the data. In
EAGER
mode, the mapping is carried out as soon as a snapshot is sampled and
only the tensorial features are stored. For large lattices and modest ranks,
this provides potentially huge savings in storage and memory. On the flipside,
the user commits themselves to a fixed choice of the tensorial kernel (TK) rank
and cluster choice.
In LAZY
mode, the mapping to monomial features is deferred until the learning
stage and the full spin configurations are stored in the *.clone.h5
file. This
may be desirable when the Monte Carlo simulation is particularly hard and
requires long intermediate updates to generate independent configurations. In
this scenario, large lattices and huge amounts of samples are typically
unreachable anyway, so storing the full spin configurations remains feasible.
The lazy mode then helps to avoid having to run the costly simulations
separately for each choice of the TK rank or spin cluster; rather, the
corresponding parameters can be overridden when launching the *-learn
program.
Note that the choice of the CONFIG_MAPPING
build option modifies the internal
representation used in the *.clone.h5
files, so files generated using one
option are incompatible to programs built using the other.
Learning: SVM optimization
The *-learn
binary takes the samples stored in one or multiple *.clone.h5
files and performs the actual SVM optimization, saving the result in a file
ending in *.out.h5
.
$ gauge-learn Td-hyperplane.clone.h5
Labeling samples: Within the *.clone.h5
files, the samples are stored
along with the full point in parameter space they have been sampled at. Before
the data can be fed to the SVM, those points have to be mapped to discrete
(one-dimensional) labels. This labeling is done by a component called the
→ classifier. Since many different strategies for labeling are
conceivable, corresponding to varying degrees of physical prior knowledge on the
system or knowledge gained through previous analyses, the classifier policy
can be selected and configured from a range of options at runtime. For example,
the parameter file
gauge/params/Td-hyperplane.ini
uses the
hyperplane
classifier, mapping the phase diagram points to a binary label,
indicating if samples are from the ordered or disordered phase, respectively.
For options on the choice and configuration of the different classifiers, refer
to the dedicated section.
Multiclassification: In contrast, when using the D2h symmetry group,
$ mpirun -n 4 gauge-sample D2h-phase-diagram.ini
$ gauge-learn D2h-phase-diagram.clone.h5
the phase diagram features three phases called O(3) (isotropic),
D∞h (uniaxial), and D2h (biaxial). The parameter file
gauge/params/D2h-phase-diagram.ini
is
intended to be used with the phase_diagram
classifier which maps the
(J1, J3) points to three labels, O(3)
, Dinfh
, and
D2h
, according to the known phase diagram. The SVM will then solve three
classification problems, yielding one decision function to distinguish between
each pair of the labels. Given that the labels indeed correspond to the physical
phases and the rank is chosen appropriately (as is the case for the
D2h example), these three decision functions can be related to the
order parameter(s).
Note that the same samples can also be used with a different classifier by
overriding the classifier.policy
parameter. For example, when selecting the
fixed_from_sweep
classifier, the samples will be labeled according to the grid
points in phase diagram space from which they have been sampled:
$ gauge-learn D2h-phase-diagram.clone.h5 \
--nu=0.1 \
--classifier.policy=fixed_from_sweep \
--outputfile=D2h-grid.out.h5
The phase diagram has previously been sampled on a 10 x 10 grid, hence, the
samples are now assigned 100 distinct labels, resulting in 100 * 99 / 2 = 4950
decision functions. The result is written to the alternative location
D2h-grid.out.h5
to avoid overwriting the previous result. It is obviously not
feasible to interpret each of the decision functions individually, but we can
use their bias parameters for the
→ graph analysis. To that end, a
moderate regularization level is appropriate; hence ν = 0.1 is chosen.
Merging additional clones: One may also perform the SVM optimization on a
collection of samples drawn from a variety of different runs of the *-sample
program with different, potentially overlapping sweeps through the parameter
space. This is possible; one can provide additional *.clone.h5
files in a
colon-separated string to the merge
parameter:
$ gauge-learn first.clone.h5 --merge=second.clone.h5:third.clone.h5
Note that first.clone.h5
is special in that its parameters are used
subsequently, specifically they determine the selected classifier used to label
all samples, including those contributed by second.clone.h5
and
third.clone.h5
.
Command line options
Long flag | Short | Description |
---|---|---|
--help |
-h |
Display ALPSCore help message (lists parameters) |
--merge=<clone-list> |
Specify a colon-separated list of additional *.clone.h5 files whose samples should be included in the analysis |
|
--infinite-temperature |
-i |
Include sweep.samples fictitious samples as obtained from Simulation::random_configuration() as a control group |
--statistics-only |
Collect all samples, label them by the classifer and print their statistics, but forego the actual SVM optimization |
Note that additionally runtime parameters may also be overridden using command line arguments.
Testing: measuring decision function and observables
$ mpirun -n 4 gauge-test Td-hyperplane.out.h5
will read the SVM model that was previously learned and measure its decision
function as an observable in independent Monte Carlo simulations. Additionally,
any observables that are ordinarily measured in the simulation's measure()
function will be measured as well. It is also possible to invoke the *-test
program on an *.ini
file, e.g.
$ mpirun -n 4 gauge-test Td-hyperplane.ini
which will measure only those ordinary observables and not attempt to evaluate any SVM model. This may be useful then debugging the simulation itself.
The set of parameter points at which measurements are to be taken is again
specified by a sweep policy in much the same way this is done for the
*-sample
program. Refer to the detailed explanation of the testing stage
parameters below. The example
Td-hyperplane.ini
sweeps a number of points
on a line segment between two end points in the parameter space.
Output format: The results are summarized in a text file
Td-hyperplane.test.txt
and full observables are stored in
Td-hyperplane.test.h5
. The *.test.txt
is self-documenting (generating
comments indicating the column layout at the beginning of the file) and suitable
for plotting using gnuplot
. A plot of the decision function can be created
using:
gnuplot> plot 'Td.test.txt' using 2:5:6 with yerrorlines
Variances: For any observable Obs
that is measured by the simulation's
measure()
function and for which there exists another observable named
Obs^2
, the *.test.txt
will additionally contain columns giving the variance
of Obs
, as calculated using the formula <Obs^2> - <Obs>^2
, including an
error estimate on that variance, as obtained from jackknife resampling. This
extends to the SVM decision function(s) whose squares are also always measured.
The physical utility of the variables of decision functions, particularly some
analogies to generalized susceptibilities, are discussed in
arXiv:1907.12322.
Note: Variances are also calculated for observables of higher powers or
absolute values, e.g. when measuring all of Obs
, |Obs|
, Obs^2
, Obs^4
,
and Obs^8
, the variances of Obs
, |Obs|
, Obs^2
, and Obs^4
are all
calculated.
Command line options
Long flag | Short | Description |
---|---|---|
--help |
-h |
Display ALPSCore help message (lists parameters) |
--rescale |
-r |
Shift decision functions by the bias and rescale to unit interval |
Note that additionally runtime parameters may also be overridden using command line arguments.
Extracting the coefficient matrix
$ gauge-coeffs Td-hyperplane.out.h5 -u
reads the SVM model and contracts over the support vectors to extract the coefficient matrix, block structure, and performs miscellaneous analyses and processing steps, such as fitting and removing self-contractions or comparing to the exact result.
The three panels in Fig. 2 of arXiV:1804.08557 can be reproduced by the following invocations:
$ gauge-coeffs Td-hyperplane.out.h5 -u
$ gauge-coeffs Td-hyperplane.out.h5 -u --block=[lmn:lmn]
$ gauge-coeffs Td-hyperplane.out.h5 -u --block=[lmn:lmn] --remove-self-contractions
These options are detailed below.
Index arrangement
-
-u | --unsymmetrize
: coefficients involving redundant monomial (cf. Supplementary Materials) are reconstructed, i.e. the coefficient value of the corresponding non-redundant coefficient is copied. -
-r | --raw
: the indicies are not rearranged. By default, the indices are reshuffled in the form (α1, ..., αn, a1, ..., an) and lexicographically ordered such that the color-block structure becomes apparent. Specifying this flag disables this, and indices are arranged as (α1, a1, ..., αn, an).
Extracting individual blocks
The program has two major modes: extraction of the full coefficient matrix and
single-block mode where a single block is targeted and extracted exclusively.
This is done by specifying the command line parameter --block=<block-spec>
where <block-spec>
identifies a block by its color indices. E.g. in the case
of the tetrahedral (rank-3) order, the non-trival block can be found by
specifying --block=[lmn:lmn]
.
Contraction analysis
-
-c | --contraction-weights
: for each block (or only the single block), perform a least-squares fit of contraction "masks" to obtain the coefficients with which each contraction contributes and output those. -
-s | --remove-self-contractions
: perform the same analysis as above, but consequently isolate the contributions due to self-contractions and subtract them from the full coefficient matrix.
Note that "blocks" refers to the representation where indices have been
reshuffled. When the full coefficient matrix is extracted, the contraction
analysis (-c
, -s
) operates on the symmetrized representation where redundant
element have not been reinstated and "blocks" are actually non-local. In
single-block mode however, this would not work because the equivalents of
redundant elements are actually part of different blocks. Thus, in single-block
mode, the contraction-analysis is performed on the "unsymmetrized" (redundant)
representation. Consequently, the contraction coefficients will be different
from those obtained for the same block in full-matrix mode, as now equivalent
contractions will contribute evenly, as opposed to only those contractions that
are compatible with the symmetrization.
Multiclassification
When analyzing the result of a multiclassification problem (see the D2h example above), the analysis is performed for each decision function separately. Because of the quadratic growth of the number of decision functions, this can get messy for big multiclassification problems. Instead, one can list the decision functions and their biases without extracting the coefficient matrices using:
$ gauge-coeffs D2h-hyperplane.out.h5 --list
0: O(3) -- Dinfh rho = -1.0052
1: O(3) -- D2h rho = -1.00266
2: Dinfh -- D2h rho = -0.659734
In this lists, the decision functions are given a running number. This number can be used to select a single decision function and extract only its coefficient matrix. For example, to limit the analysis to the O(3) / D2h transition:
$ gauge-coeffs D2h-hyperplane.out.h5 -u --transition=1
Command line options
Long flag | Short | Description |
---|---|---|
--help |
-h |
Display ALPSCore help message (lists parameters) |
--verbose |
-v |
Print information on currently running step |
--list |
-l |
List the transitions between labels and their bias value, but do not write any files |
--transition=<num> |
-t |
Rather than extracting coefficients for all transitions in a multiclassification setting, only consider the transition numbered <num> (as listed using the --list option) |
--unsymmetrize |
-u |
Assign the extracted coefficient to all equivalent monomials, rather than only to a single representative |
--raw |
-r |
Indicies are not rearranged; redundant monomials are not included in the output pattern |
--blocks-only |
-b |
Skip output of the full coefficient matrix and only output block structure; unavailable (and pointless) in single-block mode. |
--block=<block-spec> |
Only output the coefficient matrix block identified by <block-spec>
|
|
--contraction-weights |
-c |
For each block (or only the single block), perform a least-squares fit of contraction "masks" to obtain the contributions of each contraction and output those. |
--remove-self-contractions |
-s |
Perform the same analysis as above, but consequently isolate the contributions due to self-contractions and subtract them from the full coefficient matrix. |
--exact |
-e |
Calculate and output the exact result (as specified by --result ) for the coefficient matrix (if available) |
--diff |
-d |
Output the deviation from the exact result (as specified by --result ) for the coefficient matrix and print total relative deviation δ
|
--result=<name> |
The name of the exact result to look up; hardcoded values: Cinfv , Dinfh , D2h ; required when using --exact or --diff
|
Short flags may be combined into one multi-flag, as is customary for POSIX-compatible programs.
Spectral graph partitioning analysis
The *-segregate-phases
program is intended to be used in conjunction with
massive multiclassification programs (such as the previously constructed file
D2h-grid.out.h5
) to infer the topology of the phase diagram based on the bias
criterion.
The rationale is, that if the bias of a decision function between two points is far from one (the ideal value expected for a physical transitions), those two points are either in the same phase, or the phase transition could not be captured at the rank considered. We then go ahead and construct an undirected simple graph by adding an edge connecting those two points (vertices). On the other hand, if the bias is "close" to one, we assume there is a phase transition taking place between the points and we don't put an edge. The threshold value ρc can be tuned to control the inclusivity of the criterion: an edge is included in the graph if the deviation of corresponding unoriented bias |ρ| from unity exceeds ρc: ||ρ| - 1| > ρc.
The resulting graph can be subjected to a spectral partitioning analysis. To that end, the eigendecomposition of the Laplacian of the graph is calculated. The degeneracy of the lowest eigenvalue (zero) indicates the number of connected components within the graph. We recommend to choose ρc such that the graph consists of a single connected component, i.e. a slightly larger value of ρc would lead to a degeneracy of eigenvalue zero. For the above example of the D2h symmetry, this sweet spot happens to be achieved at around ρc=0.9, but your mileage may vary:
$ gauge-segregate-phases D2h-grid.out.h5 --rhoc=1.0
Degeneracy of smallest eval: 2
$ gauge-segregate-phases D2h-grid.out.h5 --rhoc=0.9
Degeneracy of smallest eval: 1
Then, the second largest eigenvalue is the so-called algebraic connectivity and its corresponding eigenvector, the Fiedler vector can be used to infer the phase diagram: if the graph features regions which are strongly intraconnected, but weakly interconnected, the algebraic connectivity will be small and the entries of the Fiedler vector will cluster such that the entries corresponding to intraconnected regions will have a similar value.
Thus, the Fiedler vector can be visualized and will be reminiscient of the phase
diagram. To that end, the *-segregate-phases
program outputs two files,
edges.txt
and phases.txt
. The former lists pairs of points corresponding to
the edges of the graph, the latter gives all the eigenvectors of the Laplacian
matrix such that the second dataset (index 1) is the Fiedler vector. Using
gnuplot, we can see that we achieve a decent approximation of the
D2h phase diagram:
gnuplot> plot 'edges.txt' using 2:1 with lines, \
'phases.txt' index 1 using 2:1:3 with points pt 7 lc palette
Weighted graphs
Rather than making a binary decision on whether or not to include an edge into
the graph, one can in some situations obtain better results by constructing a
weighted graph whose edges have weights which are calculated from their biases
by means of a weighting function w(ρ). Three options are provided which may be
selected using the --weight
flag:
<weight-func-name> |
Equation w(ρ) | Default value of ρc |
---|---|---|
box (default) |
![]() |
1 |
gaussian |
![]() |
standard deviation of the distribution of ρ values |
lorentzian |
![]() |
half interquartile spacing of the distribution of ρ values |
In each of these, the parameter ρc constitutes a characteristic
scale. The default box
option corresponds to the strategy using unweighted
graphs outlined above.
Masks
In phase diagrams where multiple different transitions occur at different ranks,
it may be desirable to exclude a region of the parameter space which has been
identified to exhibit order at a certain rank from the subsequent analysis of
the remaining parameter space at higher ranks. An example of this can be found
in the easy-plane antiferromagnetic (AF) phase in the pyrochlore XXZ model
studied in arXiv:1907.12322. To this end, the *-segregate-phases
program
allows for the creation and application of "masks". Sticking with the example of
the AF phase, one first performs a graph analysis of the whole parameter space
at rank 1. The resulting Fiedler vector clearly reflects the bipartition of the
parameter space into the AF phase and "the rest". By specifying the option
--threshold=0
, a file mask.txt
will be written which identifies all those
parameter points whose associated Fiedler vector entry is larger than 0
; in
the example, this would correspond to the AF phase. Since we want to further
analyze the remainder of the parameter space, the flag --invert-mask
is
provided. Now mask.txt
identifies all those points which are not in the AF
phase.
Subsequently, we perform a graph analysis at rank 2, passing the flag
--mask=mask.txt
on the command line. This way, the previously created mask is
applied and the parameter points belonging to the AF phase are ignored in the
construction of the graph. They are also omitted from the output files
phases.txt
and edges.txt
. Particularly when visualizing the Fiedler vector
(phases.txt
) using e.g. gnuplot
this can be a problem as each row of the
"image" may have a different number of pixels due to the omitted points. As a
remedy, the additional flag --masked-value=42
will cause these ignored points
to be included yet again in the output but with the value 42
in lieu of their
Fiedler vector element. When rendering the result, the color palette may be
chosen such that pixels with the value 42
are shown in gray, indicating their
absence from the analysis.
Command line options
Long flag | Short | Description |
---|---|---|
--help |
-h |
Display ALPSCore help message (lists parameters) |
--verbose |
-v |
Print information on currently running step |
--weight=<weight-func-name> |
-w |
Specify the weighting function used to mapping biases to edge weights; <weight-func-name must be one of: box (default), gaussian , lorentzian (see below) |
--rhoc=<number> |
-r |
Specify the characteristics bias scale ρc used to rescale the weighting function; default value depends on choice of weighting function (see below) |
--radius=<max-distance> |
-R |
Specify a cutoff radius: edges between parameters points farther than <max-distance> apart from oneanother will not be included in the graph; defaults to infinity |
--threshold=<tval> |
-t |
Trigger the output of a mask.txt file, identifying all parameter points with Fiedler vector elements larger or equal to <tval>
|
--invert-mask |
Inverts the behavior of the --threshold option, i.e. identifies points with Fiedler vector elements less than <tval>
|
|
--mask=<mask-filename> |
-m |
Applies a previously saved mask, i.e. ignores all parameter points not included in the mask from the graph analysis |
--masked-value=<mval> |
When outputting the results of the graph analysis, specifically phases.txt , replace the would-be Fiedler vector element of those parameter points which have been ignored due to the --mask with <mval> ; otherwise omit the point entirely |
Runtime parameters
This section lists all parameters that can be used to control the behavior of the simulations with respect to runtime duration and output, regularization of the SVM optimization, the way the parameter space is swept, and the way samples are labeled (classified) before being fed to the SVM. These parameters are shared by all client codes.
Client codes may define parameters of their own, conventionally pertaining to
the physical parameters of the simulated model (lattice size, coupling, etc.),
the MC update scheme (number of sweeps, update attempt probabilities, etc.),
and the choice of the tensorial kernel that is to be used (rank, cluster choice,
etc.).
For the three client codes shipped as part of this repository, those parameters
are documented in the READMEs of each of the corresponding subdirectories.
Note that some of these parameters appear in each of these cases (e.g. rank
,
total_sweeps
); however no assumption is made on part of the shared framework
that any of them need to be present.
Most parameters feature a hierarchical name where levels in the hierarchy are
separated by dots. These parameters can also be grouped by hierarchical prefix.
The following three specification of the parameter sweep.grid.N1
are
equivalent:
sweep.grid.N1 = 42
[sweep]
grid.N1 = 42
[sweep.grid]
N1 = 42
Note that the hierarchical prefix has to be specified in full and cannot be nested. I.e., the following version does not work: