MADNESS  0.10.1
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
madness::BaseTensor Class Reference

The base class for tensors defines generic capabilities. More...

#include <basetensor.h>

Inheritance diagram for madness::BaseTensor:
Inheritance graph
[legend]

Public Member Functions

 BaseTensor ()
 
virtual ~BaseTensor ()
 
bool conforms (const BaseTensor *t) const
 Returns true if this and *t are the same shape and size. More...
 
long dim (int i) const
 Returns the size of dimension i. More...
 
const long * dims () const
 Returns the array of tensor dimensions. More...
 
long id () const
 Returns the typeid of the tensor (c.f., TensorTypeData<T> ) More...
 
bool iscontiguous () const
 Returns true if the tensor refers to contiguous memory locations. More...
 
long ndim () const
 Returns the number of dimensions in the tensor. More...
 
long size () const
 Returns the number of elements in the tensor. More...
 
long stride (int i) const
 Returns the stride associated with dimension i. More...
 
const long * strides () const
 Returns the array of tensor strides. More...
 

Static Public Member Functions

static int get_instance_count ()
 Returns the count of all current instances of tensors & slice tensors of all types. More...
 

Protected Member Functions

void cycledim_inplace (long shift, long start, long end)
 Cyclic shift of dimensions. More...
 
void flat_inplace ()
 Reshapes the tensor inplace into 1D. More...
 
void fusedim_inplace (long i)
 Fuses dimensions i and i+1. More...
 
void mapdim_inplace (const std::vector< long > &map)
 General permutation of dimensions. More...
 
void reshape_inplace (const std::vector< long > &d)
 Reshapes the tensor inplace. More...
 
void reshape_inplace (int ndimnew, const long *d)
 Reshapes the tensor inplace. More...
 
void set_dims_and_size (long nd, const long d[])
 
void splitdim_inplace (long i, long dimi0, long dimi1)
 Splits dimension i. More...
 
void swapdim_inplace (long i, long j)
 Swaps the dimensions. More...
 

Protected Attributes

long _dim [TENSOR_MAXDIM]
 Size of each dimension. More...
 
long _id
 Id from TensorTypeData<T> in type_data.h. More...
 
long _ndim
 Number of dimensions (-1=invalid; 0=no supported; >0=tensor) More...
 
long _size
 Number of elements in the tensor. More...
 
long _stride [TENSOR_MAXDIM]
 Increment between elements in each dimension. More...
 

Detailed Description

The base class for tensors defines generic capabilities.

The base class manages the size, dimension and stride information, and provides operations to manipulate them.

It also provides methods for type-safe operation on tensors using just the base class pointers. This interface is primarily useful only to the interface to Python, since Python is largely neutral to (and ignorant of) the type. These are still being re-implemented after the big clean up.

Since the base tensor class is virtual, you cannot have an instance of it. Thus, in addition to methods that return information or perform checks, there are two types of base tensor operations.

Constructor & Destructor Documentation

◆ BaseTensor()

madness::BaseTensor::BaseTensor ( )
inline

◆ ~BaseTensor()

virtual madness::BaseTensor::~BaseTensor ( )
inlinevirtual

Member Function Documentation

◆ conforms()

bool madness::BaseTensor::conforms ( const BaseTensor t) const
inline

Returns true if this and *t are the same shape and size.

References _dim, and _ndim.

Referenced by madness::Tensor< T >::conforms().

◆ cycledim_inplace()

void madness::BaseTensor::cycledim_inplace ( long  shift,
long  start,
long  end 
)
protected

Cyclic shift of dimensions.

Cyclic shift by nshift places of the inclusive range of dimensions [start,....,end].

References _dim, _ndim, _stride, TENSOR_ASSERT, and TENSOR_MAXDIM.

Referenced by madness::Tensor< T >::cycledim().

◆ dim()

long madness::BaseTensor::dim ( int  i) const
inline

Returns the size of dimension i.

References _dim.

Referenced by madness::ConvolutionData1D< Q >::ConvolutionData1D(), MolecularNuclearPotentialFunctor2::MolecularNuclearPotentialFunctor2(), madness::TensorIterator< T, Q, R >::TensorIterator(), madness::TensorTrain< T >::TensorTrain(), madness::BSHApply< T, NDIM >::add_coupling_and_levelshift(), madness::SeparatedConvolution< Q, NDIM >::apply(), madness::apply(), KPeriodicBSHOperator::apply(), madness::GTHPseudopotential< Q >::apply_potential(), madness::GTHPseudopotential< Q >::apply_potential_simple(), madness::projector_irrep::apply_symmetry_operators(), madness::SRConf< T >::check_dimensions(), madness::LowRankFunction< T, NDIM, LDIM >::check_orthonormality(), madness::PNO::check_orthonormality(), madness::SRConf< T >::check_right_orthonormality(), madness::Zcis::compute_fock_pt(), madness::CC2::compute_local_coupling(), madness::OEP::compute_ocep_correction(), madness::RandomizedMatrixDecomposition< T >::compute_range(), Calculation::compute_residuals(), madness::GenTensor< T >::convert_inplace(), madness::SRConf< T >::copy_slice(), madness::SCF::diag_fock_matrix(), madness::GenTensor< T >::dim(), madness::distributed_localize_PM(), madness::RandomizedMatrixDecomposition< T >::do_compute_range(), madness::FunctionImpl< T, NDIM >::do_print_grid(), madness::FunctionImpl< T, NDIM >::do_print_plane(), madness::TensorTrain< T >::emul(), madness::Function< T, NDIM >::eval_cube(), madness::fcube(), fixphases(), madness::SRConf< T >::flat_vector_with_weights(), fock_trace(), madness::TensorTrain< T >::fusedim(), madness::TensorTrain< T >::gaxpy(), madness::GradBSHOperator(), GradBSHOperator_Joel(), madness::GradCoulombOperator(), Test2::gradient(), madness::GradSlaterOperator(), madness::SRConf< T >::has_structure(), madness::QuasiNewton::hessian_update_bfgs(), madness::Solver< T, NDIM >::initial_guess(), madness::inner(), madness::FunctionImpl< T, NDIM >::inner_local(), madness::inner_result(), madness::Znemo::iterate(), madness::OEP::iterate(), madness::Localizer::localize(), madness::Localizer::localize_new(), madness::ConvolutionData1D< Q >::make_approx(), MiniDFT::make_bsh_operators(), make_bsh_operators(), madness::SCF::make_bsh_operators(), madness::TDHF::make_cis_matrix(), madness::Zcis::make_guess(), madness::TDHF::make_guess_from_initial_diagonalization(), madness::SRConf< T >::make_vector_dimensions(), madness::SRConf< T >::make_vector_with_weights(), madness::Solver< T, NDIM >::matrix_exponential(), matrix_exponential(), op(), madness::FunctionImpl< T, NDIM >::pointwise_multiplier< LDIM >::operator()(), madness::BSHApply< T, NDIM >::operator()(), madness::operator<<(), madness::TensorTrain< T >::operator=(), QuasiNewton::optimize(), madness::QuasiNewton::optimize(), madness::MolecularOptimizer::optimize_quasi_newton(), madness::ortho3(), madness::ortho5(), madness::outer(), madness::FunctionImpl< T, NDIM >::partial_inner_contract(), madness::PeriodicBSHOp(), madness::PeriodicBSHOpPtr(), madness::PeriodicCoulombOp(), madness::PeriodicCoulombOpPtr(), print(), madness::Solver< T, NDIM >::print_potential_matrix_eigs(), madness::Solver< T, NDIM >::print_tensor2d(), madness::TensorTrain< T >::ranks(), Plotter::read(), madness::FunctionImpl< T, NDIM >::read_grid(), madness::FunctionImpl< T, NDIM >::read_grid2(), madness::FunctionDefaults< NDIM >::recompute_cell_info(), madness::TensorTrain< T >::reconstruct(), madness::SCF::restart_aos(), madness::Solver< T, NDIM >::solve(), madness::Nemo::solve(), madness::TensorTrain< T >::splitdim(), Test1(), OptimizationTargetInterface::test_gradient(), madness::OptimizationTargetInterface::test_gradient(), test_orthogonalization(), madness::FunctionImpl< T, NDIM >::tnorm(), madness::TensorTrain< T >::trace(), madness::TensorTrain< T >::truncate(), madness::TensorTrain< T >::two_mode_representation(), madness::Localizer::undo_degenerate_rotations(), madness::Localizer::undo_reordering(), Test2::value(), madness::TensorTrain< T >::verify(), and madness::TensorTrain< T >::zero_me().

◆ dims()

const long* madness::BaseTensor::dims ( ) const
inline

Returns the array of tensor dimensions.

References _dim.

Referenced by madness::TensorTrain< T >::TensorTrain(), madness::abs(), madness::apply(), madness::arg(), madness::SVDTensor< T >::compute_randomized_svd(), madness::SVDTensor< T >::compute_svd(), madness::SVDTensor< T >::compute_svd_from_range(), madness::conj(), madness::convert(), madness::copy(), madness::TensorTrain< T >::decompose(), madness::GenTensor< T >::dims(), madness::SRConf< T >::emul(), madness::TensorTrain< T >::emul(), madness::SRConf< T >::general_transform(), madness::general_transform(), madness::SRConf< T >::get_configs(), madness::imag(), madness::SRConf< T >::inplace_add(), madness::SRConf< T >::make_vector_with_weights(), madness::TensorTrain< T >::operator TensorTrain< Q >(), madness::guessfactory::ExopUnaryOpStructure::operator()(), function_real2complex_op< Q, NDIM >::operator()(), madness::abs_square_op< Q, NDIM >::operator()(), madness::real_op< Q, NDIM >::operator()(), madness::imag_op< Q, NDIM >::operator()(), madness::abs_op< Q, NDIM >::operator()(), madness::conj_op< Q, NDIM >::operator()(), madness::function_real2complex_op< Q, NDIM >::operator()(), madness::BinaryOpStructure< NDIM >::operator()(), madness::UnaryOpStructure< NDIM >::operator()(), madness::SRConf< T >::operator=(), madness::TensorTrain< T >::operator=(), madness::real(), madness::SRConf< T >::reconstruct(), madness::TensorTrain< T >::reconstruct(), madness::SRConf< T >::set_size_and_dim(), madness::TensorTrain< T >::set_size_and_dim(), madness::TensorTrain< T >::splitdim(), madness::archive::ArchiveStoreImpl< Archive, Tensor< T > >::store(), madness::tensor_real2complex(), madness::SRConf< T >::transform(), madness::transform(), madness::SRConf< T >::transform_dir(), madness::transform_dir(), madness::SystolicMatrixAlgorithm< T >::unshuffle(), and madness::TensorTrain< T >::zero_me().

◆ flat_inplace()

void madness::BaseTensor::flat_inplace ( )
protected

Reshapes the tensor inplace into 1D.

Reshape the current tensor to be the same size and 1-d. It must be contiguous.

References _size, d(), iscontiguous(), set_dims_and_size(), and TENSOR_ASSERT.

◆ fusedim_inplace()

void madness::BaseTensor::fusedim_inplace ( long  i)
protected

Fuses dimensions i and i+1.

Fuse the contiguous dimensions i and i+1.

References _dim, _ndim, _stride, and TENSOR_ASSERT.

Referenced by madness::TensorTrain< T >::fusedim(), and madness::Tensor< T >::fusedim().

◆ get_instance_count()

static int madness::BaseTensor::get_instance_count ( )
inlinestatic

Returns the count of all current instances of tensors & slice tensors of all types.

Referenced by main().

◆ id()

long madness::BaseTensor::id ( ) const
inline

Returns the typeid of the tensor (c.f., TensorTypeData<T> )

References _id.

Referenced by madness::archive::ArchiveLoadImpl< Archive, Tensor< T > >::load(), madness::archive::ArchiveStoreImpl< Archive, Tensor< T > >::store(), and Test1().

◆ iscontiguous()

bool madness::BaseTensor::iscontiguous ( ) const
inline

◆ mapdim_inplace()

void madness::BaseTensor::mapdim_inplace ( const std::vector< long > &  map)
protected

General permutation of dimensions.

General permuation of the dimensions.

References _dim, _ndim, _stride, TENSOR_ASSERT, and TENSOR_MAXDIM.

Referenced by madness::Tensor< T >::mapdim().

◆ ndim()

long madness::BaseTensor::ndim ( ) const
inline

Returns the number of dimensions in the tensor.

References _ndim.

Referenced by madness::SRConf< T >::SRConf(), madness::TensorIterator< T, Q, R >::TensorIterator(), madness::TensorTrain< T >::TensorTrain(), madness::abs(), madness::BSHApply< T, NDIM >::add_coupling_and_levelshift(), madness::SRConf< T >::append(), madness::SeparatedConvolution< Q, NDIM >::apply(), madness::apply(), madness::arg(), madness::SRConf< T >::check_dimensions(), madness::SVDTensor< T >::compute_randomized_svd(), madness::RandomizedMatrixDecomposition< T >::compute_range(), madness::SVDTensor< T >::compute_svd(), madness::SVDTensor< T >::compute_svd_from_range(), madness::SVDTensor< T >::concatenate(), madness::Zcis::concatenate(), madness::conj(), madness::conj_transpose(), madness::convert(), madness::GenTensor< T >::convert_inplace(), madness::copy(), madness::SRConf< T >::copy_slice(), madness::Localizer::determine_frozen_orbitals(), madness::SRConf< T >::emul(), madness::TensorTrain< T >::emul(), madness::fast_transform(), fock_trace(), madness::TensorTrain< T >::gaxpy(), madness::SRConf< T >::general_transform(), madness::general_transform(), madness::SRConf< T >::get_configs(), madness::TensorTrain< T >::get_core(), madness::imag(), madness::inner(), madness::inner_result(), madness::SRConf< T >::inplace_add(), madness::SRConf< T >::is_flat(), madness::TensorTrain< T >::is_tensor(), madness::SRConf< T >::make_empty_vectors_and_weights(), madness::TensorTrain< T >::make_operator(), madness::SRConf< T >::make_slices(), madness::SRConf< T >::make_structure(), madness::TensorTrain< T >::make_tensor(), madness::SRConf< T >::make_vector_dimensions(), madness::SRConf< T >::make_vector_with_weights(), madness::GenTensor< T >::ndim(), madness::SRConf< T >::normf(), madness::guessfactory::ExopUnaryOpStructure::operator()(), function_real2complex_op< Q, NDIM >::operator()(), madness::abs_square_op< Q, NDIM >::operator()(), madness::real_op< Q, NDIM >::operator()(), madness::imag_op< Q, NDIM >::operator()(), madness::abs_op< Q, NDIM >::operator()(), madness::conj_op< Q, NDIM >::operator()(), madness::function_real2complex_op< Q, NDIM >::operator()(), madness::BinaryOpStructure< NDIM >::operator()(), madness::UnaryOpStructure< NDIM >::operator()(), madness::Tensor< T >::operator()(), madness::BSHApply< T, NDIM >::operator()(), madness::operator<<(), madness::SRConf< T >::operator=(), madness::TensorTrain< T >::operator=(), madness::outer(), madness::outer_result(), madness::FunctionImpl< T, NDIM >::partial_inner_contract(), madness::TensorTrain< T >::ranks(), madness::real(), madness::FunctionDefaults< NDIM >::recompute_cell_info(), madness::SRConf< T >::reconstruct(), madness::TensorTrain< T >::reconstruct(), madness::SRConf< T >::set_size_and_dim(), madness::TensorTrain< T >::splitdim(), madness::archive::ArchiveStoreImpl< Archive, Tensor< T > >::store(), madness::tensor_real2complex(), Test1(), madness::SRConf< T >::trace(), madness::TensorTrain< T >::trace(), madness::SRConf< T >::transform(), madness::transform(), madness::SRConf< T >::transform_dir(), madness::transform_dir(), madness::transpose(), madness::TensorTrain< T >::truncate(), madness::TensorTrain< T >::two_mode_representation(), madness::TensorTrain< T >::verify(), and madness::TensorTrain< T >::zero_me().

◆ reshape_inplace() [1/2]

void madness::BaseTensor::reshape_inplace ( const std::vector< long > &  d)
protected

Reshapes the tensor inplace.

Reshape the size and number of dimensions.

Modifies the current tensor to have the number and size of dimensions as described in the vector d . The total number of elements must be the same before and after, and the current tensor must be contiguous.

References d().

Referenced by madness::Tensor< T >::reshape().

◆ reshape_inplace() [2/2]

void madness::BaseTensor::reshape_inplace ( int  nd,
const long *  d 
)
protected

Reshapes the tensor inplace.

Reshape the size and number of dimensions.

Modifies the current tensor to have the number and size of dimensions as described in the vector d . The total number of elements must be the same before and after, and the current tensor must be contiguous.

References _size, d(), iscontiguous(), set_dims_and_size(), and TENSOR_ASSERT.

◆ set_dims_and_size()

void madness::BaseTensor::set_dims_and_size ( long  nd,
const long  d[] 
)
inlineprotected

◆ size()

long madness::BaseTensor::size ( ) const
inline

Returns the number of elements in the tensor.

References _size.

Referenced by madness::GenTensor< T >::GenTensor(), madness::SliceTensor< T >::SliceTensor(), madness::TensorTrain< T >::TensorTrain(), madness::SCF::analyze_vectors(), madness::apply(), madness::GTHPseudopotential< Q >::apply_potential(), madness::FunctionImpl< T, NDIM >::binaryXXa(), madness::LowRankFunction< T, NDIM, LDIM >::check_orthonormality(), madness::PNO::check_orthonormality(), madness::SRConf< T >::check_right_orthonormality(), madness::PNO::compute_cispd_correction_gs(), madness::PNO::compute_cispd_f12_correction_gs(), madness::Localizer::compute_core_valence_separation_transformation_matrix(), madness::OEP::compute_E_zeroth(), compute_frequencies(), madness::Zcis::compute_potentials(), madness::SVDTensor< T >::compute_randomized_svd(), madness::RandomizedMatrixDecomposition< T >::compute_range(), madness::Znemo::compute_residuals(), madness::SVDTensor< T >::compute_svd_from_range(), madness::SVDTensor< T >::concatenate(), madness::Zcis::concatenate(), madness::convert(), madness::GenTensor< T >::convert_inplace(), madness::copy(), madness::TensorTrain< T >::decompose(), madness::SCF::derivatives(), madness::Localizer::determine_frozen_orbitals(), madness::SCF::do_plots(), Ebar(), errsq(), madness::FunctionImpl< T, NDIM >::eval_plot_cube(), madness::XCfunctional::exc(), madness::exchange_anchor_test(), exchange_anchor_test(), madness::GFit< T, NDIM >::f12_fit(), madness::GFit< T, NDIM >::f12sq_fit(), madness::Tensor< T >::fillrandom(), madness::Localizer::find_degenerate_blocks(), fit(), madness::SRConf< T >::flat_vector_with_weights(), madness::FunctionImpl< T, NDIM >::gaxpy_ext_recursive(), madness::geqp3_result(), get_fock_transformation(), madness::SCF::get_fock_transformation(), madness::MolecularOrbitals< T, NDIM >::get_subset(), madness::SRConf< T >::has_data(), madness::Tensor< T >::has_data(), AtomicBasis::has_guess_info(), madness::AtomicBasis::has_guess_info(), madness::AtomicBasis::has_guesspsp_info(), madness::Nemo::hessian(), initial_guess(), madness::inner_result(), madness::OEP::iterate(), madness::archive::ArchiveLoadImpl< Archive, Tensor< T > >::load(), madness::Localizer::localize_boys(), madness::Localizer::localize_new(), madness::lq_result(), main(), madness::PNO::make_bsh_operators(), madness::Zcis::make_CIS_matrix(), madness::SCF::make_fock_matrix(), Fred::make_g(), madness::Zcis::make_guess(), madness::TDHF::make_guess_from_initial_diagonalization(), madness::TDHF::make_perturbed_fock_matrix(), make_pw_matrix(), madness::SRConf< T >::make_vector_with_weights(), madness::AtomicBasisSet::modify_dmat_psp(), madness::FunctionImpl< T, NDIM >::mulXXa(), madness::FunctionImpl< T, NDIM >::mulXXveca(), madness::Zcis::noct(), madness::norm2s_T(), madness::xc_lda_potential::operator()(), myunaryop_square< T, NDIM >::operator()(), madness::BSHApply< T, NDIM >::operator()(), BoysLocalization::operator()(), madness::operator<<(), madness::QuasiNewton::optimize(), optimize_coeffs(), madness::MolecularOptimizer::optimize_conjugate_gradients(), madness::MolecularOptimizer::optimize_quasi_newton(), madness::orgqr(), madness::ortho3(), madness::TDHF::orthonormalize(), madness::outer(), madness::FunctionImpl< T, NDIM >::partial_inner_contract(), madness::plotdx(), madness::OEP::print_orbens(), madness::FunctionImpl< T, NDIM >::print_stats(), madness::GFit< T, NDIM >::prune_small_coefficients(), madness::rank_revealing_decompose(), Plotter::read(), madness::InitParameters::readnw(), madness::MolecularOrbitals< T, NDIM >::save_restartaodata(), madness::Localizer::separate_core_valence(), madness::MolecularOrbitals< T, NDIM >::serialize(), madness::GenTensor< T >::size(), madness::Nemo::solve(), madness::Zcis::split(), madness::archive::ArchiveStoreImpl< Archive, Tensor< T > >::store(), madness::svd_result(), tensor2vec(), madness::guessfactory::tensor_to_coord(), Test1(), test_asymmetric(), test_hermiticity(), test_integration(), madness::test_inverse(), test_orthogonalization(), madness::FunctionImpl< T, NDIM >::unaryXXa(), madness::Localizer::undo_degenerate_rotations(), madness::SystolicMatrixAlgorithm< T >::unshuffle(), and Fred::value_and_gradient().

◆ splitdim_inplace()

void madness::BaseTensor::splitdim_inplace ( long  i,
long  dimi0,
long  dimi1 
)
protected

Splits dimension i.

Split dimension i in two ... the product of the new dimensions must match the old.

References _dim, _ndim, _stride, TENSOR_ASSERT, and TENSOR_MAXDIM.

Referenced by madness::Tensor< T >::splitdim(), and madness::TensorTrain< T >::splitdim().

◆ stride()

long madness::BaseTensor::stride ( int  i) const
inline

Returns the stride associated with dimension i.

References _stride.

Referenced by madness::TensorIterator< T, Q, R >::TensorIterator(), madness::inner_result(), and Test1().

◆ strides()

const long* madness::BaseTensor::strides ( ) const
inline

Returns the array of tensor strides.

References _stride.

◆ swapdim_inplace()

void madness::BaseTensor::swapdim_inplace ( long  i,
long  j 
)
protected

Swaps the dimensions.

Swap the dimensions i and j.

References _dim, _ndim, _stride, and TENSOR_ASSERT.

Referenced by madness::Tensor< T >::swapdim().

Member Data Documentation

◆ _dim

long madness::BaseTensor::_dim[TENSOR_MAXDIM]
protected

◆ _id

long madness::BaseTensor::_id
protected

◆ _ndim

long madness::BaseTensor::_ndim
protected

◆ _size

long madness::BaseTensor::_size
protected

◆ _stride

long madness::BaseTensor::_stride[TENSOR_MAXDIM]
protected

The documentation for this class was generated from the following files: