3.6.4. horton.matrix.dense
– Dense matrix implementations¶
3.6.4.1. Naming scheme for the contract, slice or expand methods¶
The name of contract
, slice
and expand
methods is as follows:
[iadd_]{contract,slice,expand}[_X][_Y][_to_Z]
where each part between square brackets is optional. X
, Y
and Z
can be any of one
, two
, three
or four
. The names
contract
, slice
and expand
have the following meaning:
contract
- Some actual contraction takes place, i.e. summing of (products) of matrix elements.
slice
- A subset of the elements is merely selected, without any addition of
elements. This is never combined with
iadd_
. expand
- Products of elements are computed but these products are not added. Similar to an outer product but more general.
When iadd_
is used as a prefix, the result of the contraction is added
in-place to the self object. In this case, the _to_Z
part is never
present. A contraction of all input arguments is made. The dimensionality
of the input arguments is indicated by _X
and _Y
.
When _to_Z
is present, the contraction involves self and possibly other
arguments whose dimensionality is indicated by _X
and _Y
. In this
case, iadd_
can not be present. The result is written to an output
argument. If the output argument is not provided, fresh memory is allocated
to compute the contraction and the result is returned. (This allocation is
provided for convenience but should not be used in critical situations.)
Some remarks:
- Similar conventions apply to an
expand
method. - All
contract
andexpand
methods are implemented with the driver methodDenseLinalgFactory.einsum
. However, other implementations than Dense are free to implement things differently. - All
contract
andexpand
methods never touch the internals of higher-index objects.
For more specific information, read the documentation of the individual classes and methods below.
3.6.4.2. Handling of index symmetry¶
The dense matrix classes do not exploit matrix symmetry to reduce memory
needs. Instead they will happily store non-symmetric data if need be. There
are however a few methods in the DenseTwoIndex
and
DenseFourIndex
classes below that take a symmetry
argument to
check or enforce a certain index symmetry.
The symmetry argument is always an integer that corresponds to the redundancy of the off-diagonal matrix elements in the dense storage. In practice this means the following:
DenseTwoIndex
symmetry=1
: Nothing is checked/enforcedsymmetry=2
: Hermitian index symmetry is checked/enforced (default), i.e. \(\langle i \vert A \vert j \rangle = ` :math:\)langle j vert A vert i rangle`
DenseFourIndex
symmetry=1
: Nothing is checked/enforcedsymmetry=2
: Dummy index symmetry is checked/enforced, i.e. \(\langle ij \vert B \vert kl \rangle =\) \(\langle ji \vert B \vert lk \rangle\)symmetry=4
: Hermitian and real index symmetry are checked/enforced, i.e. \(\langle ij \vert B \vert kl \rangle =\) \(\langle kl \vert B \vert ij \rangle =\) \(\langle kj \vert B \vert il \rangle =\) \(\langle il \vert B \vert kj \rangle\). (This only makes sense because the basis functions are assumed to be real.)symmetry=8
: All possible symmetries are checked/enforced, i.e. \(\langle ij \vert B \vert kl \rangle =\) \(\langle kl \vert B \vert ij \rangle =\) \(\langle kj \vert B \vert il \rangle =\) \(\langle il \vert B \vert kj \rangle =\) \(\langle ji \vert B \vert lk \rangle =\) \(\langle lk \vert B \vert ji \rangle =\) \(\langle jk \vert B \vert li \rangle =\) \(\langle li \vert B \vert jk \rangle\). (This only makes sense because the basis functions are assumed to be real.)
3.6.4.3. Dense matrix classes¶
-
class
horton.matrix.dense.
DenseLinalgFactory
(default_nbasis=None)¶ Bases:
horton.matrix.base.LinalgFactory
Optional arguments:
- default_nbasis
- The default basis size when constructing new operators/expansions.
-
__init__
(default_nbasis=None)¶ Optional arguments:
- default_nbasis
- The default basis size when constructing new operators/expansions.
-
create_expansion
(nbasis=None, nfn=None)¶ Create a DenseExpansion with defaults from the LinalgFactory
Optional arguments:
- nbasis
- The number of basis functions. When not given, the default_nbasis value of the DenseLinalgFactory instance will be used.
- nfn
- The number of orbitals. When not given, the default_nbasis value of the DenseLinalgFactory instance will be used.
-
create_four_index
(nbasis=None, nbasis1=None, nbasis2=None, nbasis3=None)¶ Create a DenseFourIndex with defaults from the LinalgFactory
Optional arguments:
- nbasis
- The number of basis functions. When not given, the default_nbasis value of the DenseLinalgFactory instance will be used.
-
create_one_index
(nbasis=None)¶ Create a DenseOneIndex with defaults from the LinalgFactory
Optional arguments:
- nbasis
- The number of basis functions. When not given, the default_nbasis value of the DenseLinalgFactory instance will be used.
-
create_three_index
(nbasis=None, nbasis1=None, nbasis2=None)¶ Create a DenseThreeIndex with defaults from the LinalgFactory
Optional arguments:
- nbasis
- The number of basis functions. When not given, the default_nbasis value of the DenseLinalgFactory instance will be used.
-
create_two_index
(nbasis=None, nbasis1=None)¶ Create a DenseTwoIndex with defaults from the LinalgFactory
Optional arguments:
- nbasis
- The number of basis functions. When not given, the default_nbasis value of the DenseLinalgFactory instance will be used.
- nbasis1
- The number of basis functions for the second axis if it differes
from
nbasis
.
-
static
einsum
(subscripts, out, factor, clear, *operands)¶ einsum wrapper for DenseNBody objects
Arguments:
- subscripts
- The usual subscripts argument for the einsum function, except
that the format is slightly limited. All indexes must be
explicitly stated. If the output argument is a scalar,
->
may not be present. In all other cases,->
must be present. - out
- An output DenseNBody object. This argument is not allowed when the result is a scalar.
- factor
- Multiply the contraction by a scalar factor.
- clear
- When set to False, the output argument is not zeroed before adding the contraction.
- operands
- A list where each element is (i) a DenseNBody object or (ii) a two-tuple of DenseNBody object and ranges to select for slicing. The ranges must be a tuple of integers of the form (begin0, end0, begin1, end1, …) where the number of (begin, end) pairs depends on the dimensionality of the corresponding tensor.
Returns: the out object. When the out argument is given, this is returned with in-place modifications.
-
from_hdf5
(grp)¶ Construct an instance from data previously stored in an h5py.Group.
Arguments:
- grp
- An h5py.Group object.
-
set_default_nbasis
(nbasis)¶
-
static
tensordot
(a, b, axes, out=None, factor=1.0, clear=True)¶ tensordot wrapper for dense n-index objects.
Arguments:
- a, b, axes
- See documentation of numpy.tensordot
Optional arguments:
- out
- The output argument, an instance of one of the Dense?Index classes.
-
to_hdf5
(grp)¶ Write a LinalgFactory to an HDF5 group
Argument:
- grp
- A h5py.Group instance to write to.
-
class
horton.matrix.dense.
DenseOneIndex
(nbasis)¶ Bases:
horton.matrix.base.OneIndex
Arguments:
- nbasis
- The number of basis functions.
-
__init__
(nbasis)¶ Arguments:
- nbasis
- The number of basis functions.
-
assign
(other, begin0=0, end0=None)¶ Assign with the contents of another object
Arguments:
- other
- Another DenseOneIndex object or a scalar.
Optional arguments:
- begin0, end0
- If specified, only a sublock of self is assigned.
-
change_basis_signs
(signs)¶ Correct for different sign conventions of the basis functions.
Arguments:
- signs
- A numpy array with sign changes indicated by +1 and -1.
-
clear
()¶ Reset all elements to zero.
-
copy
(begin=0, end=None)¶ Return a copy of (a part of) the object
Optional arguments:
- begin, end
- Can be used to select a subblock of the object. When not given, the full range is used.
-
divide
(other, factor=1.0, out=None)¶ Divide two DenseOneIndex object, multiplied by factor, and return output
Arguments:
- other
- A DenseOneIndex instance to be divided
Optional arguments:
- factor
- A scalar factor
- out
- The output DenseOneIndex
-
dot
(other, factor=1.0)¶ Dot product with other DenseOneIndex object, multiplied by factor
Arguments:
- other
- A DenseOneIndex instance to be added
Optional arguments:
- factor
- A scalar factor
-
classmethod
from_hdf5
(grp)¶ Construct an instance from data previously stored in an h5py.Group.
Arguments:
- grp
- An h5py.Group object.
-
get_element
(i)¶ Return a matrix element
-
get_max
(abs=True)¶ Return maximum (absolute) element
-
iadd
(other, factor=1.0)¶ Add another DenseOneIndex object in-place, multiplied by factor
Arguments:
- other
- A DenseOneIndex instance to be added
Optional arguments:
- factor
- A scalar factor
-
iscale
(factor)¶ In-place multiplication with a scalar
Arguments:
- factor
- A scalar factor.
-
mult
(other, out=None, factor=1.0)¶ Muliply with other DenseOneIndex object, multiplied by factor
Arguments:
- other
- A DenseOneIndex instance to be added
Optional arguments:
- out
- The output argument (DenseOneIndex with proper size).
- factor
- A scalar factor
-
new
()¶ Return a new one-index object with the same nbasis
-
norm
()¶ Calculate L2 norm
-
permute_basis
(permutation)¶ Reorder the coefficients for a given permutation of basis functions.
Arguments:
- permutation
- An integer numpy array that defines the new order of the basis functions.
-
randomize
()¶ Fill with random normal data
-
set_element
(i, value)¶ Set a matrix element
-
sort_indices
(reverse=False)¶ Return indices of sorted arguments in decreasing order
Optional arguements
- reverse
- If True search order is reversed to increasing arguements
-
to_hdf5
(grp)¶ Dump this object in an h5py.Group
Arguments:
- grp
- An h5py.Group object.
-
trace
(begin0=0, end0=None)¶ Calculate trace
Optional arguments:
- begin0, end0
- Can be used to select a subblock of the object. When not given, the full range is used.
-
nbasis
¶ The number of basis functions
-
ndim
¶ The number of axes in the N-index object.
-
shape
¶ The shape of the object
-
class
horton.matrix.dense.
DenseExpansion
(nbasis, nfn=None)¶ Bases:
horton.matrix.base.Expansion
Arguments:
- nbasis
- The number of basis functions.
Optional arguments:
- nfn
- The number of functions to store. Defaults to nbasis.
- do_energies
- Also allocate an array to store an energy corresponding to each function.
-
__init__
(nbasis, nfn=None)¶ Arguments:
- nbasis
- The number of basis functions.
Optional arguments:
- nfn
- The number of functions to store. Defaults to nbasis.
- do_energies
- Also allocate an array to store an energy corresponding to each function.
-
assign
(other)¶ Assign with the contents of another object
Arguments:
- other
- Another DenseExpansion or DenseTwoIndex object. If DenseTwoIndex, energies and occupations are set to zero.
-
assign_dot
(left, right)¶ Dot product of orbitals in a DenseExpansion and TwoIndex object
Arguments:
- left, right:
- An expansion and a two-index object, or a two-index and an expansion object.
The transformed orbitals are stored in self.
-
assign_occupations
(occupation)¶ Assign orbital occupations
Arguments:
- occupation
- The orbital occupations to be updated. An OneIndex instance
-
change_basis_signs
(signs)¶ Correct for different sign conventions of the basis functions.
Arguments:
- signs
- A numpy array with sign changes indicated by +1 and -1.
-
check_normalization
(overlap, eps=0.0001)¶ Check that the occupied orbitals are normalized.
Arguments:
- overlap
- The overlap two-index operator
Optional arguments:
- eps
- The allowed deviation from unity, very loose by default.
-
check_orthonormality
(overlap, eps=0.0001)¶ Check that the occupied orbitals are orthogonal and normalized.
Arguments:
- overlap
- The overlap two-index operator
Optional arguments:
- eps
- The allowed deviation from unity, very loose by default.
-
clear
()¶ Reset all elements to zero.
-
copy
()¶ Return a copy of the object
-
derive_naturals
(dm, overlap)¶ Arguments:
- dm
- A DenseTwoIndex object with the density matrix
- overlap
- A DenseTwoIndex object with the overlap matrix
Optional arguments:
-
error_eigen
(fock, overlap)¶ Compute the error of the orbitals with respect to the eigenproblem
Arguments:
- fock
- A DenseTwoIndex Hamiltonian (or Fock) operator.
- overlap
- A DenseTwoIndex overlap operator.
Returns: the RMSD error on the orbital energies
-
from_fock
(fock, overlap)¶ Diagonalize a Fock matrix to obtain orbitals and energies
This method updated the attributes
coeffs
andenergies
in-place.Arguments:
- fock
- The fock matrix, an instance of DenseTwoIndex.
- overlap
- The overlap matrix, an instance of DenseTwoIndex.
-
from_fock_and_dm
(fock, dm, overlap, epstol=1e-08)¶ Combined Diagonalization of a Fock and a density matrix
This routine first diagonalizes the Fock matrix to obtain orbitals and orbital energies. Then, using first order (degenerate) perturbation theory, the occupation numbers are computed and, if needed, the the degeneracies of the Fock orbitals are lifted. It is assumed that the Fock and the density matrices commute. This method updated the attributes
coeffs
,energies
andoccupations
in-place.Arguments:
- fock
- The fock matrix, an instance of DenseTwoIndex.
- dm
- The density matrix, an instance of DenseTwoIndex.
- overlap
- The overlap matrix, an instance of DenseTwoIndex.
Optional arguments:
- epstol
- The threshold for recognizing degenerate energy levels. When two
subsequent energy levels are separated by an orbital energy less
than
epstol
, they are considered to be degenerate. When a series of energy levels have an orbital energy spacing between subsequent levels that is smaller thanepstol
, they are all considered to be part of the same degenerate group. For every degenerate set of orbitals, the density matrix is used to (try to) lift the degeneracy.
-
classmethod
from_hdf5
(grp)¶ Construct an instance from data previously stored in an h5py.Group.
Arguments:
- grp
- An h5py.Group object.
-
get_homo_energy
(offset=0)¶ Return a homo energy
Optional arguments:
- offset
- By default, the (highest) homo energy is returned. When this index is above zero, the corresponding lower homo energy is returned.
-
get_homo_index
(offset=0)¶ Return the index of a HOMO orbital.
Optional arguments:
- offset
- By default, the (highest) homo energy is returned. When this index is above zero, the corresponding lower homo energy is returned.
-
get_lumo_energy
(offset=0)¶ Return a lumo energy
Optional arguments:
- offset
- By default, the (lowest) lumo energy is returned. When this index is above zero, the corresponding higher homo energy is returned.
-
get_lumo_index
(offset=0)¶ Return the index of a LUMO orbital.
Optional arguments:
- offset
- By default, the (lowest) lumo energy is returned. When this index is above zero, the corresponding higher homo energy is returned.
-
imul
(other)¶ Inplace multiplication with other DenseOneIndex.
The attributes
energies
andoccupations
are not altered.Arguments:
- other
- A DenseOneIndex object.
-
itranspose
()¶ In-place transpose
-
new
()¶ Return a new expansion object with the same nbasis and nfn
-
permute_basis
(permutation)¶ Reorder the coefficients for a given permutation of basis functions (rows).
Arguments:
- permutation
- An integer numpy array that defines the new order of the basis functions.
-
permute_orbitals
(permutation)¶ Reorder the coefficients for a given permutation of orbitals (columns).
Arguments:
- permutation
- An integer numpy array that defines the new order of the orbitals.
-
randomize
()¶ Fill with random normal data
-
rotate_2orbitals
(angle=0.7853981633974483, index0=None, index1=None)¶ Rotate two orbitals
Optional arguments:
- angle
- The rotation angle, defaults to 45 deg.
- index0, index1
- The orbitals to rotate, defaults to HOMO and LUMO,
The attributes
energies
andoccupations
are not altered.
-
rotate_random
()¶ Apply random unitary transformation distributed with Haar measure
The attributes
energies
andoccupations
are not altered.
-
swap_orbitals
(swaps)¶ Change the order of the orbitals using pair-exchange
Arguments:
- swaps
- An integer numpy array with two columns where every row corresponds to one swap.
The attributes
energies
andoccupations
are also reordered.
-
to_dm
(out=None, factor=1.0, clear=True, other=None)¶ Compute the density matrix
Optional arguments:
- out
- An output density matrix (DenseTwoIndex instance).
- factor
- The density matrix is multiplied by the given scalar.
- clear
- When set to False, the output density matrix is not zeroed first.
- other
- Another DenseExpansion object to construct a transfer-density matrix.
-
to_hdf5
(grp)¶ Dump this object in an h5py.Group
Arguments:
- grp
- An h5py.Group object.
-
coeffs
¶ The matrix with the expansion coefficients
-
energies
¶ The orbital energies
-
homo_energy
¶ Return a homo energy
Optional arguments:
- offset
- By default, the (highest) homo energy is returned. When this index is above zero, the corresponding lower homo energy is returned.
-
lumo_energy
¶ Return a lumo energy
Optional arguments:
- offset
- By default, the (lowest) lumo energy is returned. When this index is above zero, the corresponding higher homo energy is returned.
-
nbasis
¶ The number of basis functions
-
nfn
¶ The number of orbitals (or functions in general)
-
occupations
¶ The orbital occupations
-
class
horton.matrix.dense.
DenseTwoIndex
(nbasis, nbasis1=None)¶ Bases:
horton.matrix.base.TwoIndex
Arguments:
- nbasis
- The number of basis functions. (Number of rows. Also number of columns, unless nbasis1 is given.)
Optional arguments:
- nbasis1
- When given, this is the number of columns (second index).
Note that by default the two-index object is assumed to be Hermitian. Only when nbasis1 is given, this assumption is dropped.
-
__init__
(nbasis, nbasis1=None)¶ Arguments:
- nbasis
- The number of basis functions. (Number of rows. Also number of columns, unless nbasis1 is given.)
Optional arguments:
- nbasis1
- When given, this is the number of columns (second index).
Note that by default the two-index object is assumed to be Hermitian. Only when nbasis1 is given, this assumption is dropped.
-
assign
(other, ind=None)¶ Assign a new contents to the two-index object
Arguments:
- other
- The new data, may be DenseTwoIndex, DenseOneIndex, a scalar value, or an ndarray.
- ind
- If provided, only these elements (of DenseOneIndex) are assigned
-
assign_diagonal
(value)¶ Set diagonal elements to value
Arguments:
- value
- Either a scalar or a DenseOneIndex object
-
assign_dot
(other, tf2)¶ Dot product of orbitals in a DenseExpansion and TwoIndex object
Arguments:
- other
- An expansion object with input orbitals
- tf2
- A two-index object
The transformed array is stored in self.
-
assign_two_index_transform
(ao_integrals, exp0, exp1=None)¶ Perform two index transformation:
exp0.T * ao_integrals * exp0
Arguments:
- ao_integrals
- Something computed in the atomic orbital basis
- exp0
- The molecular orbitals.
Optional arguments:
- exp1
- When given, the following is computed:
exp0.T * ao_integrals * exp1
-
change_basis_signs
(signs)¶ Correct for different sign conventions of the basis functions.
Arguments:
- signs
- A numpy array with sign changes indicated by +1 and -1.
-
clear
()¶ Reset all elements to zero.
-
contract_to_one
(subscripts, out=None, factor=1.0, clear=True, begin0=0, end0=None, begin1=0, end1=None)¶ Contract self to OneIndex.
Arguments:
- subscripts
ab->b
: contract first index,ab->a
: contract second index.
Optional arguments
- out, factor, clear
- See
DenseLinalgFactory.einsum()
- begin0, end0, begin1, end1
- Can be used to contract only a part of the two-index object
Returns: the contracted one-index object.
-
contract_two
(subscripts, two, begin0=0, end0=None, begin1=0, end1=None)¶ Compute the trace using with other two-index objects
Arguments:
- subscripts
ab,ba
: trace of matrix product.ab,ab
: trace after element-wise multiplication (expectation_value).- two
- A DenseTwoIndex instance
- begin0, end0, begin1, end1
- Can be used to select a subblock of the (other) object. When not given, the full range is used.
-
contract_two_to_one
(subscripts, other, out=None, factor=1.0, clear=True, begin0=0, end0=None, begin1=0, end1=None)¶ Contract two TwoIndex objects to a one OneIndex.
Arguments:
- subscripts:
ab,ab->b
: contract first index.ab,ab->a
: contract second index.- other
- The second DenseTwoIndex object.
Optional arguments:
- out, factor, clear
- See
DenseLinalgFactory.einsum()
- begin0, end0, begin1, end1
- Can be used to contract only a part of the other two-index object. When not given, the full range is taken.
Returns: the contracted one-index object.
-
contract_two_to_two
(subscripts, other, out=None, factor=1.0, clear=True, begin0=0, end0=None, begin1=0, end1=None)¶ Contract two TwoIndex objects to a one twoIndex.
Arguments:
- subscripts:
ab,cb->ac
,ab,ca->cb
,ab,bc->ac
,ab,ac->cb
,ab,ac->bc
,ab,cb->ca
- other
- The second DenseTwoIndex object.
Optional arguments:
- out, factor, clear
- See
DenseLinalgFactory.einsum()
- begin0, end0, begin1, end1
- Can be used to select a subblock of the (other) object. When not given, the full range is used.
Returns: the contracted two-index object.
-
copy
(begin0=0, end0=None, begin1=0, end1=None)¶ Return a copy of (a part of) the object
Optional arguments:
- begin0, end0, begin1, end1
- Can be used to select a subblock of the object. When not given, the full range is used. Only when the ranges are equal for both axes, an Hermitian two-index may be returned.
-
copy_diagonal
(out=None, begin=0, end=None)¶ Copy (part of) the diagonal of the two-index object
Optional arguments:
- out
- The output argument (DenseOneIndex with proper size).
- begin, end
- Can be used to select a range of the diagonal. If not given, then the entire diagonal is copied.
-
copy_slice
(ind, out=None)¶ Copy slice of the two-index object into one-index object
Arguments:
- ind
- 2-Tuple of 1-dim arrays with row and column indices of TwoIndex object to be copied.
Optional arguments:
- out
- The output argument (DenseOneIndex with proper size).
-
diagonalize
()¶ Return eigenvalues of two-index object
-
distance_inf
(other)¶ The infinity norm distance between self and other
Arguments:
- other
- A DenseTwoIndex instance.
-
classmethod
from_hdf5
(grp)¶ Construct an instance from data previously stored in an h5py.Group.
Arguments:
- grp
- An h5py.Group object.
-
get_element
(i, j)¶ Return a matrix element
-
iabs
()¶ In-place absolute values
-
iadd
(other, factor=1.0, begin0=0, end0=None, begin1=0, end1=None, transpose=False)¶ Add another DenseTwoIndex object in-place, multiplied by factor. If begin0, end0, begin1, end1 are specified, other is added to the selected range.
Arguments:
- other
- A DenseTwoIndex, DenseOneIndex instance or float to be added
Optional arguments:
- factor
- A scalar factor
- begin0, end0, begin1, end1
- When given, specify the ranges where the contribution will be added. When not given, the full range is used.
-
iadd_contract_two_one
(subscripts, two, one, factor=1.0, begin0=0, end0=None, begin1=0, end1=None, begin2=0, end2=None)¶ In-place addition of a contraction of a two- and one-index object
Arguments:
- two, one
- The two- and one-index objects, respectively. Instances of
- subscripts
ab,b->ab
: contract with the first index of the two-index object.ab,a->ab
: contract with the second index of the two-index.ba,a->ab
: contract with the second index and transpose. object.- factor
- The term added is scaled by this factor.
- begin0, end0, begin1, end1
- Can be used to select a subblock of the two-index object (two). When not given, the full range is used.
- begin2, end2
- Can be used to select a subblock of the one-index object (one). When not given, the full range is used.
-
iadd_dot
(other0, other1, factor=1.0, begin0=0, end0=None, begin1=0, end1=None, begin2=0, end2=None, begin3=0, end3=None)¶ In-place addition of dot product:
other0 * other1
Arguments:
- other0, other1
- Two-index objects that go into the kronecker product.
Optional arguments:
- factor
- The term added is scaled by this factor.
- begin0, end0, begin1, end1
- Can be used to select a subblock of the other0 object. When not given, the full range is used.
- begin2, end2, begin3, end3
- Can be used to select a subblock of the other1 object. When not given, the full range is used.
-
iadd_dott
(other0, other1, factor=1.0)¶ In-place addition of dot product:
other0 * other1.T
Arguments:
- other0, other1
- Two-index objects that go into the kronecker product.
Optional arguments:
- factor
- The term added is scaled by this factor.
-
iadd_kron
(other0, other1, factor=1.0)¶ In-place addition of kronecker product of two other DenseTwoIndex
Arguments:
- other0, other1
- Two-index objects that go into the kronecker product.
Optional arguments:
- factor
- The term added is scaled by this factor.
-
iadd_mult
(other0, other1, factor=1.0, begin0=0, end0=None, begin1=0, end1=None, begin2=0, end2=None, begin3=0, end3=None, transpose0=False, transpose1=False)¶ In-place addition of multiplication:
other0 * other1
Arguments:
- other0, other1
- Two-index objects that go into the product.
Optional arguments:
- factor
- The term added is scaled by this factor.
- begin0, end0, begin1, end1
- Can be used to select a subblock of the other1 object. When not given, the full range is used.
- begin2, end2, begin3, end3
- Can be used to select a subblock of the other0 object. When not given, the full range is used.
- transpose0, transpose1
- Can be used to select transpose of other0, other1 objects
-
iadd_one_mult
(other0, other1, factor=1.0, transpose0=False, transpose1=False)¶ In-place addition of multiplication:
other0 * other1
Arguments:
- other0, other1
- One-index objects that go into the product.
Optional arguments:
- factor
- The term added is scaled by this factor.
- transpose0, transpose1
- Can be used to select transpose of one-index objects.
-
iadd_outer
(other0, other1, factor=1.0)¶ In-place addition of outer product of two other DenseTwoIndex
Arguments:
- other0, other1
- Two-index objects that go into the outer product. They are raveled before taking the outer product.
Optional arguments:
- factor
- The term added is scaled by this factor.
-
iadd_shift
(lshift)¶ Add positive shift to elements. If negative subtract shift
Arguments:
- lshift
- A scalar used to augment the matrix elements.
-
iadd_slice
(other, factor=1.0, begin0=0, end0=None, begin1=0, end1=None)¶ Add slice of another DenseTwoIndex object in-place, multiplied by factor. The slice of other is added to the full range of the array.
Arguments:
- other
- A DenseTwoIndex instance to be added
Optional arguments:
- factor
- A scalar factor
- begin0, end0, begin1, end1
- When given, specify the ranges of contribution to be added. When not given, the full range is used.
-
iadd_t
(other, factor=1.0, begin0=0, end0=None, begin1=0, end1=None)¶ See
DenseTwoIndex.iadd()
, transpose=True
-
iadd_tdot
(other0, other1, factor=1.0)¶ In-place addition of dot product:
other0.T * other1
Arguments:
- other0, other1
- Two-index objects that go into the kronecker product.
Optional arguments:
- factor
- The term added is scaled by this factor.
-
idot
(other)¶ In-place dot product: self = self * other
Arguments:
- other
- The other array
-
imul
(other, factor=1.0, begin0=0, end0=None, begin1=0, end1=None)¶ In-place element-wise multiplication:
self *= other * factor
Arguments:
- other
- A DenseTwoIndex instance.
Optional arguments:
- factor
- The two-index object is scaled by this factor.
- begin0, end0, begin1, end1
- Can be used to select a subblock of the (other) object. When not given, the full range is used.
-
imul_t
(other, factor=1.0, begin0=0, end0=None, begin1=0, end1=None)¶ In-place element-wise multiplication:
self *= other.T * factor
Arguments:
- other
- A DenseTwoIndex instance.
Optional arguments:
- factor
- The two-index object is scaled by this factor.
- begin0, end0, begin1, end1
- Can be used to select a subblock of the (other) object. When not given, the full range is used.
-
inner
(vec0, vec1)¶ Compute an inner product of two vectors using the two-index as a metric
Arguments:
- vec0, vec1
- The vectors, either DenseOneIndex or numpy arrays.
-
inverse
()¶ Return the inverse of two-index object
-
iortho
()¶ In-place orthogonalization
Arguments:
-
is_shape_symmetric
(symmetry)¶ Check whether the symmetry argument matches the shape
-
is_symmetric
(symmetry=2, rtol=1e-05, atol=1e-08)¶ Check the symmetry of the array.
Optional arguments:
- symmetry
- The symmetry to check. See Handling of index symmetry for more details.
- rtol and atol
- relative and absolute tolerance. See to
np.allclose
.
-
iscale
(factor)¶ In-place multiplication with a scalar
Arguments:
- factor
- A scalar factor.
-
itranspose
()¶ In-place transpose
-
new
()¶ Return a new two-index object with the same nbasis (and nbasis1)
-
permute_basis
(permutation)¶ Reorder the coefficients for a given permutation of basis functions.
The same permutation is applied to all indexes.
Arguments:
- permutation
- An integer numpy array that defines the new order of the basis functions.
-
randomize
()¶ Fill with random normal data
-
set_element
(i, j, value, symmetry=2)¶ Set a matrix element
Arguments:
- i, j
- The matrix indexes to be set
- value
- The value to be assigned to the matrix element.
Optional arguments:
- symmetry
- When 2 (the default), the element (j,i) is set to the same value. When set to 1 the opposite off-diagonal is not set. See Handling of index symmetry for more details.
-
sqrt
()¶ Return the real part of the square root of two-index object
-
sum
(begin0=0, end0=None, begin1=0, end1=None)¶ Return the sum of all elements (in the selected range)
Optional arguments:
- begin0, end0, begin1, end1
- Can be used to select a subblock of the object to be contracted.
-
symmetrize
(symmetry=2)¶ Symmetrize in-place
Optional arguments:
- symmetry
- The symmetry to impose. See Handling of index symmetry for more details.
-
to_hdf5
(grp)¶ Dump this object in an h5py.Group
Arguments:
- grp
- An h5py.Group object.
-
trace
(begin0=0, end0=None, begin1=0, end1=None)¶ Return the trace of the two-index object.
Optional arguments:
- begin0, end0, begin1, end1
- Can be used to select a subblock of the object to be contracted.
-
nbasis
¶ The number of basis functions
-
nbasis1
¶ The other size of the two-index object
-
ndim
¶ The number of axes in the N-index object.
-
shape
¶ The shape of the object
-
class
horton.matrix.dense.
DenseThreeIndex
(nbasis, nbasis1=None, nbasis2=None)¶ Bases:
horton.matrix.base.ThreeIndex
Arguments:
- nbasis
- The number of basis functions.
-
__init__
(nbasis, nbasis1=None, nbasis2=None)¶ Arguments:
- nbasis
- The number of basis functions.
-
assign
(other)¶ Assign with the contents of another object
Arguments:
- other
- Another DenseThreeIndex object or a scalar.
-
change_basis_signs
(signs)¶ Correct for different sign conventions of the basis functions.
Arguments:
- signs
- A numpy array with sign changes indicated by +1 and -1.
-
clear
()¶ Reset all elements to zero.
-
contract_to_two
(subscripts, out=None, factor=1.0, clear=True, begin0=0, end0=None, begin1=0, end1=None, begin2=0, end2=None)¶ Contract self to TwoIndex.
Arguments:
- subscripts
abc->ac
: contract first index.
Optional arguments
- out, factor, clear
- See
DenseLinalgFactory.einsum()
- begin0, end0, begin1, end1, begin2, end2
- Can be used to contract only a part of the three-index object. When not given, the full range is used.
Returns: the contracted two-index object.
-
contract_two_to_three
(subscripts, two, out=None, factor=1.0, clear=True, begin0=0, end0=None, begin1=0, end1=None, begin2=0, end2=None)¶ Contracts with two-index to obtain three-index.
Arguments:
- subscripts
- One of
abc,ad->cdb
,abc,ad->cbd
,abc,da->dbc
,abc,da->bdc
- inp
- A DenseTwoIndex input object.
Optional arguments:
- out, factor, clear
- See
DenseLinalgFactory.einsum()
- begin0, end0, begin1, end1, begin2, end2
- Can be used to contract only a part of the three-index object
-
contract_two_to_two
(subscripts, inp, out=None, factor=1.0, clear=True)¶ Contracts self with two-index to obtain two-index.
Arguments:
- subscripts
- One of
abc,ab->ac
,abc,ab->ca
,abc,bc->ba
,abc,bc->ab
,abc,cb->ac
. - inp
- A DenseTwoIndex input object.
Optional arguments:
- out, factor, clear
- See
DenseLinalgFactory.einsum()
-
copy
(begin0=0, end0=None, begin1=0, end1=None, begin2=0, end2=None)¶ Return a copy of (a part of) the object
Optional arguments:
- begin0, end0, begin1, end1, begin2, end2
- Can be used to select a subblock of the object. When not given, the full range is used.
-
classmethod
from_hdf5
(grp)¶ Construct an instance from data previously stored in an h5py.Group.
Arguments:
- grp
- An h5py.Group object.
-
get_element
(i, j, k)¶ Return a matrix element
-
iadd
(other, factor=1.0)¶ Add another DenseThreeIndex object in-place, multiplied by factor
Arguments:
- other
- A DenseThreeIndex instance to be added
Optional arguments:
- factor
- The added term is scaled by this factor.
-
iadd_expand_two_one
(subscripts, two, one, factor=1.0)¶ In-place addition of expanded two-index and one-index to three-index.
Arguments:
- subscripts
- Contraction type:
ab,c->cab
orab,c->acb
. - two
- A DenseTwoIndex object.
- one
- A DenseOneIndex object.
Optional arguments:
- factor
- The added term is scaled by this factor.
-
iadd_expand_two_two
(subscripts, two0, two1, factor=1.0, begin0=0, end0=None, begin1=0, end1=None, begin2=0, end2=None, begin3=0, end3=None)¶ In-place addition of expanded two-index and two-index to three-index.
Arguments:
- subscripts
- Contraction type:
ac,bc->abc
,ab,bc->abc
,ab,ac->acb
,cb,ac->acb
,ac,ab->abc
,ab,ac->abc
- two0, two1
- A DenseTwoIndex object.
Optional arguments:
- factor
- The added term is scaled by this factor.
- begin0, end0, begin1, end1
- Can be used to contract only a part of the two1 object. When not given, the full range is used.
- begin2, end2, begin3, end3
- Can be used to contract only a part of the two0 object. When not given, the full range is used.
-
iadd_slice
(other, factor=1.0, begin0=0, end0=None, begin1=0, end1=None, begin2=0, end2=None)¶ Add slice of another DenseThreeIndex object in-place, multiplied by factor
Arguments:
- other
- A DenseThreeIndex instance to be added
Optional arguments:
- factor
- The added term is scaled by this factor.
- begin0, end0, begin1, end1, begin2, end2
- Can be used to add only a part of the three-index object
-
iscale
(factor)¶ In-place multiplication with a scalar
Arguments:
- factor
- A scalar factor.
-
new
()¶ Return a new three-index object with the same nbasis
-
permute_basis
(permutation)¶ Reorder the coefficients for a given permutation of basis functions.
Arguments:
- permutation
- An integer numpy array that defines the new order of the basis functions.
-
randomize
()¶ Fill with random normal data
-
set_element
(i, j, k, value)¶ Set a matrix element
-
to_hdf5
(grp)¶ Dump this object in an h5py.Group
Arguments:
- grp
- An h5py.Group object.
-
nbasis
¶ The number of basis functions
-
nbasis1
¶ The number of basis functions
-
nbasis2
¶ The number of basis functions
-
ndim
¶ The number of axes in the N-index object.
-
shape
¶ The shape of the object
-
class
horton.matrix.dense.
DenseFourIndex
(nbasis, nbasis1=None, nbasis2=None, nbasis3=None)¶ Bases:
horton.matrix.base.FourIndex
Arguments:
- nbasis
- The number of basis functions.
Optional arguments:
nbasis1, nbasis2, nbasis3
-
__init__
(nbasis, nbasis1=None, nbasis2=None, nbasis3=None)¶ Arguments:
- nbasis
- The number of basis functions.
Optional arguments:
nbasis1, nbasis2, nbasis3
-
assign
(other)¶ Assign with the contents of another object
Arguments:
- other
- Another DenseFourIndex object or a np ndarrray.
-
assign_four_index_transform
(ao_integrals, exp0, exp1=None, exp2=None, exp3=None, method=’tensordot’)¶ Perform four index transformation.
Arguments:
- oa_integrals
- A DenseFourIndex with integrals in atomic orbitals.
- exp0
- A DenseExpansion object with molecular orbitals
Optional arguments:
- exp1, exp2, exp3
- Can be provided to transform each index differently.
- method
- Either
einsum
ortensordot
(default).
-
change_basis_signs
(signs)¶ Correct for different sign conventions of the basis functions.
Arguments:
- signs
- A numpy array with sign changes indicated by +1 and -1.
-
clear
()¶ Reset all elements to zero.
-
contract_four
(subscripts, other, factor=1.0, begin0=0, end0=None, begin1=0, end1=None, begin2=0, end2=None, begin3=0, end3=None)¶ Compute the trace with other four-index object
Arguments:
- subscripts
abcd,abcd
,abcd,adcb
,abcd,acbd
,abcd,acdb
- other
- A DenseFourIndex instance
Optional arguments:
- factor
- A scalar factor
- begin0, end0, begin1, end1, begin2, end2, begin3, end3
- Can be used to select a subblock of the other object. When not given, the full range is used.
-
contract_four_to_four
(subscripts, four, out=None, factor=1.0, clear=True, begin0=0, end0=None, begin1=0, end1=None, begin2=0, end2=None, begin3=0, end3=None)¶ Contracts with a four-index object to obtain a four-index object.
Arguments:
- subscripts
- Any of
abcd,cedf->abfe
,abcd,cefd->abfe
,abcd,cedf->feab
,abcd,cefd->feab
,abcd,eadf->fbce
,abcd,eadf->cefb
,abcd,aedf->cbfe
,abcd,aedf->fecb
,abcd,acef->ebfd
,abcd,bdef->aecf
,abcd,cedf->abef
,abcd,cefd->abef
,abcd,cedf->efab
,abcd,cefd->efab
,abcd,eadf->ebcf
,abcd,eadf->cfeb
,abcd,eadf->cbef
,abcd,eadf->efcb
,abcd,eafd->cbef
,abcd,eafd->efcb
- four
- An instance of DenseFourIndex.
Optional arguments:
- out, factor, clear
- See
DenseLinalgFactory.einsum()
- begin0, end0, begin1, end1, begin2, end2, begin3, end3
- Can be used to select a subblock of the other four-index object (four). When not given, the full range is used.
-
contract_four_to_two
(subscripts, four, out=None, factor=1.0, clear=True, begin0=0, end0=None, begin1=0, end1=None, begin2=0, end2=None, begin3=0, end3=None)¶ Contracts with a four-index object to obtain a two-index object.
Arguments:
- subscripts
- Any of
abcd,aced->be
,abcd,acde->be
,abcd,aced->eb
,abcd,acde->eb
,abcd,ecbd->ae
,abcd,ecdb->ae
,abcd,ecdb->ea
,abcd,ecbd->ea
,abcd,acbe->ed
,abcd,aceb->ed
,abcd,aebd->ce
,abcd,aedb->ce
- four
- An instance of DenseFourIndex.
Optional arguments:
- out, factor, clear
- See
DenseLinalgFactory.einsum()
- begin0, end0, begin1, end1, begin2, end2, begin3, end3
- Can be used to select a subblock of the other four-index object (four). When not given, the full range is used.
-
contract_three_to_three
(subscripts, three, out=None, factor=1.0, clear=True, begin0=0, end0=None, begin1=0, end1=None, begin2=0, end2=None)¶ Contracts with a three-index object to obtain a three-index object.
Arguments:
- subscripts
- Any of
abcd,ace->ebd
,abcd,ebd->ace
- three
- An instance of DenseThreeIndex.
Optional arguments:
- out, factor, clear
- See
DenseLinalgFactory.einsum()
- begin0, end0, begin1, end1, begin2, end2
- Can be used to select a subblock of the three-index object (three). When not given, the full range is used.
-
contract_to_two
(subscripts, out=None, factor=1.0, clear=True, begin0=0, end0=None, begin1=0, end1=None, begin2=0, end2=None, begin3=0, end3=None)¶ Contract four-index object to two-index object
Arguments:
- subscripts
abcb->ac
,abbc->ac
Optional arguments:
- out, factor, clear
- See
DenseLinalgFactory.einsum()
- begin0, end0, begin1, end1, begin2, end2, begin3, end3
- Can be used to select a subblock of the object. When not given, the full range is used.
-
contract_two
(subscripts, other, factor=1.0, begin0=0, end0=None, begin1=0, end1=None, begin2=0, end2=None, begin3=0, end3=None)¶ Compute diagonal trace with two-index objects
Arguments:
- subscripts
aabb,ab
.- other
- A DenseTwoIndex instance
Optional arguments:
- factor
- A scalar factor. See
DenseLinalgFactory.einsum()
- begin0, end0, begin1, end1, begin2, end2, begin3, end3
- Can be used to select a subblock of the four-index object. When not given, the full range is used.
-
contract_two_to_four
(subscripts, two, out=None, factor=1.0, clear=True, begin0=0, end0=None, begin1=0, end1=None, begin2=0, end2=None, begin3=0, end3=None, begin4=0, end4=None, begin5=0, end5=None)¶ Contracts with a two-index object to obtain a four-index object.
Arguments:
- subscripts
- Any of
abcd,cd->acbd
,abcd,cd->acdb
,abcd,cb->acdb
,abcd,cb->acbd
,abcd,ab->acbd
,abcd,ab->acdb
,abcd,ad->acbd
,abcd,ad->acdb
,abcd,ad->abcd
,abcd,ad->abdc
,abcd,bd->abcd
,abcd,bd->abdc
,abcd,bc->abdc
,abcd,bc->abcd
,abcd,ac->abcd
,abcd,ac->abdc
,abcd,bd->cabd
,abcd,bc->dabc
,abcd,ca->cabd
,abcd,da->dabc
,abcd,dc->dabc
,abcd,ba->dabc
,aabc,dc->adbc
,aabc,db->adbc
,abcc,ad->abcd
,abcc,bd->abcd
,abcd,bc->acbd
,abcd,eb->aecd
,abcd,eb->cdae
,abcd,ed->ceab
,abcd,ed->abce
,abcd,ae->cdeb
,abcd,ae->ebcd
,abcd,ce->edab
,abcd,ce->abed
,abcd,ab->abcd
,abcd,cb->abcd
,abcd,ec->eadb
,abcd,ec->dbea
,abcd,ae->cedb
,abcd,ae->dbce
,abcd,ec->ebad
,abcd,ed->ebac
,abcd,ec->adeb
,abcd,ed->aceb
,abcd,ae->debc
,abcd,ae->cebd
,abcd,ae->bcde
,abcd,ae->bdce
- two
- An instance of DenseTwoIndex.
Optional arguments:
- out, factor, clear
- See
DenseLinalgFactory.einsum()
- begin0, end0, begin1, end1, begin2, end2, begin3, end3
- Can be used to select a subblock of the four-index object. When not given, the full range is used.
- begin4, end4, begin5, end5
- Can be used to select a subblock of the two-index object (two). When not given, the full range is used.
-
contract_two_to_three
(subscripts, two, out=None, factor=1.0, clear=True, begin0=0, end0=None, begin1=0, end1=None, begin2=0, end2=None, begin3=0, end3=None)¶ Contracts with a two-index object to obtain a three-index object.
Arguments:
- subscripts
- Any of
aabc,ad->bcd
,'abcd,ac->bdc
,abcd,ad->bcd
,abcd,bc->abd
,abcd,ac->abd
,abcc,dc->dab
,abcd,ac->bcd
,abcc,cd->dab
,abcc,dc->abd
,abcb,db->adc
,abcc,dc->adb
,abcb,db->dac
,abcc,cd->abd
- two
- An instance of DenseTwoIndex.
Optional arguments:
- out, factor, clear
- See
DenseLinalgFactory.einsum()
- begin0, end0, begin1, end1, begin2, end2, begin3, end3
- Can be used to select a subblock of the four-index object. When not given, the full range is used.
-
contract_two_to_two
(subscripts, two, out=None, factor=1.0, clear=True, begin0=0, end0=None, begin1=0, end1=None, begin2=0, end2=None, begin3=0, end3=None, begin4=0, end4=None, begin5=0, end5=None)¶ Contract self with a two-index to obtain a two-index.
Arguments:
- subscripts
- Any of
abcd,bd->ac
(direct),abcd,cb->ad
(exchange),aabb,cb->ac
,abcc,bc->ab
,aabc,ab->bc
,aabc,ac->bc
,abcc,ac->ab
,abcb,cb->ac
,abcb,ab->ac
,abcc,bc->ba
,abcc,bc->ab
,abcd,ac->db
,abcd,ad->cb
,abcd,ac->bd
,abcd,ad->bc
,abcd,ab->cd
- two
- The input two-index object. (DenseTwoIndex)
Optional arguments:
- out, factor, clear
- See
DenseLinalgFactory.einsum()
- begin0, end0, begin1, end1, begin2, end2, begin3, end3
- Can be used to select a subblock of the four-index object. When not given, the full range is used.
- begin4, end4, begin5, end5
- Can be used to select a subblock of the two-index object (two). When not given, the full range is used.
-
copy
(begin0=0, end0=None, begin1=0, end1=None, begin2=0, end2=None, begin3=0, end3=None)¶ Return a copy of (a part of) the object
Optional arguments:
- begin0, end0, begin1, end1, begin2, end2, begin3, end3
- Can be used to select a subblock of the object. When not given, the full range is used.
-
classmethod
from_hdf5
(grp)¶ Construct an instance from data previously stored in an h5py.Group.
Arguments:
- grp
- An h5py.Group object.
-
get_element
(i, j, k, l)¶ Return a matrix element
-
iadd
(other, factor=1.0)¶ Add another DenseFourIndex object in-place, multiplied by factor
Arguments:
- other
- A DenseFourIndex instance to be added
Optional arguments:
- factor
- The added term is scaled by this factor.
-
iadd_exchange
()¶ In-place addition of its own exchange contribution
-
iadd_expand_three_to_four
(axis, three, factor=1.0)¶ Expand three-index object along one axis and add it to four-index object.
Arguments:
- axis
- Any of
1-3-1-2
(abac += abc),1-2-1-2
(abca += abc),0-2-0-2
(abcb += abc),0-3-0-2
(abbc += abc),0-2-0-1
(abcb += acb). The expansion is performed for repeated indices in the four-index object - three
- An instance of DenseThreeIndex.
Optional arguments:
- factor
- See
DenseLinalgFactory.einsum()
-
iadd_expand_two_to_four
(axis, two, factor=1.0, begin0=0, end0=None, begin1=0, end1=None)¶ Expand two-index object along two axes and add it to four-index object.
Arguments:
- axis
- Any of
1-3
(abac += bc),0-2
(abcb += ac),0-3
(abbc += ac),1-2
(abca += bc),diag
(abab += ab). The expansion is performed for the repeated indices in the four-index object - two
- An instance of DenseTwoIndex.
Optional arguments:
- factor
- See
DenseLinalgFactory.einsum()
- begin0, end0, begin1, end1
- Can be used to select a subblock of the two-index object (two). When not given, the full range is used.
-
imul
(other, factor=1.0)¶ In-place element-wise multiplication:
self *= other * factor
Arguments:
- other
- A DenseFourIndex instance.
Optional arguments:
- factor
- The four-index object is scaled by this factor.
-
is_shape_symmetric
(symmetry)¶ Check whether the symmetry argument matches the shape
-
is_symmetric
(symmetry=8, rtol=1e-05, atol=1e-08)¶ Check the symmetry of the array.
Optional arguments:
- symmetry
- The symmetry to check. See Handling of index symmetry for more details. In addition to 1, 2, 4, 8, also ‘cdab’ is supported.
- rtol and atol
- relative and absolute tolerance. See to
np.allclose
.
-
iscale
(factor)¶ In-place multiplication with a scalar
Arguments:
- factor
- A scalar factor.
-
itranspose
()¶ In-place transpose:
0,1,2,3 -> 1,0,3,2
-
new
()¶ Return a new four-index object with the same nbasis
-
permute_basis
(permutation)¶ Reorder the coefficients for a given permutation of basis functions.
Arguments:
- permutation
- An integer numpy array that defines the new order of the basis functions.
-
randomize
()¶ Fill with random normal data
-
reshape
(shape)¶ Reshape array
Optional arguments:
- shape
- List containing the new dimension of each axis.
-
set_element
(i, j, k, l, value, symmetry=8)¶ Set a matrix element
Arguments:
- i, j, k, l
- The matrix indexes to be set
- value
- The value to be assigned to the matrix element.
Optional arguments:
- symmetry
- The level of symmetry to be enforced when setting the matrix element. See Handling of index symmetry for more details.
-
slice_to_four
(subscripts, out=None, factor=1.0, clear=True, begin0=0, end0=None, begin1=0, end1=None, begin2=0, end2=None, begin3=0, end3=None)¶ Returns a four-index contraction of the four-index object.
Arguments:
- superscripts
- Any of
abcd->abcd
Optional arguments:
- out, factor, clear
- See
DenseLinalgFactory.einsum()
- begin0, end0, begin1, end1, begin2, end2, begin3, end3
- Can be used to select a subblock of the object. When not given, the full range is used.
-
slice_to_three
(subscripts, out=None, factor=1.0, clear=True)¶ Returns a three-index contraction of the four-index object.
Arguments:
- superscripts
- Any of
abcc->bac
,abcc->abc
,abcb->abc
,abbc->abc
Optional arguments:
- out, factor, clear
- See
DenseLinalgFactory.einsum()
-
slice_to_two
(subscripts, out=None, factor=1.0, clear=True, begin0=0, end0=None, begin1=0, end1=None, begin2=0, end2=None, begin3=0, end3=None)¶ Returns a two-index contraction of the four-index object.
Arguments:
- superscripts
- Any of
aabb->ab
,abab->ab
,abba->ab
Optional arguments:
- out, factor, clear
- See
DenseLinalgFactory.einsum()
- begin0, end0, begin1, end1, begin2, end2, begin3, end3
- Can be used to select a subblock of the object. When not given, the full range is used.
-
sum
()¶ Return the sum of all elements
-
symmetrize
(symmetry=8)¶ Symmetrize in-place
Optional arguments:
- symmetry
- The symmetry to impose. See Handling of index symmetry for more details.
-
to_hdf5
(grp)¶ Dump this object in an h5py.Group
Arguments:
- grp
- An h5py.Group object.
-
nbasis
¶ The number of basis functions
-
nbasis1
¶ The number of basis functions
-
nbasis2
¶ The number of basis functions
-
nbasis3
¶ The number of basis functions
-
ndim
¶ The number of axes in the N-index object.
-
shape
¶ The shape of the object