API reference¶
Note
Work in progress. Many function do not yet follow the numpy standard in their documentation!
This page contains a structured overview of the python API of adcc. See also the Index.
adcc module¶
Functions¶
|
Run an ADC calculation. |
|
Evaluate Einstein summation convention for the operands similar to numpy’s einsum function. |
|
Form a single, einsum-like contraction, that is contract tensor a and be to form out via a contraction defined by the first argument string, e.g. |
|
Return a copy of the input tensor. |
|
Form the scalar product between two tensors. |
|
Return an empty tensor of the same shape and symmetry as the input tensor. |
|
Force full evaluation of a tensor expression |
|
Form a linear combination from a list of tensors. |
|
Return tensor of the same shape, but without the symmetry setup of the input tensor. |
|
Return tensor of the same shape and symmetry as the input tensor, but initialised to 1, that is the canonical blocks are 1 and the other ones are symmetry-equivalent (-1 or 0) |
|
Return the transpose of a tensor as a copy. |
|
|
|
Return a zero tensor of the same shape and symmetry as the input tensor. |
|
|
|
Set the number of running worker threads used by adcc |
|
Get the number of running worker threads used by adcc. |
|
Obtain guesses for computing singlet states by inspecting the passed ADC matrix. |
|
Obtain guesses for computing triplet states by inspecting the passed ADC matrix. |
|
Obtain guesses by inspecting a block of the diagonal of the passed ADC matrix. |
|
Return guess symmetry objects (one for each AmplitudeVector block) such that the specified requirements on the guesses are satisfied. |
|
Obtain guesses for computing spin-flip states by inspecting the passed ADC matrix. |
|
Return an AmplitudeVector object filled with zeros, but where the symmetry has been properly set up to meet the specified requirements on the guesses. |
|
Run an ADC calculation. |
|
Run an ADC calculation. |
|
Run an ADC calculation. |
|
Run an ADC calculation. |
|
Run an ADC calculation. |
|
Run an ADC calculation. |
|
Run an ADC calculation. |
|
Run an ADC calculation. |
|
Run an ADC calculation. |
|
Run an ADC calculation. |
|
Run an ADC calculation. |
|
Return a nice banner describing adcc and its components |
Classes¶
Exception thrown during the validation stage of the arguments passed to :py:`run_adc` to signal that an input is not valid. |
|
|
|
|
|
|
|
|
|
|
|
|
|
Abstract class defining the interface for passing data from the host program to adcc. |
|
|
|
|
|
|
|
|
|
|
|
|
|
|
libadcc: Python bindings¶
The libadcc Python module contains python bindings of libadcc: C++ library. They are generated directly from the C++ source code using pybind11 and the sources contained in the extension subfolder of the adcc github repository.
It is not recommended to call these functions directly, but instead resort to the higher-level functionality from the adcc module.
Functions¶
|
Apply the spin symmetrisation required to make the doubles and higher parts of an amplitude vector consist of components for a particular spin kind only. |
|
|
|
|
|
Fill the passed vector of doubles blocks with doubles guesses using the delta-Fock matrices df02 and df13, which are the two delta-Fock matrices involved in the doubles block. |
Get the number of running worker threads used by adcc. |
|
Get the total number of threads (running and sleeping) used by adcc. |
|
|
|
|
Return the Symmetry object like it would be set up for the passed subspace of the electron-repulsion tensor. |
|
Return the Symmetry object for an orbital subspace block of a one-particle operator |
|
Return the symmetry object for an operator in the AO basis. The object will represent a block-diagonal matrix of the form ( M 0 ) ( 0 M ). where M is an n_bas x n_bas block and is indentical in upper-left and lower-right. |
|
Return the Symmetry object like it would be set up for the passed subspace of the orbital coefficients tensor. |
|
Return the Symmetry object like it would be set up for the passed subspace of the orbital energies tensor. |
|
Set the number of running worker threads used by adcc |
|
Set the total number of threads (running and sleeping) used by adcc. |
|
Overloaded function. |
|
Overloaded function. |
Classes¶
Class controlling the memory allocations for adcc ADC calculations. |
|
Abstract class defining the interface for passing data from the host program to adcc. |
|
Interface class representing the data expected in adcc from an interfacing HF / SCF program. |
|
Helper object to extract information from indices into orbitals subspaces and to map them between different indexing conventions (full MO space, MO subspaces, indexing convention in the HF Provider / SCF host program, . |
|
Class setting up the molecular orbital index spaces and subspaces and exposing information about them. |
|
Class representing a one-particle operator. |
|
Class representing information about the reference state for adcc. |
|
Container for Tensor symmetry information |
|
Class representing the Tensor objects used for computations in adcc |
libadcc: C++ library¶
A reference of the C++ part of libadcc and its classes and functions can be found in the following. The functions and classes discussed here are exposed to Python as the libadcc: Python bindings python module.
Reference state¶
This category lists the libadcc functionality,
which imports the data from the libadcc::HartreeFockSolution_i
interface into the libadcc::ReferenceState
for internal use by the library.
See Connecting host programs to adcc for details how to connect
host programs to adcc.
Important classes in the process are libadcc::MoSpaces
,
which collects information about the occupied and virtual
orbital spaces, and libadcc::MoIndexTranslation
,
which maps orbitals indices between the ordering used by libadcc
and the one used by the SCF program.
-
std::vector<size_t>
construct_blocks
(std::vector<size_t> block_starts_crude, size_t length, size_t max_block_size)¶ Construct a list of Tensor blocks (by identifying their starting indices) from a list of indices that identify places where a block definitely has to start and a maximal block size.
The function checks the block sizes resulting from block_starts_crude and if they are beyond max_block_size tries to split them evenly.
- Parameters
block_starts_crude
: list of crude block startslength
: Total number of indicesmax_block_size
: Maximal size of a block
-
struct
libadcc
::
AxisRange
: public std::pair<size_t, size_t>¶ - #include <MoIndexTranslation.hh>
Very simple structure holding the range of indices for a single axis
- Note
This structure is not intended for wide use throughout the adcc code, much rather it should help to make the MoIndexTranslation code more clear. Outside this class it can be well thought of as a std::pair<size_t, size_t>.
-
struct
libadcc
::
SimpleRange
: public std::vector<std::pair<size_t, size_t>>¶ - #include <MoIndexTranslation.hh>
Very simple structure holding an index range.
The ranges are given as lists of ranges in each axis, where for each axis the first number is the first index of the range to be taken and the second number is the first index past the range, i.e. the index range in each axis are given as half-open intervals [start, end).
- Note
This structure is not intended for wide use throughout the adcc code, much rather it should help to make the MoIndexTranslation code more clear. Outside this class it can be well thought of as a std::vector<std::pair<size_t, size_t>>.
Public Functions
-
SimpleRange
(std::vector<size_t> starts, std::vector<size_t> ends)¶ Build a Range object by specifying all starts and all ends
-
SimpleRange
(std::vector<std::pair<size_t, size_t>> in)¶ Build a Range object by providing a vector of [start, end) pairs
-
SimpleRange
()¶ Build an empty SimpleRange object
-
size_t
ndim
() const¶ Return the dimensionality, i.e. the size of the vector
-
void
push_axis
(std::pair<size_t, size_t> new_axis)¶ Add an axis by specifying its range
-
std::vector<size_t>
starts
() const¶ Get the vector of all range starts
-
std::vector<size_t>
lasts
() const¶ Get the vectors of the last indices in all axes
-
std::vector<size_t>
ends
() const¶ Get the vector of all range ends
-
struct
libadcc
::
RangeMapping
: public std::pair<SimpleRange, SimpleRange>¶ - #include <MoIndexTranslation.hh>
Datastructure to describe a mapping of one range of indices to another range of indices.
- Note
This structure is not intended for wide use throughout the adcc code, much rather it should help to make the MoIndexTranslation code more clear. Outside this class it can be well thought of as a pair of above range objects, which are itself effectively std::vector<std::pair<size_t, size_t>>.
Public Functions
-
SimpleRange &
from
()¶ The range from which we map
-
SimpleRange &
to
()¶ The range to which we map
-
class
libadcc
::
MoIndexTranslation
¶ - #include <MoIndexTranslation.hh>
Translation object, which helps to translate various representations of a tensorial index over orbitals subspaces. E.g. it helps to identify the spin block, reduce indices to spatial indices, translate indices or index ranges to the original ordering used in the SCF program and similar tasks.
Public Functions
Create an MoIndexTranslation object for the provided space string.
Create an MoIndexTranslation object for a list of subspace identifiers
-
std::shared_ptr<const MoSpaces>
mospaces_ptr
() const¶ Return the MoSpaces to which this MoIndexTranslation is set up
-
std::string
space
() const¶ The space used to initialise the object
-
const std::vector<std::string> &
subspaces
() const¶ The space splitup into subspaces along each dimension
-
size_t
ndim
() const¶ Number of dimensions
-
const std::vector<size_t> &
shape
() const¶ Shape of each dimension
-
std::vector<size_t>
full_index_of
(const std::vector<size_t> &index) const¶ Map an index given in the space passed upon construction to the corresponding index in the range of all MO indices (the ffff space)
-
std::vector<size_t>
block_index_of
(const std::vector<size_t> &index) const¶ Get the block index of an index
-
std::vector<size_t>
block_index_spatial_of
(const std::vector<size_t> &index) const¶ Get the spatial block index of an index
The spatial block index is the result of block_index_of modulo the spin blocks, i.e. it maps an index onto the index of the spatial blocks only, such that the resulting value is idential for two index where the MOs only differ by spin. For example the 1st core alpha and the 1st core beta orbital will map to the same value upon a call of this function.
-
std::vector<size_t>
inblock_index_of
(const std::vector<size_t> &index) const¶ Get the in-block index, i.e. the index within the tensor block
-
std::string
spin_of
(const std::vector<size_t> &index) const¶ Get the spin block of each of the index components
-
std::pair<std::vector<size_t>, std::vector<size_t>>
split
(const std::vector<size_t> &index) const¶ Split an index into block index and in-block index and return the result
-
std::tuple<std::string, std::vector<size_t>, std::vector<size_t>>
split_spin
(const std::vector<size_t> &index) const¶ Split an index into a spin block descriptor, a spatial block index and an in-block index.
-
std::vector<size_t>
combine
(const std::vector<size_t> &block_index, const std::vector<size_t> &inblock_index) const¶ Combine a block index and an in-block index into a subspace index
-
std::vector<size_t>
combine
(const std::string &spin_block, const std::vector<size_t> &block_index_spatial, const std::vector<size_t> &inblock_index) const¶ Combine a descriptor for the spin block, a spatial-only block index and an in-block index into a subspace index.
-
std::vector<size_t>
hf_provider_index_of
(const std::vector<size_t> &index) const¶ Map an index given in the space passed upon construction to the corresponding index in the HF provider, which was used to start off the adcc calculation
-
std::vector<RangeMapping>
map_range_to_hf_provider
(const SimpleRange &range) const¶ Map a range of indices to host program indices, i.e. the indexing convention used in the HfProvider / HartreeFockSolution_i classes.
Since the mapping between subspace and host program indices might not be contiguous, a list of range pairs is returned. In each pair the first entry represents a range of indices (indexed in the MO subspace) and the second entry represents the equivalent range of indices in the Hartree-Fock provider these are mapped to.
-
struct
libadcc
::
MoSpaces
¶ - #include <MoSpaces.hh>
Information about the molecular orbitals spaces and subspaces, i.e. the occupied versus virtual orbitals, active versus frozen orbitals, core-occupied versus valence-occupied orbitals.
Public Functions
-
size_t
n_orbs
(const std::string &space = "f") const¶ The number of orbitals inside a particular orbital space
-
size_t
n_orbs_alpha
(const std::string &space = "f") const¶ The number of alpha orbitals inside a particular orbital space
-
size_t
n_orbs_beta
(const std::string &space = "f") const¶ The number of beta orbitals inside a particular orbital space
-
std::string
subspace_name
(const std::string &space = "f") const¶ Obtain a pretty name for the provided subspace
-
bool
has_core_occupied_space
() const¶ Does this object have a core-occupied space (i.e. is ready for core-valence separation)
-
bool
has_subspace
(const std::string &space) const¶ Does this object have a particular subspace
-
std::string
irrep_totsym
() const¶ Return the totally symmetric irreducible representation.
Construct an MoSpaces object from a HartreeFockSolution_i, a pointer to an AdcMemory object.
- Parameters
hf
: HartreeFockSolution_i object to useadcmem_ptr
: ADC memory keep-alive object to be used in all Tensors constructed using this MoSpaces object.core_orbitals
: List of orbitals indices (in the full fock space, original ordering of the hf object), which defines the orbitals to be put into the core space, if any. The same number of alpha and beta orbitals should be selected. These will be forcibly occupied.frozen_core_orbitals
: List of orbital indices, which define the frozen core, i.e. those occupied orbitals, which do not take part in the ADC calculation. The same number of alpha and beta orbitals has to be selected.frozen_virtuals
: List of orbital indices, which the frozen virtuals, i.e. those virtual orbitals, which do not take part in the ADC calculation. The same number of alpha and beta orbitals has to be selected.
-
std::shared_ptr<const AdcMemory>
adcmem_ptr
() const¶ Return a pointer to the memory keep-alive object
-
size_t
libtensor_irrep_index
(std::string irrep) const¶ Return the libtensor index used for the irrep string
- Note
this function is for internal usage only.
Public Members
-
std::string
point_group
¶ The name of the point group for which the data in this class has been set up.
-
std::vector<std::string>
irreps
¶ The list of all irreducible representations in this point group. The first entry will always be the totally symmetric irreducible representation.
-
bool
restricted
¶ Are the orbitals resulting from a restricted calculation, such that alpha and beta electron share the same spatial part
-
std::vector<std::string>
subspaces_occupied
¶ The list of occupied subspaces known to this object
-
std::vector<std::string>
subspaces_virtual
¶ List of virtual subspaces known to this object
-
std::vector<std::string>
subspaces
¶ List of all subspaces known to this object
-
std::map<std::string, std::vector<size_t>>
map_index_hf_provider
¶ Contains for each orbital space (e.g. f, o1) a mapping from the indices used inside the respective space to the molecular orbital index convention used in the host program, i.e. to the ordering in which the molecular orbitals have been exposed in the HartreeFockSolution_i or HfProvider object.
-
std::map<std::string, std::vector<size_t>>
map_block_start
¶ Contains for each orbital space the indices at which a new block starts. Thus this list contains at least on index, namely 0.
-
std::map<std::string, std::vector<std::string>>
map_block_irrep
¶ Contains for each orbital space the mapping from each block used inside the space to the irreducible representation it correspond to.
-
std::map<std::string, std::vector<char>>
map_block_spin
¶ Contains for each orbital space the mapping from each block used inside the space to the spin it correspond to (‘a’ is alpha and ‘b’ is beta)
-
size_t
-
class
libadcc
::
ReferenceState
¶ - #include <ReferenceState.hh>
Struct containing data about the reference state on top of which the ADC calculation is performed
Public Functions
-
bool
restricted
() const¶ Return whether the reference is restricted or not
-
size_t
spin_multiplicity
() const¶ Return the spin multiplicity of the reference state. 0 indicates that the spin cannot be determined or is not integer (e.g. UHF)
-
bool
has_core_occupied_space
() const¶ Flag to indicate whether this is a cvs reference state
-
std::string
irreducible_representation
() const¶ Ground state irreducible representation
-
std::shared_ptr<const MoSpaces>
mospaces_ptr
() const¶ Return the MoSpaces to which this ReferenceState is set up
-
size_t
n_orbs
() const¶ Return the number of orbitals
-
size_t
n_orbs_alpha
() const¶ Return the number of alpha orbitals
-
size_t
n_orbs_beta
() const¶ Return the number of beta orbitals
-
size_t
n_alpha
() const¶ Return the number of alpha electrons
-
size_t
n_beta
() const¶ Return the number of beta electrons
-
std::vector<scalar_type>
nuclear_multipole
(size_t order) const¶ Return the nuclear contribution to the cartesian multipole moment (in standard ordering, i.e. xx, xy, xz, yy, yz, zz) of the given order.
-
double
conv_tol
() const¶ Return the SCF convergence tolerance
-
double
energy_scf
() const¶ Final total SCF energy
-
std::string
backend
() const¶ String identifying the backend used for the computation
-
std::shared_ptr<Tensor>
orbital_energies
(const std::string &space) const¶ Return the orbital energies corresponding to the provided space
-
std::shared_ptr<Tensor>
orbital_coefficients
(const std::string &space) const¶ Return the molecular orbital coefficients corresponding to the provided space (alpha and beta coefficients are returned)
-
std::shared_ptr<Tensor>
orbital_coefficients_alpha
(const std::string &space) const¶ Return the alpha molecular orbital coefficients corresponding to the provided space
-
std::shared_ptr<Tensor>
orbital_coefficients_beta
(const std::string &space) const¶ Return the beta molecular orbital coefficients corresponding to the provided space
-
std::shared_ptr<Tensor>
fock
(const std::string &space) const¶ Return the fock object of the context corresponding to the provided space.
-
std::shared_ptr<Tensor>
eri
(const std::string &space) const¶ Return the ERI (electron-repulsion integrals) contained in the context corresponding to the provided space.
-
void
import_all
() const¶ Normally the class only imports the Fock matrix blocks and electron-repulsion integrals of a particular space combination when this is requested by a call to above fock() or eri() functions.
This function call, however, instructs the class to immediately import all ERI tensors.
- Note
This typically does not need to be called explicitly.
-
std::vector<std::string>
cached_fock_blocks
() const¶ Return the list of momentarily cached fock matrix blocks
-
void
set_cached_fock_blocks
(std::vector<std::string> newlist)¶ Set the list of momentarily cached fock matrix blocks
- Note
This function is both capable to drop blocks (to save memory) and to request extra caching.
-
std::vector<std::string>
cached_eri_blocks
() const¶ Return the list of momentarily cached eri tensor blocks
-
void
set_cached_eri_blocks
(std::vector<std::string> newlist)¶ Set the list of momentarily cached ERI tensor blocks
- Note
This function is both capable to drop blocks (to save memory) and to request extra caching.
-
void
flush_hf_cache
() const¶ Tell the contained HartreeFockSolution_i object, that a larger amount of ERI data has just been computed and that the next request to further ERI tensor objects will potentially take some time, such that intermediate caches can now be flushed to save memory or other resources.
-
bool
ADC guess setup¶
Apply the spin projection required to make the doubles part of an amplitude vector consist of components for a singlet state only.
- Note
This function assumes a restricted reference with a closed-shell ground state (e.g. RHF).
- Parameters
tensor
: The tensor to apply the symmetrisation toblock
: The block of an amplitude this tensor representsspin_kind
: The kind of spin to enforce
One-particle operators¶
Warning
doxygengroup: Cannot find group “Properties” in doxygen xml output for project “libadcc” from directory: /home/mfh/Dokumente/Forschung_Wissenschaft/ADC/adcc/docs/../build/libadcc_docs/xml
Tensor interface¶
The generalised libadcc::Tensor
interface
used by adcc and libadcc to perform tensor operations.
-
TensorBackend
tensor_backend
()¶ Get some info about libtensor
Return the Symmetry object for the orbital energies tensor for the passed orbital subspace.
Return the Symmetry object for the orbital coefficient tensor for the passed orbital subspace.
- Parameters
mospaces_ptr
: MoSpaces pointerspace
: Space to use (e.g. fb or o1b)n_bas
: Number of AO basis functionsblocks
: Which blocks of the coefficients to return. Valid values are “ab” to return a tensor for both alpha and beta block as a block-diagonal tensor, “a” to only return a tensor for only alpha block and “b” for only the beta block. The fourth option is “abstack”, which returns the symmetry for a alpha block stacked on top of a beta block (This is the convention used by adcman).
Return the Symmetry object for the (antisymmetrised) electron-repulsion integral tensor in the physicists’ indexing convention for the passed orbital subspace.
Return the Symmetry object for the (non-antisymmetrised) electron-repulsion integral tensor in the physicists’ indexing convention for the passed orbital subspace.
- Note
This is not the symmetry of the eri objects in adcc. This is only needed for some special routines importing the ERI tensor into adcc. Use it only if you know you need it.
Return the Symmetry object for a MO spaces block of a one-particle operator
Note: This can be used for the Fock matrix with operator_class = “1”
- Parameters
mospaces_ptr
: MoSpaces pointerspace
: Space to use (e.g. o1o1 or o1v1)symmetric
: Is the tensor symmetric (only in effect if both space axes identical). false disables a setup of permutational symmetry.cartesian_transformation
: The cartesian function accordung to which the operator transforms. Valid values are: “1” Totally symmetric (default) “x”, “y”, “z” Coordinate axis “xx”, “xy”, “yz” … Products of two coordinate axis “Rx”, “Ry”, “Rz” Rotations about the coordinate axis
Return the symmetry object for an operator in the AO basis. The object will represent a block-diagonal matrix of the form ( M 0 ) ( 0 M ). where M is an n_bas x n_bas block and is indentical in upper-left and lower-right.
- Parameters
mospaces_ptr
: MoSpaces pointern_bas
: Number of AO basis functionssymmetric
: Is the tensor symmetric (only in effect if both space axes identical). false disables a setup of permutational symmetry.
Construct an uninitialised tensor using a Symmetry object
Construct an uninitialised tensor using no symmetry at all
Construct a tensor initialised to zero using a Symmetry object
-
struct
TensorBackend
¶ - #include <backend.hh>
Structure to hold information about the tensor backends
-
class
libadcc
::
Symmetry
¶ - #include <Symmetry.hh>
Container for Tensor symmetry
Public Functions
Create an empty symmetry object for a Tensor with the provided space string
Create an empty symmetry object for a Tensor with the provided space string, optionally using a map to supply the number of orbitals for some additional axes.
For the additional axis the pair contains either two numbers (for the number of alpha and beta orbitals in this axis) or only one number and a zero (for an axis, which as only one spin kind, alpha or beta).
-
const std::vector<AxisInfo>
axes
() const¶ Return the info about each of the contained tensor axes in order
-
std::string
space
() const¶ The space used to initialise the object
-
const std::vector<std::string> &
subspaces
() const¶ The space splitup into subspaces along each dimension
-
size_t
ndim
() const¶ Number of dimensions
-
std::vector<size_t>
shape
() const¶ Shape of each dimension
-
void
clear
()¶ Clear all symmetry elements stored in the datastructure
-
bool
empty
() const¶ Is the datastructure empty, i.e. are there no symmetry elements stored in here
-
void
set_irreps_allowed
(std::vector<std::string> irreps)¶ Set the list of irreducible representations, for which the symmetry shall be non-zero. If this is not set all irreps will be allowed
-
std::vector<std::string>
irreps_allowed
() const¶ Get the list of allowed irreducible representations. If this empty, all irreps are allowed.
-
bool
has_irreps_allowed
() const¶ Is a restriction to irreps set
-
void
clear_irreps_allowed
()¶ Clear the restriction to particular irreps, i.e. allow all irreps
-
void
set_permutations
(std::vector<std::string> permutations)¶ Set a list of index permutations, which do not change the tensor. A minus may be used to indicate anti-symmetric permutations with respect to the first (reference) permutation.
For example the vector {“ij”, “ji”} defines a symmetric matrix and {“ijkl”, “-jikl”, “-ijlk”, “klij”} the symmetry of the ERI tensor. Not all permutations need to be given to fully describe the symmetry. The check for errors and conflicts is only rudimentary at the moment, however.
-
std::vector<std::string>
permutations
() const¶ Return the list of equivalent permutations
-
bool
has_permutations
() const¶ Are there equivalent permutations set up
-
void
clear_permutations
()¶ Clear the equivalent permutations
-
std::vector<std::pair<std::vector<size_t>, scalar_type>>
permutations_parsed
() const¶ Return the list of equivalent permutations in a parsed internal format.
- Note
Internal function. May go or change at any time
-
void
set_spin_block_maps
(std::vector<std::tuple<std::string, std::string, double>> spin_maps)¶ Set lists of tuples of the form {“aaaa”, “bbbb”, -1.0}, i.e. two spin blocks followed by a factor. This maps the second onto the first with a factor of -1.0 between them.
-
std::vector<std::tuple<std::string, std::string, double>>
spin_block_maps
() const¶ Return the list of equivalent spin blocks
-
bool
has_spin_block_maps
() const¶ Are there equivalent spin blocks set up?
-
void
clear_spin_block_maps
()¶ CLear the equivalent spin blocks
-
void
set_spin_blocks_forbidden
(std::vector<std::string> forbidden)¶ Mark spin-blocks as forbidden (i.e. enforce them to stay zero). Blocks are give as a list in the letters ‘a’ and ‘b’, e.g. {“aaaa”, “abba”}
-
std::vector<std::string>
spin_blocks_forbidden
() const¶ Return the list of currently forbidden spin blocks
-
bool
has_spin_blocks_forbidden
() const¶ Are there forbidden spin blocks set up?
-
void
clear_spin_blocks_forbidden
()¶ Clear the forbidden spin blocks, i.e. allow all blocks
-
class
libadcc
::
Tensor
¶ - #include <Tensor.hh>
The tensor interface used by libadcc
Subclassed by libadcc::TensorImpl< N >
Interface of a tensor, which is exposed to python.
-
size_t
ndim
() const¶ Number of dimensions
-
const std::vector<size_t> &
shape
() const¶ Shape of each dimension
-
size_t
size
() const¶ Number of elements
-
std::vector<std::string>
subspaces
() const¶ Subspace labels identifying each axis
-
std::string
space
() const¶ Space to which the tensor is initialised
-
const std::vector<std::string> &
flags
() const¶ Set some flags for playing around
-
std::shared_ptr<Tensor>
empty_like
() const = 0¶ Return a new tensor with the same dimensionality and symmetry as the passed tensor. The elements have undefined values.
-
std::shared_ptr<Tensor>
nosym_like
() const = 0¶ Return a new tensor with same dimensionality and shape (and bispace) but without any symmetry copied over
-
std::shared_ptr<Tensor>
zeros_like
() const¶ Return a new tensor with the same dimensionality and symmetry as the passed tensor and all elements set to zero.
-
std::shared_ptr<Tensor>
ones_like
() const¶ Return a new tensor with the same dimensionality and symmetry as the passed tensor and all elements set to one.
-
void
set_mask
(std::string mask, scalar_type value) = 0¶ Set a mask into a tensor to a specific value. The mask to set is defined by the mask string. Repetitive indices define the values to be set. E.g. for a 6D tensor the mask string “iijkli” would set all elements T_{iijkli} for all ijkl to the given value.
-
void
fill
(scalar_type value)¶ Shift the elements matched by a mask into a tensor by a specific value. The mask to set is defined by the mask string. Repetitive indices define the values to be shifted. E.g. for a 6D tensor the mask string “iijkli” would shift all elements T_{iijkli} for all ijkl by the given value. Fill the tensor with a particular value, adhereing to symmetry.
-
std::shared_ptr<Tensor>
diagonal
(std::vector<size_t> axes) = 0¶ Extract a generalised diagonal from this Tensor.
axes
defines the axes over which to extract. The resulting new axis is appended to the remaining axis. For example if we have a Tensor fourth ordert
then diagonal(1, 2) returnsd_{iax} = t_{ixxa}
and diagonal(0, 1, 2) returnsd_{ax} = t_{xxxa}
.
-
void
set_random
() = 0¶ Set the tensor to random data preserving symmetry
-
std::shared_ptr<Tensor>
scale
(scalar_type c) const = 0¶ Scale tensor by a scalar value and return the result
Add another tensor and return the result
Add a linear combination to the represented tensor and store the result. Note that unlike
add
andscale
this operation is not lazy and is provided for optimising on cases where a lazy linear combination introduces a too large overhead.
Multiply another tensor elementwise and return the result
Divide another tensor elementwise and return the result
-
std::shared_ptr<Tensor>
transpose
(std::vector<size_t> axes) const = 0¶ Return a transpose form of this tensor
- Parameters
axes
: New permutation of the dimensions of this tensor (e.g. (1,0,2,3) will permute dimension 0 and 1 and keep the others) If missing, just reverse the dimension order.
Return the tensordot of this tensor with another
- Parameters
axes
: Axes to contract over. The first vector refers to the axes of this tensor, the second vector to the axes of the other tensor.
Compute the Frobenius or l2 inner product with a list of other tensors, returning the results as a scalar array
-
double
trace
(std::string contraction) const = 0¶ Return the full trace of a tensor
- Parameters
contraction
: Indices to contract over using Einstein summation convention, e.g. “abab” traces over the 1st and 3rd and 2nd and 4th axis.
Compute the direct sum of this tensor with another tensor. All axes are used in the order of appearance (this and then other).
-
std::shared_ptr<Tensor>
symmetrise
(const std::vector<std::vector<size_t>> &permutations) const = 0¶ Symmetrise with respect to the given index permutations by adding the elements resulting from appropriate index permutations.
Examples for permutations. Take a rank-4 tensor T_{ijkl} {{0,1}} Permute the first two indices, i.e. form 0.5 * (T_{ijkl} + T_{jikl}) {{0,1}, {2,3}} Form 0.5 * (T_{ijkl} + T_{jilk})
- Parameters
permutations
: The list of permutations to be applied simultaneously
- Note
Unlike the implementation in libtensor, this function includes the prefactor 0.5
-
std::shared_ptr<Tensor>
antisymmetrise
(const std::vector<std::vector<size_t>> &permutations) const = 0¶ Antisymmetrise with respect to the given index permutations and return the result. For details with respect to the format of the permutations, see symmetrise. In this case the output tensor will be antisymmetric with respect to these permutations.
- Note
Unlike the implementation in libtensor, this function includes the prefactor 0.5
-
void
evaluate
() const = 0¶ Make sure to evaluate this tensor in case this has not yet happened
-
bool
needs_evaluation
() const = 0¶ Check weather the object represents an evaluated expression or not.
-
void
set_immutable
() = 0¶ Flag the tensor as immutable, allows some optimisations to be performed
-
bool
is_mutable
() const = 0¶ Is the tensor mutable
-
void
set_element
(const std::vector<size_t> &idx, scalar_type value) = 0¶ Set the value of a single tensor element
- Note
This is a slow function and should be avoided.
-
scalar_type
get_element
(const std::vector<size_t> &tidx) const = 0¶ Get the value of a single tensor element
- Note
This is a slow function and should be avoided.
-
bool
is_element_allowed
(const std::vector<size_t> &tidx) const = 0¶ Return whether the element referenced by tidx is allowed (non-zero) by the symmetry of the Tensor or not.
-
std::vector<std::pair<std::vector<size_t>, scalar_type>>
select_n_absmax
(size_t n, bool unique_by_symmetry = true) const = 0¶ Get the n absolute largest elements along with their values
- Parameters
n
: Number of elements to selectunique_by_symmetry
: By default the returned elements are made unique by symmetry of the tensor. Setting this to false disables this feature.
-
std::vector<std::pair<std::vector<size_t>, scalar_type>>
select_n_absmin
(size_t n, bool unique_by_symmetry = true) const = 0¶ Get the n absolute smallest elements along with their values
- Parameters
n
: Number of elements to selectunique_by_symmetry
: By default the returned elements are made unique by symmetry of the tensor. Setting this to false disables this feature.
-
std::vector<std::pair<std::vector<size_t>, scalar_type>>
select_n_max
(size_t n, bool unique_by_symmetry = true) const = 0¶ Get the n largest elements along with their values
- Parameters
n
: Number of elements to selectunique_by_symmetry
: By default the returned elements are made unique by symmetry of the tensor. Setting this to false disables this feature.
-
std::vector<std::pair<std::vector<size_t>, scalar_type>>
select_n_min
(size_t n, bool unique_by_symmetry = true) const = 0¶ Get the n smallest elements along with their values
- Parameters
n
: Number of elements to selectunique_by_symmetry
: By default the returned elements are made unique by symmetry of the tensor. Setting this to false disables this feature.
-
void
export_to
(scalar_type *memptr, size_t size) const = 0¶ Extract the tensor to plain memory provided by the given pointer.
- Note
This will return a full, dense tensor. At least size elements of space are assumed at the provided memory location. The data is stored in row-major (C-like) format
-
void
export_to
(std::vector<scalar_type> &output) const¶ Extract the tensor into a std::vector, which will be resized to fit the data.
- Note
All data is stored in row-major (C-like) format
-
void
import_from
(const scalar_type *memptr, size_t size, scalar_type tolerance = 0, bool symmetry_check = true) = 0¶ Import the tensor from plain memory provided by the given pointer. The memory will be copied and all existing data overwritten. If symmetry_check is true, the process will check that the data has the required symmetry to fit into the tensor. This requires a slower algorithm to be chosen.
- Note
This function requires a full, dense tensor with the data stored in row-major (C-like) format.
- Parameters
memptr
: Full, dense memory pointer to the tensor data to be imported.size
: Size of the dense memory.tolerance
: Threshold to account for numerical inconsistencies when checking the symmetry or for determining zero blocks.symmetry_check
: Should symmetry be explicitly checked during the import.
-
void
import_from
(const std::vector<scalar_type> &input, scalar_type tolerance = 0, bool symmetry_check = true)¶ Import the tensor from plain memory provided by the given vector. The memory will be copied and all existing data overwritten. If symmetry_check is true, the process will check that the data has the required symmetry to fit into the tensor. This requires a slower algorithm to be chosen.
- Note
This function requires a full, dense tensor with the data stored in row-major (C-like) format.
- Parameters
input
: Input data in a linearised vectortolerance
: Threshold to account for numerical inconsistencies when checking the symmetry or for determining zero blocks.symmetry_check
: Should symmetry be explicitly checked during the import.
-
void
import_from
(std::function<void(const std::vector<std::pair<size_t, size_t>>&, scalar_type*)> generator, scalar_type tolerance = 0, bool symmetry_check = true, ) = 0¶ Import the tensor from a generator functor. All existing data will be overwritten. If symmetry_check is true, the process will check that the data has the required symmetry to fit into the tensor. This requires a slower algorithm to be chosen.
- Note
The generator is required to produce data in row-major (C-like) format at a designated memory location.
- Parameters
generator
: Generator functor. The functor is called with a list of ranges for each dimension. The ranges are half-open, left-inclusive and right-exclusive. The corresponding data should be written to the passed pointer into raw memory. This is an advanced functionality. Use only if you know what you are doing.tolerance
: Threshold to account for numerical inconsistencies when checking the symmetry or for determining zero blocks.symmetry_check
: Should symmetry be explicitly checked during the import.
-
std::string
describe_symmetry
() const = 0¶ Return a std::string providing hopefully helpful information about the symmetries stored inside the tensor object.
-
std::string
describe_expression
(std::string stage = "unoptimised") const = 0¶ Return a std::string providing hopefully helpful information about the expression tree stored inside the tensor object. Valid values for stage are “unoptimised”, “optimised”, “evaluation” returning the unoptimised expression tree as stored, an optimised form of it or the tree as used actually for evaluating the expression.
-
std::shared_ptr<const AdcMemory>
adcmem_ptr
() const¶ Return a pointer to the memory keep-alive object
-
const std::vector<AxisInfo> &
axes
() const¶ Return the axes info for this tensor
Public Functions
Construct a tensor interface class.
- Parameters
adcmem_ptr
: Memory to use for the tensor.axes
: Information about the tensor axes.
-
size_t
Utilities¶
Some random things to set up shop.
-
typedef double
real_type
¶ The real type to use
-
typedef double
scalar_type
¶ The scalar type to use
-
std::string
shape_to_string
(const std::vector<size_t> &shape)¶ Convert a std::vector<size_t> representing a tensor shape to a string
-
class
libadcc
::
AdcMemory
¶ - #include <AdcMemory.hh>
Manages memory for adc calculations.
This class has to remain in scope for as long as the allocated memory is needed. Normally this is achieved by keeping some sort of keep-alive shared pointer in the private context of all classes which require access to this memory. All changes to the environment are undone when the class goes out of scope.
- Note
Undefined behaviour may result if this class is ever initialised twice.
Public Functions
-
AdcMemory
(std::string pagefile_directory, size_t max_memory, size_t tbs_param = 16)¶ Setup the class and directly call initialise with the default allocator
-
std::string
allocator
() const¶ Return the allocator to which the class is initialised.
-
size_t
max_memory
() const¶ Return the max_memory parameter value to which the class was initialised.
- Note
This value is only a meaningful upper bound if allocator() != “standard”
-
std::string
pagefile_directory
() const¶ Return the pagefileprefix value
- Note
This value is only meaningful if allocator() != “standard”
-
size_t
tbs_param
() const¶ Return the tbs_param value
-
size_t
contraction_batch_size
() const¶ Get the contraction batch size, this is the number of tensor blocks, which are processed in a batch in a tensor contraction.
-
void
set_contraction_batch_size
(size_t bsize)¶ Set the contraction batch size
-
void
initialise
(std::string pagefile_directory, size_t max_memory, size_t tbs_param = 16, std::string allocator = "default")¶ Setup the environment for the memory management.
- Note
For some allocators
pagefile_directory
andmax_memory
are not supported and thus ignored.- Note
”libxm” and “libvmm” are extra features, which are not available in a default setup.
- Parameters
pagefile_directory
: File prefix for page filesmax_memory
: The maximal memory adc makes use of (in bytes).tbs_param
: The tensor block size parameter. This parameter roughly has the meaning of how many indices are handled together on operations. A good value seems to be 16.allocator
: The allocator to be used. Valid values are “libxm”, “libvmm”, “standard” and “default”, where “default” uses a default chosen from the first three.
-
struct
dimension_mismatch
: public invalid_argument¶
-
struct
not_implemented_error
: public runtime_error¶
-
class
libadcc
::
ThreadPool
¶ - #include <ThreadPool.hh>
Pool managing how many threads may be used for calculations.
Public Functions
-
ThreadPool
(size_t n_running, size_t n_total)¶ Initialise the thread pool.
- Parameters
n_running
: The total number of running threads to employn_total
: The total number of worker threads to use (Either running or ready)
-
void
reinit
(size_t n_running, size_t n_total)¶ Reinitialise the thread pool.
- Parameters
n_running
: The total number of running threads to employn_total
: The total number of worker threads to use (Either running or idle)
-
ThreadPool
()¶ Initialise a thread pool without parallelisation
-
size_t
n_running
() const¶ Return the number of running threads.
-
size_t
n_total
() const¶ Return the total number of worker threads (running or ready)
-
-
class
libadcc
::
Timer
¶ - #include <Timer.hh>
Minimalistic timer class: Just records the start and stop times for tasks, i.e. the intervals they ran
Public Functions
-
Timer
()¶ Construct an empty timer object
-
void
start
(const std::string &task)¶ Start a particular task
-
double
stop
(const std::string &task)¶ Stop the task again
Public Members
-
double
time_construction
¶ Time when this class was constructed
-
std::map<std::string, std::vector<std::pair<double, double>>>
intervals
¶ The intervals stored for a particular key
-
std::map<std::string, double>
start_times
¶ The start time stored for each key
Public Static Functions
-
double
now
()¶ A float representing the current time on the clock used for measurement. For consistency, this function should always be used to obtain a representation of “now”. The unit is seconds.
-
-
class
libadcc
::
RecordTime
¶ - #include <Timer.hh>
RAI-like class for timing tasks. Construction starts the timer, destruction ends it automatically. It is assumed that the passed task has a longer lifetime than the RecordTime object
Public Functions
Tensor implementation using libtensor¶
This section describes the implementation of the
Tensor functionality of libadcc::Tensor
using the libtensor tensor library.
-
std::map<size_t, std::string>
setup_point_group_table
(libtensor::product_table_container &ptc, const std::string &point_group)¶ Setup point group symmetry table inside libtensor and return the irrep mapping as libtensor needs it.
-
template<size_t
N
>
libtensor::bispace<N>as_bispace
(const std::vector<AxisInfo> &axes)¶ Form a bispace from a list of AxisInfo representing the axes of a tensor
Translate a adcc::Symmetry object into a shared pointer of the libtensor symmetry object
Construction of the tensor implementation class
- Note
This constructor is for experts only. Non-experts should use the make_tensor method instead.
-
std::shared_ptr<Tensor>
empty_like
() const override¶ Return a new tensor with the same dimensionality and symmetry as the passed tensor. The elements have undefined values.
-
std::shared_ptr<Tensor>
nosym_like
() const override¶ Return a new tensor with same dimensionality and shape (and bispace) but without any symmetry copied over
-
void
set_mask
(std::string mask, scalar_type value) override¶ Set a mask into a tensor to a specific value. The mask to set is defined by the mask string. Repetitive indices define the values to be set. E.g. for a 6D tensor the mask string “iijkli” would set all elements T_{iijkli} for all ijkl to the given value.
-
std::shared_ptr<Tensor>
diagonal
(std::vector<size_t> axes) override¶ Extract a generalised diagonal from this Tensor.
axes
defines the axes over which to extract. The resulting new axis is appended to the remaining axis. For example if we have a Tensor fourth ordert
then diagonal(1, 2) returnsd_{iax} = t_{ixxa}
and diagonal(0, 1, 2) returnsd_{ax} = t_{xxxa}
.
-
void
set_random
() override¶ Set the tensor to random data preserving symmetry
-
std::shared_ptr<Tensor>
scale
(scalar_type c) const override¶ Scale tensor by a scalar value and return the result
Add another tensor and return the result
Multiply another tensor elementwise and return the result
Divide another tensor elementwise and return the result
Add a linear combination to the represented tensor and store the result. Note that unlike
add
andscale
this operation is not lazy and is provided for optimising on cases where a lazy linear combination introduces a too large overhead.
-
std::shared_ptr<Tensor>
copy
() const override¶ Return a deep copy of this tensor.
-
std::shared_ptr<Tensor>
transpose
(std::vector<size_t> axes) const override¶ Return a transpose form of this tensor
- Parameters
axes
: New permutation of the dimensions of this tensor (e.g. (1,0,2,3) will permute dimension 0 and 1 and keep the others) If missing, just reverse the dimension order.
Compute the direct sum of this tensor with another tensor. All axes are used in the order of appearance (this and then other).
Return the tensordot of this tensor with another
- Parameters
axes
: Axes to contract over. The first vector refers to the axes of this tensor, the second vector to the axes of the other tensor.
Compute the Frobenius or l2 inner product with a list of other tensors, returning the results as a scalar array
-
std::shared_ptr<Tensor>
symmetrise
(const std::vector<std::vector<size_t>> &permutations) const override¶ Symmetrise with respect to the given index permutations by adding the elements resulting from appropriate index permutations.
Examples for permutations. Take a rank-4 tensor T_{ijkl} {{0,1}} Permute the first two indices, i.e. form 0.5 * (T_{ijkl} + T_{jikl}) {{0,1}, {2,3}} Form 0.5 * (T_{ijkl} + T_{jilk})
- Parameters
permutations
: The list of permutations to be applied simultaneously
- Note
Unlike the implementation in libtensor, this function includes the prefactor 0.5
-
std::shared_ptr<Tensor>
antisymmetrise
(const std::vector<std::vector<size_t>> &permutations) const override¶ Antisymmetrise with respect to the given index permutations and return the result. For details with respect to the format of the permutations, see symmetrise. In this case the output tensor will be antisymmetric with respect to these permutations.
- Note
Unlike the implementation in libtensor, this function includes the prefactor 0.5
-
double
trace
(std::string contraction) const override¶ Return the full trace of a tensor
- Parameters
contraction
: Indices to contract over using Einstein summation convention, e.g. “abab” traces over the 1st and 3rd and 2nd and 4th axis.
-
void
evaluate
() const override¶ Make sure to evaluate this tensor in case this has not yet happened
-
bool
needs_evaluation
() const override¶ Check weather the object represents an evaluated expression or not.
-
void
set_immutable
() override¶ Flag the tensor as immutable, allows some optimisations to be performed
-
bool
is_mutable
() const override¶ Is the tensor mutable
-
void
set_element
(const std::vector<size_t> &tidx, scalar_type value) override¶ Set the value of a single tensor element
- Note
This is a slow function and should be avoided.
-
scalar_type
get_element
(const std::vector<size_t> &tidx) const override¶ Get the value of a single tensor element
- Note
This is a slow function and should be avoided.
-
bool
is_element_allowed
(const std::vector<size_t> &tidx) const override¶ Return whether the element referenced by tidx is allowed (non-zero) by the symmetry of the Tensor or not.
-
std::vector<std::pair<std::vector<size_t>, scalar_type>>
select_n_absmax
(size_t n, bool unique_by_symmetry) const override¶ Get the n absolute largest elements along with their values
- Parameters
n
: Number of elements to selectunique_by_symmetry
: By default the returned elements are made unique by symmetry of the tensor. Setting this to false disables this feature.
-
std::vector<std::pair<std::vector<size_t>, scalar_type>>
select_n_absmin
(size_t n, bool unique_by_symmetry) const override¶ Get the n absolute smallest elements along with their values
- Parameters
n
: Number of elements to selectunique_by_symmetry
: By default the returned elements are made unique by symmetry of the tensor. Setting this to false disables this feature.
-
std::vector<std::pair<std::vector<size_t>, scalar_type>>
select_n_max
(size_t n, bool unique_by_symmetry) const override¶ Get the n largest elements along with their values
- Parameters
n
: Number of elements to selectunique_by_symmetry
: By default the returned elements are made unique by symmetry of the tensor. Setting this to false disables this feature.
-
std::vector<std::pair<std::vector<size_t>, scalar_type>>
select_n_min
(size_t n, bool unique_by_symmetry) const override¶ Get the n smallest elements along with their values
- Parameters
n
: Number of elements to selectunique_by_symmetry
: By default the returned elements are made unique by symmetry of the tensor. Setting this to false disables this feature.
-
void
export_to
(scalar_type *memptr, size_t size) const override¶ Extract the tensor to plain memory provided by the given pointer.
- Note
This will return a full, dense tensor. At least size elements of space are assumed at the provided memory location. The data is stored in row-major (C-like) format
-
void
import_from
(const scalar_type *memptr, size_t size, scalar_type tolerance, bool symmetry_check) override¶ Import the tensor from plain memory provided by the given pointer. The memory will be copied and all existing data overwritten. If symmetry_check is true, the process will check that the data has the required symmetry to fit into the tensor. This requires a slower algorithm to be chosen.
- Note
This function requires a full, dense tensor with the data stored in row-major (C-like) format.
- Parameters
memptr
: Full, dense memory pointer to the tensor data to be imported.size
: Size of the dense memory.tolerance
: Threshold to account for numerical inconsistencies when checking the symmetry or for determining zero blocks.symmetry_check
: Should symmetry be explicitly checked during the import.
-
void
import_from
(std::function<void(const std::vector<std::pair<size_t, size_t>>&, scalar_type*)> generator, scalar_type tolerance, bool symmetry_check, ) override¶ Import the tensor from a generator functor. All existing data will be overwritten. If symmetry_check is true, the process will check that the data has the required symmetry to fit into the tensor. This requires a slower algorithm to be chosen.
- Note
The generator is required to produce data in row-major (C-like) format at a designated memory location.
- Parameters
generator
: Generator functor. The functor is called with a list of ranges for each dimension. The ranges are half-open, left-inclusive and right-exclusive. The corresponding data should be written to the passed pointer into raw memory. This is an advanced functionality. Use only if you know what you are doing.tolerance
: Threshold to account for numerical inconsistencies when checking the symmetry or for determining zero blocks.symmetry_check
: Should symmetry be explicitly checked during the import.
-
std::string
describe_symmetry
() const override¶ Return a std::string providing hopefully helpful information about the symmetries stored inside the tensor object.
-
std::string
describe_expression
(std::string stage = "unoptimised") const override¶ Return a std::string providing hopefully helpful information about the expression tree stored inside the tensor object. Valid values for stage are “unoptimised”, “optimised”, “evaluation” returning the unoptimised expression tree as stored, an optimised form of it or the tree as used actually for evaluating the expression.
-
operator libtensor::btensor<N, scalar_type>&
()¶ Convert object to btensor for use in libtensor functions.
-
std::shared_ptr<libtensor::btensor<N, scalar_type>>
libtensor_ptr
()¶ Return inner btensor object
- Note
This is an advanced function. Use with care, since this can corrupt the internal state of the Tensor object.
-
std::shared_ptr<ExpressionTree>
expression_ptr
() const¶ Return permutation, expression object and required keep-alives to this tensor
- Note
This is an advanced function. Use with care, since this can corrupt the internal state of the Tensor object.
Extractor function to convert the contained tensor from adcc::Tensor to an expression.
- Note
This is an internal function. Use only if you know what you are doing.
Extractor functions to convert the contained tensor of an adcc::Tensor object to an appropriate libtensor::btensor, The dimensionality N has to match
Extractor functions to convert the contained tensor of an adcc::Tensor object to an appropriate libtensor::btensor, The dimensionality N has to match Extractor function to convert the contained tensor of an adcc::Tensor object to an appropriate pointer to libtensor::btensor. The dimensionality N has to match.
-
template<size_t
N
>
classlibadcc
::
TensorImpl
: public libadcc::Tensor¶ - #include <TensorImpl.hh>
Actual implementation of the Tensor class. For details about the functions see Tensor.hh
Public Functions
Construction of the tensor implementation class
- Note
This constructor is for experts only. Non-experts should use the make_tensor method instead.
-
std::shared_ptr<Tensor>
empty_like
() const override¶ Return a new tensor with the same dimensionality and symmetry as the passed tensor. The elements have undefined values.
-
std::shared_ptr<Tensor>
nosym_like
() const override¶ Return a new tensor with same dimensionality and shape (and bispace) but without any symmetry copied over
-
void
set_mask
(std::string mask, scalar_type value) override¶ Set a mask into a tensor to a specific value. The mask to set is defined by the mask string. Repetitive indices define the values to be set. E.g. for a 6D tensor the mask string “iijkli” would set all elements T_{iijkli} for all ijkl to the given value.
-
std::shared_ptr<Tensor>
diagonal
(std::vector<size_t> axes) override¶ Extract a generalised diagonal from this Tensor.
axes
defines the axes over which to extract. The resulting new axis is appended to the remaining axis. For example if we have a Tensor fourth ordert
then diagonal(1, 2) returnsd_{iax} = t_{ixxa}
and diagonal(0, 1, 2) returnsd_{ax} = t_{xxxa}
.
-
void
set_random
() override¶ Set the tensor to random data preserving symmetry
-
std::shared_ptr<Tensor>
scale
(scalar_type c) const override¶ Scale tensor by a scalar value and return the result
Add another tensor and return the result
Multiply another tensor elementwise and return the result
Divide another tensor elementwise and return the result
Add a linear combination to the represented tensor and store the result. Note that unlike
add
andscale
this operation is not lazy and is provided for optimising on cases where a lazy linear combination introduces a too large overhead.
-
std::shared_ptr<Tensor>
transpose
(std::vector<size_t> axes) const override¶ Return a transpose form of this tensor
- Parameters
axes
: New permutation of the dimensions of this tensor (e.g. (1,0,2,3) will permute dimension 0 and 1 and keep the others) If missing, just reverse the dimension order.
Compute the direct sum of this tensor with another tensor. All axes are used in the order of appearance (this and then other).
Return the tensordot of this tensor with another
- Parameters
axes
: Axes to contract over. The first vector refers to the axes of this tensor, the second vector to the axes of the other tensor.
Compute the Frobenius or l2 inner product with a list of other tensors, returning the results as a scalar array
-
std::shared_ptr<Tensor>
symmetrise
(const std::vector<std::vector<size_t>> &permutations) const override¶ Symmetrise with respect to the given index permutations by adding the elements resulting from appropriate index permutations.
Examples for permutations. Take a rank-4 tensor T_{ijkl} {{0,1}} Permute the first two indices, i.e. form 0.5 * (T_{ijkl} + T_{jikl}) {{0,1}, {2,3}} Form 0.5 * (T_{ijkl} + T_{jilk})
- Parameters
permutations
: The list of permutations to be applied simultaneously
- Note
Unlike the implementation in libtensor, this function includes the prefactor 0.5
-
std::shared_ptr<Tensor>
antisymmetrise
(const std::vector<std::vector<size_t>> &permutations) const override¶ Antisymmetrise with respect to the given index permutations and return the result. For details with respect to the format of the permutations, see symmetrise. In this case the output tensor will be antisymmetric with respect to these permutations.
- Note
Unlike the implementation in libtensor, this function includes the prefactor 0.5
-
double
trace
(std::string contraction) const override¶ Return the full trace of a tensor
- Parameters
contraction
: Indices to contract over using Einstein summation convention, e.g. “abab” traces over the 1st and 3rd and 2nd and 4th axis.
-
void
evaluate
() const override¶ Make sure to evaluate this tensor in case this has not yet happened
-
bool
needs_evaluation
() const override¶ Check weather the object represents an evaluated expression or not.
-
void
set_immutable
() override¶ Flag the tensor as immutable, allows some optimisations to be performed
-
bool
is_mutable
() const override¶ Is the tensor mutable
-
void
set_element
(const std::vector<size_t> &tidx, scalar_type value) override¶ Set the value of a single tensor element
- Note
This is a slow function and should be avoided.
-
scalar_type
get_element
(const std::vector<size_t> &tidx) const override¶ Get the value of a single tensor element
- Note
This is a slow function and should be avoided.
-
bool
is_element_allowed
(const std::vector<size_t> &tidx) const override¶ Return whether the element referenced by tidx is allowed (non-zero) by the symmetry of the Tensor or not.
-
std::vector<std::pair<std::vector<size_t>, scalar_type>>
select_n_absmax
(size_t n, bool unique_by_symmetry) const override¶ Get the n absolute largest elements along with their values
- Parameters
n
: Number of elements to selectunique_by_symmetry
: By default the returned elements are made unique by symmetry of the tensor. Setting this to false disables this feature.
-
std::vector<std::pair<std::vector<size_t>, scalar_type>>
select_n_absmin
(size_t n, bool unique_by_symmetry) const override¶ Get the n absolute smallest elements along with their values
- Parameters
n
: Number of elements to selectunique_by_symmetry
: By default the returned elements are made unique by symmetry of the tensor. Setting this to false disables this feature.
-
std::vector<std::pair<std::vector<size_t>, scalar_type>>
select_n_max
(size_t n, bool unique_by_symmetry) const override¶ Get the n largest elements along with their values
- Parameters
n
: Number of elements to selectunique_by_symmetry
: By default the returned elements are made unique by symmetry of the tensor. Setting this to false disables this feature.
-
std::vector<std::pair<std::vector<size_t>, scalar_type>>
select_n_min
(size_t n, bool unique_by_symmetry) const override¶ Get the n smallest elements along with their values
- Parameters
n
: Number of elements to selectunique_by_symmetry
: By default the returned elements are made unique by symmetry of the tensor. Setting this to false disables this feature.
-
void
export_to
(scalar_type *memptr, size_t size) const override¶ Extract the tensor to plain memory provided by the given pointer.
- Note
This will return a full, dense tensor. At least size elements of space are assumed at the provided memory location. The data is stored in row-major (C-like) format
-
void
import_from
(const scalar_type *memptr, size_t size, scalar_type tolerance, bool symmetry_check) override¶ Import the tensor from plain memory provided by the given pointer. The memory will be copied and all existing data overwritten. If symmetry_check is true, the process will check that the data has the required symmetry to fit into the tensor. This requires a slower algorithm to be chosen.
- Note
This function requires a full, dense tensor with the data stored in row-major (C-like) format.
- Parameters
memptr
: Full, dense memory pointer to the tensor data to be imported.size
: Size of the dense memory.tolerance
: Threshold to account for numerical inconsistencies when checking the symmetry or for determining zero blocks.symmetry_check
: Should symmetry be explicitly checked during the import.
-
void
import_from
(std::function<void(const std::vector<std::pair<size_t, size_t>>&, scalar_type*)> generator, scalar_type tolerance, bool symmetry_check, ) override¶ Import the tensor from a generator functor. All existing data will be overwritten. If symmetry_check is true, the process will check that the data has the required symmetry to fit into the tensor. This requires a slower algorithm to be chosen.
- Note
The generator is required to produce data in row-major (C-like) format at a designated memory location.
- Parameters
generator
: Generator functor. The functor is called with a list of ranges for each dimension. The ranges are half-open, left-inclusive and right-exclusive. The corresponding data should be written to the passed pointer into raw memory. This is an advanced functionality. Use only if you know what you are doing.tolerance
: Threshold to account for numerical inconsistencies when checking the symmetry or for determining zero blocks.symmetry_check
: Should symmetry be explicitly checked during the import.
-
std::string
describe_symmetry
() const override¶ Return a std::string providing hopefully helpful information about the symmetries stored inside the tensor object.
-
std::string
describe_expression
(std::string stage = "unoptimised") const override¶ Return a std::string providing hopefully helpful information about the expression tree stored inside the tensor object. Valid values for stage are “unoptimised”, “optimised”, “evaluation” returning the unoptimised expression tree as stored, an optimised form of it or the tree as used actually for evaluating the expression.
-
operator libtensor::btensor<N, scalar_type>&
()¶ Convert object to btensor for use in libtensor functions.
-
std::shared_ptr<libtensor::btensor<N, scalar_type>>
libtensor_ptr
()¶ Return inner btensor object
- Note
This is an advanced function. Use with care, since this can corrupt the internal state of the Tensor object.
-
std::shared_ptr<ExpressionTree>
expression_ptr
() const¶ Return permutation, expression object and required keep-alives to this tensor
- Note
This is an advanced function. Use with care, since this can corrupt the internal state of the Tensor object.
-
std::vector<scalar_type>
dot
(std::vector<std::shared_ptr<Tensor>> tensors) const = 0 Compute the Frobenius or l2 inner product with a list of other tensors, returning the results as a scalar array
-
void
export_to
(scalar_type *memptr, size_t size) const = 0 Extract the tensor to plain memory provided by the given pointer.
- Note
This will return a full, dense tensor. At least size elements of space are assumed at the provided memory location. The data is stored in row-major (C-like) format
-
void
export_to
(std::vector<scalar_type> &output) const¶ Extract the tensor into a std::vector, which will be resized to fit the data.
- Note
All data is stored in row-major (C-like) format
-
void
import_from
(const scalar_type *memptr, size_t size, scalar_type tolerance = 0, bool symmetry_check = true) = 0 Import the tensor from plain memory provided by the given pointer. The memory will be copied and all existing data overwritten. If symmetry_check is true, the process will check that the data has the required symmetry to fit into the tensor. This requires a slower algorithm to be chosen.
- Note
This function requires a full, dense tensor with the data stored in row-major (C-like) format.
- Parameters
memptr
: Full, dense memory pointer to the tensor data to be imported.size
: Size of the dense memory.tolerance
: Threshold to account for numerical inconsistencies when checking the symmetry or for determining zero blocks.symmetry_check
: Should symmetry be explicitly checked during the import.
-
void
import_from
(const std::vector<scalar_type> &input, scalar_type tolerance = 0, bool symmetry_check = true)¶ Import the tensor from plain memory provided by the given vector. The memory will be copied and all existing data overwritten. If symmetry_check is true, the process will check that the data has the required symmetry to fit into the tensor. This requires a slower algorithm to be chosen.
- Note
This function requires a full, dense tensor with the data stored in row-major (C-like) format.
- Parameters
input
: Input data in a linearised vectortolerance
: Threshold to account for numerical inconsistencies when checking the symmetry or for determining zero blocks.symmetry_check
: Should symmetry be explicitly checked during the import.
-
void
import_from
(std::function<void(const std::vector<std::pair<size_t, size_t>>&, scalar_type*)> generator, scalar_type tolerance = 0, bool symmetry_check = true, ) = 0 Import the tensor from a generator functor. All existing data will be overwritten. If symmetry_check is true, the process will check that the data has the required symmetry to fit into the tensor. This requires a slower algorithm to be chosen.
- Note
The generator is required to produce data in row-major (C-like) format at a designated memory location.
- Parameters
generator
: Generator functor. The functor is called with a list of ranges for each dimension. The ranges are half-open, left-inclusive and right-exclusive. The corresponding data should be written to the passed pointer into raw memory. This is an advanced functionality. Use only if you know what you are doing.tolerance
: Threshold to account for numerical inconsistencies when checking the symmetry or for determining zero blocks.symmetry_check
: Should symmetry be explicitly checked during the import.