The adccore C++ layer

This page contains a reference of the adccore C++ library and its classes and functions generated automatically from the adccore source code.

Reference state

This category lists the adccore functionality, which imports the data from the adcc::HartreeFockSolution_i interface into the adcc::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 adcc::MoSpaces, which collects information about the occupied and virtual orbital spaces, and adcc::MoIndexTranslation, which maps orbitals indices between the ordering used by adccore and the one used by the SCF program.

std::vector<size_t> adcc::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 starts
  • length: Total number of indices
  • max_block_size: Maximal size of a block

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

Unnamed Group

size_t start() const

Start of the range, i.e. the first index

Unnamed Group

size_t end() const

End of the range, i.e. the index past the last index

Public Functions

size_t length() const

Length of the range, i.e. number of indices contained

struct 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

std::vector<size_t> ends() const

Get the vector of all range ends

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

Unnamed Group

SimpleRange &from()

The range from which we map

Unnamed Group

SimpleRange &to()

The range to which we map

class 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

MoIndexTranslation(std::shared_ptr<const MoSpaces> mospaces_ptr, const std::string &space)

Create an MoIndexTranslation object for the provided space string.

MoIndexTranslation(std::shared_ptr<const MoSpaces> mospaces_ptr, const std::vector<std::string> &subspaces)

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

MoSpaces(const HartreeFockSolution_i &hf, std::shared_ptr<const AdcMemory> adcmem_ptr, std::vector<size_t> core_orbitals, std::vector<size_t> frozen_core_orbitals, std::vector<size_t> frozen_virtuals)

Construct an MoSpaces object from a HartreeFockSolution_i, a pointer to an AdcMemory object.

Parameters
  • hf: HartreeFockSolution_i object to use
  • adcmem_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.

MoSpaces(const HartreeFockSolution_i &hf, std::shared_ptr<const AdcMemory> adcmem_ptr, size_t n_core_orbitals = 0, size_t n_frozen_core = 0, size_t n_frozen_virtuals = 0)

Construct an MoSpaces object from a HartreeFockSolution_i and a maximal block size to be used for the tensors. Convenience wrapper to above constructor, which selects the orbital indices automatically on the upper / lower end of the energy range.

Parameters
  • hf: HartreeFockSolution_i object to use
  • adcmem_ptr: ADC memory keep-alive object to be used in all Tensors constructed using this MoSpaces object.
  • n_core_orbitals: Number of core orbitals to use.
  • n_frozen_core: Number of orbitals to put in a frozen core.
  • n_frozen_virtuals: Number of frozen virtual orbitals.

ctx::CtxMap to_ctx() const

Serialise the object into a ctx::CtxMap to be used inside e.g. adcman

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)

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

const Timer &timer() const

Obtain the timing info contained in this object

ctx::CtxMap to_ctx() const

Serialise the object into a ctx::CtxMap to be used inside e.g. adcman

Perturbation theory

adcc::LazyMp lazily computes second-order and third-order Møller-Plesset perturbation theory on top of the reference held by a adcc::ReferenceState.

struct CachingPolicy_i
#include <CachingPolicy.hh>

Policy class determining, which tensors are cached and which are not

Subclassed by adcc::AlwaysCachePolicy

Public Functions

virtual bool should_cache(const std::string &tensor_label, const std::string &tensor_space, const std::string &leading_order_contraction) = 0

Should a tensor be stored i.e. cached in memory according to this policy object (true) or not (false).

Parameters
  • tensor_label: Label of the tensor object (e.g. t2, t2eri, …)
  • tensor_space: The spaces string of the tensor to estimate memory requirement (e.g. “o1o1v1v1” for “t2_o1o1v1v1”)
  • leading_order_contraction: The spaces involved in the most expensive contraction to compute this tensor to estimate computational scaling for computing the tensor. (e.g. “o1o1v1v1” for “t2_o1o1v1v1”)

class AlwaysCachePolicy : public adcc::CachingPolicy_i
#include <CachingPolicy.hh>

Policy class, which always caches everything

struct LazyMp
#include <LazyMp.hh>

A lazy Møller-Plesset perturbation theory class. The first time a particular functionality is requested, it is computed

Public Functions

LazyMp(std::shared_ptr<const ReferenceState> ref_ptr)

Initialise a LazyMp object using a reference state. Will not compute anything unless explicitly requested

LazyMp(std::shared_ptr<const ReferenceState> ref_ptr, std::shared_ptr<CachingPolicy_i> caching_policy_ptr)

Initialise a LazyMp object using a reference state. Will not compute anything unless explicitly requested

LazyMp(std::shared_ptr<const ReferenceState> ref_ptr, ctx::CtxMap map)

Initialise a LazyMp object using a reference state and some cached data represented by the passed context object

double energy_correction(int level) const

Obtain the MP energy correction at a particular level

std::shared_ptr<Tensor> df(std::string space) const

Delta Fock matrix

std::shared_ptr<Tensor> t2(std::string space) const

T2 amplitudes

void set_t2(std::string space, std::shared_ptr<Tensor> tensor)

Set the t2 amplitudes tensor. This invalidates all data depending on the T2 amplitudes in this class.

Note
Potential other caches, such as computed ADC intermediates are not automatically invalidated.

std::shared_ptr<Tensor> td2(std::string space) const

Return the T^D_2 term

std::shared_ptr<Tensor> mp2_diffdm(std::string space) const

Return MP2 difference densities (in MOs)

std::shared_ptr<const OneParticleOperator> mp2_diffdm_ptr() const

Return MP2 difference density relative to HF as a OneParticleOperator object

std::shared_ptr<Tensor> t2eri(const std::string &space, const std::string &contraction) const

Return the T2 tensor with ERI tensor contraction intermediates (called pi1 to pi7 in adcman)

std::shared_ptr<const ReferenceState> reference_state_ptr() const

Access to the reference state on top of which the MP calculation has been performed.

std::shared_ptr<CachingPolicy_i> caching_policy_ptr() const

Access to the caching policy used

void set_caching_policy(std::shared_ptr<CachingPolicy_i> newcpol_ptr)

Set the caching policy used from now on

ctx::CtxMap to_ctx() const

Return the data represented by this datastructure as a ctx::CtxMap

std::shared_ptr<const MoSpaces> mospaces_ptr() const

Return the MoSpaces to which this LazyMp is set up

bool has_core_occupied_space() const

Has this Mp container the core-valence separation applied

const Timer &timer() const

Obtain the timing info contained in this object

AdcMatrix and matrix cores

adcc::AdcMatrix sets up a representation of the ADC matrix for a particular method. My the means of matrix cores, which actually do the work, this allows to perform matrix-vector products or access the diagonal of such a matrix under a common interface.

std::ostream &adcc::operator<<(std::ostream &o, const AdcIntermediates &im)

Write a descriptive string for the AdcIntermediates object to the output stream.

class AdcIntermediates
#include <AdcIntermediates.hh>

Class managing intermediate expressions of the ADC matrix, which are used frequently and/or expensive to compute.

Cache functions

These functions compute the intermediates (if not already done) and cache them inside the class if the CachingPolicy says so. In this case another call to these functions will just return the pointers to the appropriate object without extra computation. In other words setting the pointers in this class to a nullptr will always trigger another computation.

std::shared_ptr<Tensor> compute_adc2_i1()

The ADC(2) i1 intermediate

std::shared_ptr<Tensor> compute_adc2_i2()

The ADC(2) i2 intermediate

std::shared_ptr<Tensor> compute_adc3_pia()

The ADC(3) N^6 intermediate A

std::shared_ptr<Tensor> compute_adc3_pib()

The ADC(3) N^6 intermediate B

std::shared_ptr<Tensor> compute_adc3_m11()

The ADC(3) m11 intermediate (== singles block of ADC(3) matrix

std::shared_ptr<Tensor> compute_cv_p_oo()

The CVS-ADC(2) cv_p_oo intermediate for properties

std::shared_ptr<Tensor> compute_cv_p_ov()

The CVS-ADC(2) cv_p_ov intermediate for properties

std::shared_ptr<Tensor> compute_cv_p_vv()

The CVS-ADC(2) cv_p_vv intermediate for properties

std::shared_ptr<Tensor> compute_cvs_adc3_m11()

The CVS-ADC(3) m11 intermediate (== singles block of CVS-ADC(3) matrix

Public Functions

std::shared_ptr<const LazyMp> ground_state_ptr() const

Access to the ground state on top of which the ADC intermediates have been computed

const Timer &timer() const

Obtain the timing info contained in this object

Public Members

std::shared_ptr<Tensor> adc2_i1_ptr

The ADC(2) i1 intermediate

std::shared_ptr<Tensor> adc2_i2_ptr

The ADC(2) i2 intermediate

std::shared_ptr<Tensor> adc3_pia_ptr

The ADC(3) N^6 intermediate A

std::shared_ptr<Tensor> adc3_pib_ptr

The ADC(3) N^6 intermediate B

std::shared_ptr<Tensor> adc3_m11_ptr

The ADC(3) m11 intermediate (== singles block of ADC(3) matrix

std::shared_ptr<Tensor> cv_p_oo_ptr

The CVS-ADC(2) cv_p_oo intermediate for properties

std::shared_ptr<Tensor> cv_p_ov_ptr

The CVS-ADC(2) cv_p_ov intermediate for properties

std::shared_ptr<Tensor> cv_p_vv_ptr

The CVS-ADC(2) cv_p_vv intermediate for properties

std::shared_ptr<Tensor> cvs_adc3_m11_ptr

The CVS-ADC(3) m11 intermediate (== singles block of CVS-ADC(3) matrix

class AdcMatrixCoreBase
#include <AdcMatrixCoreBase.hh>

AdcMatrixCoreBase provides the interface to the actual implementations of the ADC matrix functionality

Subclassed by adcc::AdcMatrixCoreSingles, adcc::AdcMatrixCoreSinglesDoubles

Public Functions

AdcMatrixCoreBase(std::string name, std::shared_ptr<const LazyMp> mp_ptr_)

Construct the AdcMatrixCoreBase object

Parameters
  • name: Name of the ADC method
  • mp_ptr_: The LazyMp pointer describing the MP perturbation theory results

virtual std::vector<size_t> shape() const = 0

Return the shape of the represented matrix

virtual bool has_block(const std::string &block) const = 0

Return true if the referenced block exists in this matrix

bool has_block(char block) const

Return true if the referenced block exists in this matrix.

virtual std::vector<std::string> blocks() const = 0

Return the list of blocks, which exist in this matrix

virtual std::vector<std::string> block_spaces(const std::string &block) const = 0

Return the list of spaces involved for a particular block of the matrix

virtual void compute_apply_ss(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const = 0

Compute the application of the singles-singles block to the vector in, returning the result in the vector out.

virtual void compute_apply_sd(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const = 0

Compute the application of the singles-doubles block to the vector in, returning the result in the vector out.

Note
Might not be implemented for all matrix cores.

virtual void compute_apply_ds(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const = 0

Compute the application of the doubles-singles block to the vector in, returning the result in the vector out.

Note
Might not be implemented for all matrix cores.

virtual void compute_apply_dd(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const = 0

Compute the application of the doubles-doubles block to the vector in, returning the result in the vector out.

Note
Might not be implemented for all matrix cores.

void compute_apply(std::string block, const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the given block of the ADC matrix.

virtual std::shared_ptr<Tensor> diagonal_s() const = 0

Compute the singles diagonal of the ADC matrix and return it.

virtual std::shared_ptr<Tensor> diagonal_d() const = 0

Compute the doubles diagonal of the ADC matrix and return it.

Note
Might not be implemented for all matrix cores.

std::shared_ptr<Tensor> diagonal(std::string block) const

Compute the block-wise diagonal of ADC matrix.

virtual void compute_matvec(const AmplitudeVector &ins, const AmplitudeVector &outs) const

Compute the application of the ADC matrix to the amplitude vector supplied on ins and return the result in outs.

std::shared_ptr<AdcIntermediates> get_intermediates_ptr() const

Get the pointer to the AdcIntermediates cache for the precomputed intermediates

void set_intermediates_ptr(std::shared_ptr<AdcIntermediates> im_ptr)

Replace the currently held AdcIntermediates cache object by the supplied one. This allows to reuse intermediates from a previous calculation.

std::shared_ptr<const ReferenceState> reference_state_ptr()

Access to the reference state on top of which the ADC calculation has been started

std::shared_ptr<const LazyMp> ground_state_ptr()

Obtain the pointer to the MP ground state upon which this ADC matrix is based.

class AdcMatrixCoreSingles : public adcc::AdcMatrixCoreBase
#include <AdcMatrixCoreSingles.hh>

Specialisation of AdcMatrixCoreBase for ADC matrices with a singles part only

Subclassed by adcc::CvsAdcMatrixCoreSingles, adcc::GenericAdcMatrixCoreSingles

Public Functions

bool has_block(const std::string &block) const

Return true if the referenced block exists in this matrix

std::vector<std::string> blocks() const

Return the list of blocks, which exist in this matrix

void compute_apply_sd(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the singles-doubles block to the vector in, returning the result in the vector out.

Note
Might not be implemented for all matrix cores.

void compute_apply_ds(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the doubles-singles block to the vector in, returning the result in the vector out.

Note
Might not be implemented for all matrix cores.

void compute_apply_dd(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the doubles-doubles block to the vector in, returning the result in the vector out.

Note
Might not be implemented for all matrix cores.

std::shared_ptr<Tensor> diagonal_d() const

Compute the doubles diagonal of the ADC matrix and return it.

Note
Might not be implemented for all matrix cores.

void compute_matvec(const AmplitudeVector &ins, const AmplitudeVector &outs) const

Compute the application of the ADC matrix to the amplitude vector supplied on ins and return the result in outs.

class AdcMatrixCoreSinglesDoubles : public adcc::AdcMatrixCoreBase
#include <AdcMatrixCoreSinglesDoubles.hh>

Specialisation of AdcMatrixCoreBase for ADC matrices with a singles and a doubles part

Subclassed by adcc::CvsAdcMatrixCoreSinglesDoubles, adcc::GenericAdcMatrixCoreSinglesDoubles

Public Functions

bool has_block(const std::string &block) const

Return true if the referenced block exists in this matrix

std::vector<std::string> blocks() const

Return the list of blocks, which exist in this matrix

class CvsAdcMatrixCoreSingles : public adcc::AdcMatrixCoreSingles
#include <CvsAdcMatrixCores.hh>

Base class for CVS matrix cores with only a singles block

Subclassed by adcc::CvsAdc0MatrixCore, adcc::CvsAdc1MatrixCore

Public Functions

std::vector<size_t> shape() const

Return the shape of the represented matrix

std::vector<std::string> block_spaces(const std::string &block) const

Return the list of spaces involved for a particular block of the matrix

class CvsAdcMatrixCoreSinglesDoubles : public adcc::AdcMatrixCoreSinglesDoubles
#include <CvsAdcMatrixCores.hh>

Base class for CVS matrix cores with a singles and a double block

Subclassed by adcc::CvsAdc2MatrixCore

Public Functions

std::vector<size_t> shape() const

Return the shape of the represented matrix

std::vector<std::string> block_spaces(const std::string &block) const

Return the list of spaces involved for a particular block of the matrix

class CvsAdc0MatrixCore : public adcc::CvsAdcMatrixCoreSingles
#include <CvsAdcMatrixCores.hh>

Matrix core for CVS-ADC(0)

Public Functions

std::shared_ptr<Tensor> diagonal_s() const

Compute the singles diagonal of the ADC matrix and return it.

void compute_apply_ss(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the singles-singles block to the vector in, returning the result in the vector out.

class CvsAdc1MatrixCore : public adcc::CvsAdcMatrixCoreSingles
#include <CvsAdcMatrixCores.hh>

Matrix core for CVS-ADC(1)

Public Functions

std::shared_ptr<Tensor> diagonal_s() const

Compute the singles diagonal of the ADC matrix and return it.

void compute_apply_ss(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the singles-singles block to the vector in, returning the result in the vector out.

class CvsAdc2MatrixCore : public adcc::CvsAdcMatrixCoreSinglesDoubles
#include <CvsAdcMatrixCores.hh>

Matrix core for CVS-ADC(2)

Subclassed by adcc::CvsAdc2xMatrixCore

Public Functions

void compute_apply_ss(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the singles-singles block to the vector in, returning the result in the vector out.

void compute_apply_sd(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the singles-doubles block to the vector in, returning the result in the vector out.

Note
Might not be implemented for all matrix cores.

void compute_apply_ds(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the doubles-singles block to the vector in, returning the result in the vector out.

Note
Might not be implemented for all matrix cores.

void compute_apply_dd(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the doubles-doubles block to the vector in, returning the result in the vector out.

Note
Might not be implemented for all matrix cores.

std::shared_ptr<Tensor> diagonal_s() const

Compute the singles diagonal of the ADC matrix and return it.

std::shared_ptr<Tensor> diagonal_d() const

Compute the doubles diagonal of the ADC matrix and return it.

Note
Might not be implemented for all matrix cores.

void compute_matvec(const AmplitudeVector &ins, const AmplitudeVector &outs) const

Compute the application of the ADC matrix to the amplitude vector supplied on ins and return the result in outs.

class CvsAdc2xMatrixCore : public adcc::CvsAdc2MatrixCore
#include <CvsAdcMatrixCores.hh>

Matrix core for CVS-ADC(2)-x

Subclassed by adcc::CvsAdc3MatrixCore

Public Functions

void compute_apply_dd(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the doubles-doubles block to the vector in, returning the result in the vector out.

Note
Might not be implemented for all matrix cores.

std::shared_ptr<Tensor> diagonal_d() const

Compute the doubles diagonal of the ADC matrix and return it.

Note
Might not be implemented for all matrix cores.

void compute_matvec(const AmplitudeVector &ins, const AmplitudeVector &outs) const

Compute the application of the ADC matrix to the amplitude vector supplied on ins and return the result in outs.

class CvsAdc3MatrixCore : public adcc::CvsAdc2xMatrixCore
#include <CvsAdcMatrixCores.hh>

Matrix core for CVS-ADC(3)

Public Functions

void compute_apply_ss(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the singles-singles block to the vector in, returning the result in the vector out.

void compute_apply_sd(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the singles-doubles block to the vector in, returning the result in the vector out.

Note
Might not be implemented for all matrix cores.

void compute_apply_ds(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the doubles-singles block to the vector in, returning the result in the vector out.

Note
Might not be implemented for all matrix cores.

std::shared_ptr<Tensor> diagonal_s() const

Compute the singles diagonal of the ADC matrix and return it.

void compute_matvec(const AmplitudeVector &ins, const AmplitudeVector &outs) const

Compute the application of the ADC matrix to the amplitude vector supplied on ins and return the result in outs.

class GenericAdcMatrixCoreSingles : public adcc::AdcMatrixCoreSingles
#include <GenericAdcMatrixCores.hh>

Base class for generic matrix cores with only a singles block

Subclassed by adcc::Adc0MatrixCore, adcc::Adc1MatrixCore

Public Functions

std::vector<size_t> shape() const

Return the shape of the represented matrix

std::vector<std::string> block_spaces(const std::string &block) const

Return the list of spaces involved for a particular block of the matrix

class GenericAdcMatrixCoreSinglesDoubles : public adcc::AdcMatrixCoreSinglesDoubles
#include <GenericAdcMatrixCores.hh>

Base class for generic matrix cores with a singles and a doubles block

Subclassed by adcc::Adc2MatrixCore

Public Functions

std::vector<size_t> shape() const

Return the shape of the represented matrix

std::vector<std::string> block_spaces(const std::string &block) const

Return the list of spaces involved for a particular block of the matrix

class Adc0MatrixCore : public adcc::GenericAdcMatrixCoreSingles
#include <GenericAdcMatrixCores.hh>

Matrix core for ADC(0)

Public Functions

std::shared_ptr<Tensor> diagonal_s() const

Compute the singles diagonal of the ADC matrix and return it.

void compute_apply_ss(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the singles-singles block to the vector in, returning the result in the vector out.

class Adc1MatrixCore : public adcc::GenericAdcMatrixCoreSingles
#include <GenericAdcMatrixCores.hh>

Matrix core for ADC(1)

Public Functions

std::shared_ptr<Tensor> diagonal_s() const

Compute the singles diagonal of the ADC matrix and return it.

void compute_apply_ss(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the singles-singles block to the vector in, returning the result in the vector out.

class Adc2MatrixCore : public adcc::GenericAdcMatrixCoreSinglesDoubles
#include <GenericAdcMatrixCores.hh>

Matrix core for ADC(2)

Subclassed by adcc::Adc2xMatrixCore

Public Functions

void compute_apply_ss(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the singles-singles block to the vector in, returning the result in the vector out.

void compute_apply_sd(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the singles-doubles block to the vector in, returning the result in the vector out.

Note
Might not be implemented for all matrix cores.

void compute_apply_ds(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the doubles-singles block to the vector in, returning the result in the vector out.

Note
Might not be implemented for all matrix cores.

void compute_apply_dd(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the doubles-doubles block to the vector in, returning the result in the vector out.

Note
Might not be implemented for all matrix cores.

std::shared_ptr<Tensor> diagonal_s() const

Compute the singles diagonal of the ADC matrix and return it.

std::shared_ptr<Tensor> diagonal_d() const

Compute the doubles diagonal of the ADC matrix and return it.

Note
Might not be implemented for all matrix cores.

void compute_matvec(const AmplitudeVector &ins, const AmplitudeVector &outs) const

Compute the application of the ADC matrix to the amplitude vector supplied on ins and return the result in outs.

class Adc2xMatrixCore : public adcc::Adc2MatrixCore
#include <GenericAdcMatrixCores.hh>

Matrix core for ADC(2)-x

Subclassed by adcc::Adc3MatrixCore

Public Functions

void compute_apply_dd(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the doubles-doubles block to the vector in, returning the result in the vector out.

Note
Might not be implemented for all matrix cores.

std::shared_ptr<Tensor> diagonal_d() const

Compute the doubles diagonal of the ADC matrix and return it.

Note
Might not be implemented for all matrix cores.

void compute_matvec(const AmplitudeVector &ins, const AmplitudeVector &outs) const

Compute the application of the ADC matrix to the amplitude vector supplied on ins and return the result in outs.

class Adc3MatrixCore : public adcc::Adc2xMatrixCore
#include <GenericAdcMatrixCores.hh>

Matrix core for ADC(3)

Public Functions

void compute_apply_ss(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the singles-singles block to the vector in, returning the result in the vector out.

void compute_apply_sd(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the singles-doubles block to the vector in, returning the result in the vector out.

Note
Might not be implemented for all matrix cores.

void compute_apply_ds(const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Compute the application of the doubles-singles block to the vector in, returning the result in the vector out.

Note
Might not be implemented for all matrix cores.

std::shared_ptr<Tensor> diagonal_s() const

Compute the singles diagonal of the ADC matrix and return it.

void compute_matvec(const AmplitudeVector &ins, const AmplitudeVector &outs) const

Compute the application of the ADC matrix to the amplitude vector supplied on ins and return the result in outs.

class AdcMatrix
#include <AdcMatrix.hh>

Class representing an ADC matrix

The ADC matrix is tiled into the following regions,

+------------+----------------+
| singles-   |     singles-   |
| singles    |     doubles    |
+------------+----------------+
| doubles-   |     doubles-   |
| singles    |     doubles    |
+------------+----------------+

where for ADC(0) and ADC(1) only the singles-singles part is present. Calling functions, which compute or make use of the doubles part in such cases leads to a std::invalid_argument.

In case no intermediates are provided to the object prior to the first payload function call, the class will automatically generate the intermediates by itself and add them to the AdcIntermediates structure referenced by intermediates_ptr. Such structure will be created in case it does not exist.

Public Functions

AdcMatrix(std::string method, std::shared_ptr<const LazyMp> ground_state_ptr)

Construct an AdcMatrix object.

Parameters
  • method: The ADC method to use
  • ground_state_ptr: The MP ground state to build upon

std::vector<size_t> shape() const

Return the shape of the represented matrix

bool has_block(const std::string &block) const

Return true if the referenced block exists in this matrix

std::vector<std::string> blocks() const

Return the list of blocks, which exist in this matrix

std::vector<std::string> block_spaces(const std::string &block) const

Return the list of spaces involved for a particular block of the matrix

void compute_apply(std::string block, const std::shared_ptr<Tensor> &in, const std::shared_ptr<Tensor> &out) const

Apply a block of the ADC matrix (e.g. the singles-singles block or the singles-doubles block. The result is stored in a pre-allocated tensor.

Parameters
  • block: Block specification, may take the values “ss”, “sd”, “ds”, “dd”. May not take all values for all ADC variants.
  • in: Input tensor, rhs of the expression
  • out: Output tensor, result of the application.

std::shared_ptr<Tensor> diagonal(std::string block) const

Return the diagonal of the matrix

void compute_matvec(const AmplitudeVector &ins, const AmplitudeVector &outs) const

Compute the matrix-vector product of the ADC matrix with an excitation amplitude. The rhs amplitude parts are provided as std::vector, containing singles, doubles, … part and so are resulting amplitudes. The result is stored in a pre-allocated tensor.

Parameters
  • ins: Input tensors, one per amplitude part (singles, doubles)
  • outs: Output tensors as a shared pointers

std::string method() const

Access to the adc method this class represents

bool is_core_valence_separated() const

Is the core-valence separation applied in this matrix

std::shared_ptr<const ReferenceState> reference_state_ptr() const

Access to the reference state on top of which the ADC calculation has been started

std::shared_ptr<const LazyMp> ground_state_ptr() const

Obtain the pointer to the MP ground state upon which this ADC matrix is based.

std::shared_ptr<const MoSpaces> mospaces_ptr() const

Return the MoSpaces to which this LazyMp is set up

std::shared_ptr<AdcIntermediates> get_intermediates_ptr() const

Get the pointer to the AdcIntermediates cache for the precomputed intermediates

void set_intermediates_ptr(std::shared_ptr<AdcIntermediates> im_ptr)

Replace the currently held AdcIntermediates cache object by the supplied one. This allows to reuse intermediates from a previous calculation.

const Timer &timer() const

Obtain the timing info contained in this object

class AmplitudeVector : public std::vector<std::shared_ptr<Tensor>>
#include <AmplitudeVector.hh>

This class is effectively a list of Tensors, which can be thought of as blocks of the ADC amplitude vector

ADC guess setup

void adcc::amplitude_vector_enforce_spin_kind(std::shared_ptr<Tensor> tensor, std::string block, std::string spin_kind)

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 to
  • block: The block of an amplitude this tensor represents
  • spin_kind: The kind of spin to enforce

std::vector<std::shared_ptr<Symmetry>> adcc::guess_symmetries(const AdcMatrix &matrix, const AdcGuessKind &kind)

Return the symmetries for setting up guess vectors. The first entry is the singles symmetry the second the doubles symmetry and so on.

AmplitudeVector adcc::guess_zero(const AdcMatrix &matrix, const AdcGuessKind &kind)

Return an AmplitudeVector object filled with zeros, but where the symmetry has been properly set up to meet the requirements of the AdcGuessKind object

std::vector<AmplitudeVector> adcc::guesses_from_diagonal(const AdcMatrix &matrix, const AdcGuessKind &kind, std::string block, size_t n_guesses, scalar_type degeneracy_tolerance)

Obtain guesses by inspecting a block of the diagonal of the passed ADC matrix. The symmetry of the returned vectors is already set-up properly. Note that this routine may return fewer vectors than requested in case the requested number could not be found.

matrix AdcMatrix for which guesses are to be constructed kind AdcGuessKind object describing the kind of guesses to be constructed block The block of the diagonal to investigate. May be “s” (singles) or “d” (doubles). n_guesses The number of guesses to look for. degeneracy_tolerance Tolerance for two entries of the diagonal to be considered degenerate, i.e. identical.

struct AdcGuessKind
#include <guess_zero.hh>

Struct, which collects information about the kind of guess vectors to be constructed.

To make the collected quantities more relatable, here are common profiles:

  • Singlet RHF reference
    • for computing singlet excited states: spin_change == 0 and spin_block_symmetrisation == “symmetric”
    • for computing triplet excited states: spin change == 0 and spin_block_symmetrisation == “antisymmetric”
  • UHF reference
    • for computing any kind of excited state: spin_change == 0 and spin_symmetrisation == “none”
    • for spin-flip states: spin_change == -1 and spin_symmetrisation == “none”
Note

The information is collected in the sense of “what the guess vectors should

look like” as opposed to “what excitations are computed with such vectors”. The reasoning is that the latter also depends on the reference, which, however, does not influence the way the guess algorithms need to find the guess vectors.

Public Members

std::string irrep

String describing the irreducable representation to consider

float spin_change

The spin change to enforce in an excitation. Typical values are 0 and -1.

int spin_change_twice

Twice the value of spin_change as an integer

std::string spin_block_symmetrisation

Symmetrisation to enforce between equivalent spin blocks, which all yield the desired spin_change. E.g. if spin_change == 0, then both the alpha->alpha and beta->beta blocks of the singles part of the excitation vector achieve this spin change. The symmetry specified with this parameter will then be imposed between the a-a and b-b blocks. Valid values are “none”, “symmetric” and “antisymmetric”, where “none” enforces no symmetry.

ISR and one-particle densities

std::shared_ptr<OneParticleOperator> adcc::compute_gs2state_optdm(std::string method, std::shared_ptr<const LazyMp> ground_state_ptr, const AmplitudeVector &excitation_amplitude, std::shared_ptr<AdcIntermediates> intermediates_ptr)

Compute the ground-state to excited state one-particle transition-density matrix (optdm) at the ADC level of theory given by the method_ptr.

Parameters
  • method: The ADC method to use
  • ground_state_ptr: The MP ground state to build upon
  • excitation_amplitude: The excitation amplitude, given as a std::vector of blocks (singles, doubles, …)
  • intermediates_ptr: Adc intermediates to reuse during property calculation

std::shared_ptr<AmplitudeVector> adcc::compute_modified_transition_moments(std::string method, std::shared_ptr<const LazyMp> ground_state_ptr, std::shared_ptr<OneParticleOperator> op_ptr)

Compute the modified transition moments (MTM)

Parameters
  • method: The ADC method to use
  • ground_state_ptr: The MP ground state to build upon
  • op_ptr: The operator for which the MTM should be computed

std::shared_ptr<OneParticleOperator> adcc::compute_state2state_optdm(std::string method, std::shared_ptr<const LazyMp> ground_state_ptr, const AmplitudeVector &excitation_amplitude_from, const AmplitudeVector &excitation_amplitude_to, std::shared_ptr<AdcIntermediates> intermediates_ptr)

Compute the excited state to excited state one-particle transition-density matrix (optdm) at the ADC level of theory given by the method_ptr.

Parameters
  • method: The ADC method to use
  • ground_state_ptr: The MP ground state to build upon
  • excitation_amplitude_from: The excitation amplitude, given as a std::vector of blocks (singles, doubles, …) from which the electronic transition starts
  • excitation_amplitude_to: The excitation amplitude, given as a std::vector of blocks (singles, doubles, …) where the electronic transition ends.
  • intermediates_ptr: Adc intermediates to reuse during property calculation

std::shared_ptr<OneParticleOperator> adcc::compute_state_diffdm(std::string method, std::shared_ptr<const LazyMp> ground_state_ptr, const AmplitudeVector &excitation_amplitude, std::shared_ptr<AdcIntermediates> intermediates_ptr)

Compute the one-particle density matrix of a particular state at the ADC level of theory given by the method_ptr.

Parameters
  • method: The ADC method to use
  • ground_state_ptr: The MP ground state to build upon
  • excitation_amplitude: The excitation amplitude, given as a std::vector of blocks (singles, doubles, …)
  • intermediates_ptr: Adc intermediates to reuse during property calculation

struct OneParticleOperator
#include <OneParticleOperator.hh>

Container for a one-particle operator in MO representation. This container is used both for state density matrices as well as for transition-density matrices.

This object is a map from a block label (like o1o1) to the actual tensor object.

Public Functions

OneParticleOperator(std::shared_ptr<const MoSpaces> mospaces_ptr, bool is_symmetric = true, std::string cartesian_transform = "1")

Construct an empty OneParticleOperator object. Will construct an object with all blocks initialised as zero blocks.

Parameters
  • mospaces: MoSpaces object
  • is_symmetric: Is the operator symmetric
  • cartesian_transformation: The cartesian function according to which the operator transforms. See make_symmetry_operator for details.

OneParticleOperator(std::shared_ptr<const MoSpaces> mospaces_ptr, ctx::CtxMap ctx)

Construct a OneParticleOperator taking data from a context, like the one returned by the to_ctx function

size_t ndim() const

Number of dimensions

std::vector<size_t> shape() const

Shape of each dimension

size_t size() const

Number of elements

std::shared_ptr<OneParticleOperator> copy() const

Return a deep copy of this OneParticleOperator.

This can be used to conveniently add tensor corrections on top of existing OneParticleOperator objects.

std::shared_ptr<Tensor> block(std::string block) const

Get a particular block of the one-particle operator, throwing an invalid_argument in case the block is a known zero block

std::shared_ptr<Tensor> operator[](std::string label)

Get a particular block of the one-particle operator. It allocates and stores a tensor representing a zero block on the fly if required.

void set_block(std::string block, std::shared_ptr<Tensor> value)

Set a particular block of the one-particle operator

void set_zero_block(std::vector<std::string> blocks)

Set the list of specified blocks to zero without storing an explicit zero tensor for them. If previously a tensor has been stored for them it will be purged.

void set_zero_block(std::string block)

The version of set_zero_block for only a single block

bool is_zero_block(std::string block) const

Return whether the specified block is a block flagged to be exactly zero using the set_zero_block function.

bool has_block(std::string block) const

Check weather a particular block exists

std::vector<std::string> blocks() const

Return the list of all registered blocks

std::vector<std::string> blocks_nonzero() const

Return the list of all registered not-zero blocks

std::vector<std::string> orbital_subspaces() const

Return the list of all orbital subspaces contained in the operator

bool is_symmetric() const

Is the operator represented by this object symmetric If true then certain density matrices (like v1o1) might be missing, since they are equivalent to already existing blocks (such as o1v1 in our example)

std::string cartesian_transform() const

The cartesian function according to which the operator transforms. See make_symmetry_operator for details.

std::pair<std::shared_ptr<Tensor>, std::shared_ptr<Tensor>> to_ao_basis(std::shared_ptr<const ReferenceState> reference_state_ptr) const

Transform this operator object to the AO basis, which was used to obtain the provided reference state.

Return
A pair of matrices in the AO basis. The first is the matrix for alpha spin, the second for beta spin.

std::pair<std::shared_ptr<Tensor>, std::shared_ptr<Tensor>> to_ao_basis(std::map<std::string, std::shared_ptr<Tensor>> coefficient_ptrs) const

Transform this operator object to the AO basis using the coefficient tensors provided for each molecular orbital subspace. The coefficient tensors should be provided via a mapping from the subspace identifier to the actual coefficient for each spin part. E.g. “o1_a” should map to the coefficient “c_o1b” which are the alpha coefficients for the first occupied orbital subspace.

Return
A pair of matrices in the AO basis. The first is the matrix for alpha spin, the second for beta spin.

void export_to(scalar_type *memptr, size_t size) const

Extract the operator 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

ctx::CtxMap to_ctx() const

return the data represented in this object inside a ctx::CtxMap

Tensor interface

The generalised adcc::Tensor interface used by adcc and adccore to perform tensor operations.

std::shared_ptr<Symmetry> adcc::make_symmetry_orbital_energies(std::shared_ptr<const MoSpaces> mospaces_ptr, const std::string &space)

Return the Symmetry object for the orbital energies tensor for the passed orbital subspace.

std::shared_ptr<Symmetry> adcc::make_symmetry_orbital_coefficients(std::shared_ptr<const MoSpaces> mospaces_ptr, const std::string &space, size_t n_bas, const std::string &blocks = "ab")

Return the Symmetry object for the orbital coefficient tensor for the passed orbital subspace.

Parameters
  • mospaces_ptr: MoSpaces pointer
  • space: Space to use (e.g. fb or o1b)
  • n_bas: Number of AO basis functions
  • blocks: 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).

std::shared_ptr<Symmetry> adcc::make_symmetry_eri(std::shared_ptr<const MoSpaces> mospaces_ptr, const std::string &space)

Return the Symmetry object for the electron-repulsion integral tensor for the passed orbital subspace.

std::shared_ptr<Symmetry> adcc::make_symmetry_operator(std::shared_ptr<const MoSpaces> mospaces_ptr, const std::string &space, bool symmetric, const std::string &cartesian_transformation)

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 pointer
  • space: 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

std::shared_ptr<Symmetry> adcc::make_symmetry_operator_basis(std::shared_ptr<const MoSpaces> mospaces_ptr, size_t n_bas, bool symmetric)

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 pointer
  • n_bas: Number of AO basis functions
  • symmetric: Is the tensor symmetric (only in effect if both space axes identical). false disables a setup of permutational symmetry.

std::ostream &adcc::operator<<(std::ostream &out, const Tensor &tensor)

Write a string representation of the Tensor to the provided stream

std::shared_ptr<Tensor> adcc::contract(std::string contraction, std::shared_ptr<Tensor> A, std::shared_ptr<Tensor> B)

Contract tensor A with another tensor and return the result. Let this tensor be A_{abcd} and we contract it with B_{cdef} over the indices c and d to form out_{abef}. This would be achieved using A.contract(“abcd,cdef->abef”, B)

std::shared_ptr<Tensor> adcc::make_tensor_zero(std::shared_ptr<Symmetry> symmetry)

Construct a tensor initialised to zero using a Symmetry object

std::shared_ptr<Tensor> adcc::make_tensor_empty(std::shared_ptr<Symmetry> symmetry)

Construct an empty (uninitialised) tensor using a Symmetry object

struct ContractionIndices
#include <ContractionIndices.hh>

A class managing the indices involved in a contraction.

Public Functions

std::pair<size_t, size_t> contracted_tensor_axes(const char &c) const

Return the indices into the axes of first and second which are contracted by c

Exceptions
  • invalid_argument: If c is in result

bool is_contraction_index(const char &c) const

Returns whether the passed index is contracted over

bool is_result_index(const char &c) const

Returns whether the passed index is part of the result

Public Members

std::string first

The indices on the first input tensor

std::string second

The indices on the second input tensor

std::string result

The indices on the output tensor

std::string contraction

The indices over which the contraction is performed

std::string trace

Trace indices (Indices which repeat in a single tensor)

std::string all

All indices (contracted or free) which are involved in this contraction in alphabetical order.

struct AxisInfo
#include <Symmetry.hh>

Info data structure for extra axes (i.e. AO axes)

Public Functions

bool has_spin() const

Does the axis have a spin-splitting, or can no spin be identified

Public Members

size_t n_orbs_alpha

The number of (alpha) orbitals.

size_t n_orbs_beta

The number of beta orbitals (or zero for axes without spin)

std::vector<size_t> block_starts

The split indices when a new block starts.

class Symmetry
#include <Symmetry.hh>

Container for Tensor symmetry

Public Functions

Symmetry(std::shared_ptr<const MoSpaces> mospaces_ptr, const std::string &space)

Create an empty symmetry object for a Tensor with the provided space string

Symmetry(std::shared_ptr<const MoSpaces> mospaces_ptr, const std::string &space, std::map<std::string, std::pair<size_t, size_t>> extra_axes_orbs)

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

std::shared_ptr<const MoSpaces> mospaces_ptr() const

Return the MoSpaces to which this Symmetry is set up

std::shared_ptr<const AdcMemory> adcmem_ptr() const

Return the adc memory keepalive object

const std::map<std::string, AxisInfo> &extra_axes() const

Return the map for supplying the number of alpha and beta orbitals for the additional axes

Note
This is an internal function and can go at any time.

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 Tensor
#include <Tensor.hh>

The tensor interface used by adccore

Subclassed by adcc::TensorImpl< N >

Interface of a tensor, which is exposed to python.

virtual size_t ndim() const = 0

Number of dimensions

virtual std::vector<size_t> shape() const = 0

Shape of each dimension

virtual size_t size() const = 0

Number of elements

virtual std::shared_ptr<Tensor> empty_like() const = 0

Space to which the tensor is initialised Return a Symmetry object representing the Symmetry of the tensor

Note
This represents a copy of the symmetry, i.e. the tensor cannot be modified using this function.Return a new tensor with the same dimensionality and symmetry as the passed tensor. The elements have undefined values.

virtual std::shared_ptr<Tensor> zeros_like() const = 0

Return a new tensor with the same dimensionality and symmetry as the passed tensor and all elements set to zero.

virtual std::shared_ptr<Tensor> ones_like() const = 0

Return a new tensor with the same dimensionality and symmetry as the passed tensor and all elements set to one.

virtual 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

virtual std::shared_ptr<Tensor> copy() const

Return a deep copy of this tensor

virtual void copy_to(std::shared_ptr<Tensor> other) const = 0

Export a deep copy of this tensor to another tensor

virtual std::shared_ptr<Tensor> transpose(std::vector<size_t> axes) const = 0

Return a transposed form of this tensor as a copy

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.

virtual std::shared_ptr<Tensor> transpose() const

Return the transposed form of a matrix as a copy

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

virtual void scale(scalar_type c) = 0

In-place scale by a scalar value.

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

virtual void set_random() = 0

Set the tensor to random data preserving symmetry

virtual 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

virtual scalar_type dot(std::shared_ptr<Tensor> other) const

Compute the Frobenius or l2 inner product with another tensor

virtual void add(std::vector<scalar_type> scalars, std::vector<std::shared_ptr<Tensor>> tensors) = 0

Add a linear combination of tensors to this tensor

scalars List of scalar prefactors for the tensors tensors Tensors to add in

Public Functions

Tensor(std::shared_ptr<const AdcMemory> adcmem_ptr)

Construct a tensor interface class.

adcmem_ptr Memory keepalive object pointer

virtual void multiply_to(std::shared_ptr<Tensor> other, std::shared_ptr<Tensor> out) const = 0

Multiply another tensor to this tensor, that is form the expression out = this * other

virtual void divide_to(std::shared_ptr<Tensor> other, std::shared_ptr<Tensor> out) const = 0

Divide by another tensor, that is form the expression out = this / other

virtual void contract_to(std::string contraction, std::shared_ptr<Tensor> other, std::shared_ptr<Tensor> out) const = 0

Contract this tensor with another tensor and save the result in out. Let this tensor be A_{abcd} and we contract it with B_{cdef} over the indices c and d to form out_{abef}. This would be achieved using A.contract_to(“abcd,cdef->abef”, B, out)

virtual std::shared_ptr<Tensor> contract(std::string contraction, std::shared_ptr<Tensor> other) const = 0

Contract this tensor with another tensor and return the result. Let this tensor be A_{abcd} and we contract it with B_{cdef} over the indices c and d to form out_{abef}. This would be achieved using A.contract(“abcd,cdef->abef”, B)

virtual void symmetrise_to(std::shared_ptr<Tensor> out, 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
  • out: The tensor to store the result in.

Note
Unlike the implementation in libtensor, this function includes the prefactor 0.5

virtual void antisymmetrise_to(std::shared_ptr<Tensor> out, const std::vector<std::vector<size_t>> &permutations) const = 0

Antisymmetrise with respect to the given index permutations and export the result to the output tensor. For details with respect to the format of the permutations, see symmetrise_to. 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

virtual void set_immutable() = 0

Flag the tensor as immutable, allows some optimisations to be performed

virtual bool is_mutable() const = 0

Is the tensor mutable

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

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

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

virtual 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 select
  • unique_by_symmetry: By default the returned elements are made unique by symmetry of the tensor. Setting this to false disables this feature.

virtual 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 select
  • unique_by_symmetry: By default the returned elements are made unique by symmetry of the tensor. Setting this to false disables this feature.

virtual 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 select
  • unique_by_symmetry: By default the returned elements are made unique by symmetry of the tensor. Setting this to false disables this feature.

virtual 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 select
  • unique_by_symmetry: By default the returned elements are made unique by symmetry of the tensor. Setting this to false disables this feature.

virtual void print(std::ostream &o) const = 0

Print to an output stream

virtual 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

virtual 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

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

virtual 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 vector
  • 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.

virtual void import_from() = 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.

virtual std::string describe_symmetry() const = 0

Return a std::string providing hopefully helpful information about the symmetries stored inside the tensor object.

std::shared_ptr<const AdcMemory> adcmem_ptr() const

Return a pointer to the memory keep-alive object

Utilities

typedef double adcc::real_type

The real type to use

typedef double adcc::scalar_type

The scalar type to use

MethodParsed adcc::parse_method(std::string method)

Parse a method string and return a MethodParsed object containing some useful information in the form of flat data.

template<typename RandomAccessIterator, typename Compare>
std::vector<size_t> adcc::argsort(const RandomAccessIterator first, const RandomAccessIterator last, Compare cmp)

Argsort returns the index sequence, which would sort the range of elements provided.

Performs an indirect sorting of the range provided by the two iterators. The vectors of indices returned gives the ordering of the iterator range, which would give a sorted range. By the default less is used as the comparator.

E.g. It the range runs over the elements { “d”, “a”, “c”, “b” }, then the returned indices would be { 1, 3, 2, 0 } (if less is used as the comparator).

template<typename RandomAccessIterator>
std::vector<size_t> adcc::argsort(const RandomAccessIterator first, const RandomAccessIterator last)

Sort a range using the argsort function with the default comparator (less).

std::string adcc::shape_to_string(const std::vector<size_t> &shape)

Convert a std::vector<size_t> representing a tensor shape to a string

class 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

AdcMemory()

Setup the AdcMemory class using the std::allocator for memory management

~AdcMemory()

Destroy the AdcMemory class and free all the allocated memory

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 and max_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 files
  • max_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 MethodParsed
#include <parse_method.hh>

Parsed information about an ADC method

Public Members

bool extended

For ADC level 2 only: ADC(2) or ADC(2)-x

bool is_core_valence_separated

Should the CVS approximation be applied

class ThreadPool
#include <ThreadPool.hh>

Pool managing how many threads may be used for calculations.

Public Functions

ThreadPool(size_t n_cores, size_t n_threads)

Initialise the thread pool.

Parameters
  • n_cores: The number of cpu cores to employ
  • n_threads: The number of threads to use (default 2*n_cores - 1)

void reinit(size_t n_cores, size_t n_threads)

Reinitialise the thread pool.

Parameters
  • n_cores: The number of cpu cores to employ
  • n_threads: The number of threads to use

ThreadPool()

Initialise a thread pool without parallelisation

size_t n_cores() const

Return the number of cores employed for parallelisation.

size_t n_threads() const

Return the number of threads employed for parallelisation.

class 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

static 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 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

RecordTime(Timer &timer_, const std::string &task_)

Construct the class and start the timer

Public Members

Timer &timer

Timer object, which is managed.

std::string task

String describing the task.

Metadata access

These classes and functions provide access to metadate about adccore.

std::vector<std::string> adcc::__features__()

Return the list of compiled-in features

std::string adcc::__authors__()

Return the authors string

std::string adcc::__email__()

Return the email string

std::vector<Component> adcc::__components__()

Get the list of components compiled into adccore

struct version
#include <metadata.hh>

Version information

Public Static Functions

static int major_part()

Return the major part of the version

static int minor_part()

Return the minor part of the version

static int patch_part()

Return the patch part of the version

static bool is_debug()

Is the compiled version a Debug version

static std::string version_string()

Return the version as a string

struct Component
#include <metadata.hh>

Data about a library or third-party component

Public Members

std::string name

The name of the component.

std::string version

The version string.

std::string authors

The authors.

std::string description

A brief description.

std::string doi

DOI to a publication.

std::string website

Website of the upstream source.

std::string licence

A short identifier of the licence.

Tensor implementation using libtensor

This section describes the implementation of the Tensor functionality of adcc::Tensor using the libtensor tensor library.

libtensor::btensor<1, scalar_type> &adcc::asbt1(const std::shared_ptr<Tensor> &tensor)
libtensor::btensor<2, scalar_type> &adcc::asbt2(const std::shared_ptr<Tensor> &tensor)
libtensor::btensor<3, scalar_type> &adcc::asbt3(const std::shared_ptr<Tensor> &tensor)
libtensor::btensor<4, scalar_type> &adcc::asbt4(const std::shared_ptr<Tensor> &tensor)
IdedBispace<1> adcc::make_bispace(const std::string &space, std::vector<size_t> block_starts, size_t length)

Make a libtensor bispace from some block starting indices and a space string describing the space for which the bispace object should be made

std::vector<size_t> adcc::make_block_starts(const libtensor::block_index_space<1> &bis)

Extract the block starts from a libtensor block-index space

template<size_t N>
IdedBispace<2 * N> adcc::bispace_product(const IdedBispace<N> &x, const IdedBispace<N> &y)

Form higher spaces by taking the tensor products of 2 same-sized bispaces

This function takes the product between two IdedBispace objects x and y by appending the string identifiers and taking the tensor product between the bispace<N> objects.

std::map<size_t, std::string> adcc::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>
std::shared_ptr<libtensor::symmetry<N, scalar_type>> adcc::as_lt_symmetry(const Symmetry &sym)

Translate a adcc::Symmetry object into a shared pointer of the libtensor symmetry object

template<size_t DIM>
lt::permutation<DIM> adcc::construct_contraction_output_permutation(const ContractionIndices &parsed_contraction)

Construct the output tensor permutation libtensor needs to dump the output tensor data of a contraction in the right dimension order

template<size_t N>
std::shared_ptr<Tensor> adcc::wrap_libtensor(std::shared_ptr<const AdcMemory> adcmem_ptr, std::shared_ptr<libtensor::btensor<N, scalar_type>> libtensor_ptr)

Enwrap an existing libtensor object to make an adcc::Tensor

template<size_t N>
libtensor::btensor<N, scalar_type> &adcc::as_btensor(const std::shared_ptr<Tensor> &tensor)

Extractor functions to convert the contained tensor of an adcc::Tensor object to an appropriate libtensor::btensor, The dimensionality N has to match

template<size_t N>
std::shared_ptr<libtensor::btensor<N, scalar_type>> adcc::as_btensor_ptr(const std::shared_ptr<Tensor> &tensor)

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>
struct IdedBispace : public libtensor::bispace<N>
#include <bispace.hh>

Struct to logically combine a libtensor::bispace<N> and a string array of size N, which contains the ids of the bispace<1> objects which make up the bispace<N>.

Public Functions

std::string spaces_string() const

The spaces as a concatenated string

Public Members

std::array<std::string, N> spaces

The elementary spaces

template<size_t N>
class TensorImpl : public adcc::Tensor
#include <TensorImpl.hh>

Actual implementation of the Tensor class. For details about the functions see Tensor.hh

Public Functions

TensorImpl(std::shared_ptr<const AdcMemory> adcmem_ptr, std::shared_ptr<libtensor::btensor<N, scalar_type>> libtensor_ptr)

Construction of the tensor implementation class

Note
This constructor is for experts only. Non-experts should use the TensorFactory instead.
Note
We had an constructor from std::shared_ptr<libtensor::any_tensor> here, which does not work, because any_tensor is only a reference to some underlying tensor object, but indicates no ownership. This completely clashes with the idea of this class causing it to fail at random places.

size_t ndim() const

Number of dimensions

std::vector<size_t> shape() const

Shape of each dimension

size_t size() const

Number of elements

std::shared_ptr<Tensor> nosym_like() const

Return a new tensor with same dimensionality and shape (and bispace) but without any symmetry copied over

std::shared_ptr<Tensor> empty_like() const

Space to which the tensor is initialised Return a Symmetry object representing the Symmetry of the tensor

Note
This represents a copy of the symmetry, i.e. the tensor cannot be modified using this function.Return a new tensor with the same dimensionality and symmetry as the passed tensor. The elements have undefined values.

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 copy_to(std::shared_ptr<Tensor> other) const

Export a deep copy of this tensor to another tensor

std::shared_ptr<Tensor> transpose(std::vector<size_t> axes) const

Return a transposed form of this tensor as a copy

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.

double trace(std::string contraction) const

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 scale(scalar_type c)

In-place scale by a scalar value.

void set_mask(std::string mask, scalar_type value)

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 set_random()

Set the tensor to random data preserving symmetry

void add(std::vector<scalar_type> scalars, std::vector<std::shared_ptr<Tensor>> tensors)

Add a linear combination of tensors to this tensor

scalars List of scalar prefactors for the tensors tensors Tensors to add in

void multiply_to(std::shared_ptr<Tensor> other, std::shared_ptr<Tensor> out) const

Multiply another tensor to this tensor, that is form the expression out = this * other

void divide_to(std::shared_ptr<Tensor> other, std::shared_ptr<Tensor> out) const

Divide by another tensor, that is form the expression out = this / other

void contract_to(std::string contraction, std::shared_ptr<Tensor> other, std::shared_ptr<Tensor> out) const

Contract this tensor with another tensor and save the result in out. Let this tensor be A_{abcd} and we contract it with B_{cdef} over the indices c and d to form out_{abef}. This would be achieved using A.contract_to(“abcd,cdef->abef”, B, out)

std::shared_ptr<Tensor> contract(std::string contraction, std::shared_ptr<Tensor> other) const

Contract this tensor with another tensor and return the result. Let this tensor be A_{abcd} and we contract it with B_{cdef} over the indices c and d to form out_{abef}. This would be achieved using A.contract(“abcd,cdef->abef”, B)

void symmetrise_to(std::shared_ptr<Tensor> out, const std::vector<std::vector<size_t>> &permutations) const

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
  • out: The tensor to store the result in.

Note
Unlike the implementation in libtensor, this function includes the prefactor 0.5

void antisymmetrise_to(std::shared_ptr<Tensor> out, const std::vector<std::vector<size_t>> &permutations) const

Antisymmetrise with respect to the given index permutations and export the result to the output tensor. For details with respect to the format of the permutations, see symmetrise_to. 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

std::vector<scalar_type> dot(std::vector<std::shared_ptr<Tensor>> tensors) const

Compute the Frobenius or l2 inner product with a list of other tensors, returning the results as a scalar array

void set_immutable()

Flag the tensor as immutable, allows some optimisations to be performed

bool is_mutable() const

Is the tensor mutable

void set_element(const std::vector<size_t> &idx, scalar_type value)

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

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

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

Get the n absolute largest elements along with their values

Parameters
  • n: Number of elements to select
  • unique_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

Get the n absolute smallest elements along with their values

Parameters
  • n: Number of elements to select
  • unique_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

Get the n largest elements along with their values

Parameters
  • n: Number of elements to select
  • unique_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

Get the n smallest elements along with their values

Parameters
  • n: Number of elements to select
  • unique_by_symmetry: By default the returned elements are made unique by symmetry of the tensor. Setting this to false disables this feature.

void print(std::ostream &o) const

Print to an output stream

void export_to(scalar_type *memptr, size_t size) const

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)

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.

virtual void import_from()

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

Return a std::string providing hopefully helpful information about the symmetries stored inside the tensor object.

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.