mpqc::lcao::pbc::FockBuilderImpl< Tile, Policy, J_Builder, K_Builder, MA_Builder > Class Template Reference
Collaboration diagram for mpqc::lcao::pbc::FockBuilderImpl< Tile, Policy, J_Builder, K_Builder, MA_Builder >:

Documentation

template<typename Tile, typename Policy, typename J_Builder = DirectFockBuilder<Tile, Policy>, typename K_Builder = DirectFockBuilder<Tile, Policy>, typename MA_Builder = ::mpqc::pbc::ma::MA_CFFLatticeSum< pbc::gaussian::AOFactory<Tile, Policy>>>
class mpqc::lcao::pbc::FockBuilderImpl< Tile, Policy, J_Builder, K_Builder, MA_Builder >

pbc::MaJFockBuilder is an implementation of pbc::FockBuilder with multipole-accelerated 4-center Coulomb and 4-center exchange. For Coulomb, multipole approximation is used in Crystal Far Field. The builders for near-field Coulomb and exchange are customizable.

Template Parameters
Tile
Policy
J_Builderthe builder for near-field Coulomb; the default is pbc::DirectFockBuilder<Tile, Policy> which is a conventional (4-center-integral based) builder.
K_Builderthe builder for the exchange; the default is pbc::DirectFockBuilder<Tile, Policy> which is a conventional (4-center-integral based) builder.

Public Types

using factory_type = pbc::gaussian::AOFactory< Tile, Policy >
 
using array_type = typename factory_type::TArray
 
using result_type = typename FockBuilder< Tile, Policy >::result_type
 
- Public Types inherited from mpqc::lcao::pbc::FockBuilder< Tile, Policy >
using array_type = TA::DistArray< Tile, Policy >
 
using result_type = std::tuple< array_type, std::optional< double > >
 

Public Member Functions

 FockBuilderImpl (factory_type &factory)
 
 ~FockBuilderImpl ()
 
result_type operator() (array_type const &D, double target_precision, bool is_density_diagonal) override
 This computes the 2-e part of the closed-shell Fock matrix in periodic HF. More...
 
void register_fock (const array_type &fock, FormulaRegistry< array_type > &registry) const override
 
Vector3i fock_lattice_range () const override
 This returns the lattice range of Fock representation F(μ_0, ν_R) More...
 
const MA_Builder * multipole_builder () const
 
const J_Builder * coulomb_builder () const
 
const K_Builder * exchange_builder () const
 
- Public Member Functions inherited from mpqc::lcao::pbc::FockBuilder< Tile, Policy >
virtual ~FockBuilder ()
 

Static Public Attributes

static constexpr const bool have_ma = !std::is_same_v<MA_Builder, void>
 

Member Typedef Documentation

◆ array_type

template<typename Tile , typename Policy , typename J_Builder = DirectFockBuilder<Tile, Policy>, typename K_Builder = DirectFockBuilder<Tile, Policy>, typename MA_Builder = ::mpqc::pbc::ma::MA_CFFLatticeSum< pbc::gaussian::AOFactory<Tile, Policy>>>
using mpqc::lcao::pbc::FockBuilderImpl< Tile, Policy, J_Builder, K_Builder, MA_Builder >::array_type = typename factory_type::TArray

◆ factory_type

template<typename Tile , typename Policy , typename J_Builder = DirectFockBuilder<Tile, Policy>, typename K_Builder = DirectFockBuilder<Tile, Policy>, typename MA_Builder = ::mpqc::pbc::ma::MA_CFFLatticeSum< pbc::gaussian::AOFactory<Tile, Policy>>>
using mpqc::lcao::pbc::FockBuilderImpl< Tile, Policy, J_Builder, K_Builder, MA_Builder >::factory_type = pbc::gaussian::AOFactory<Tile, Policy>

◆ result_type

template<typename Tile , typename Policy , typename J_Builder = DirectFockBuilder<Tile, Policy>, typename K_Builder = DirectFockBuilder<Tile, Policy>, typename MA_Builder = ::mpqc::pbc::ma::MA_CFFLatticeSum< pbc::gaussian::AOFactory<Tile, Policy>>>
using mpqc::lcao::pbc::FockBuilderImpl< Tile, Policy, J_Builder, K_Builder, MA_Builder >::result_type = typename FockBuilder<Tile, Policy>::result_type

Constructor & Destructor Documentation

◆ FockBuilderImpl()

template<typename Tile , typename Policy , typename J_Builder = DirectFockBuilder<Tile, Policy>, typename K_Builder = DirectFockBuilder<Tile, Policy>, typename MA_Builder = ::mpqc::pbc::ma::MA_CFFLatticeSum< pbc::gaussian::AOFactory<Tile, Policy>>>
mpqc::lcao::pbc::FockBuilderImpl< Tile, Policy, J_Builder, K_Builder, MA_Builder >::FockBuilderImpl ( factory_type factory)
inline

◆ ~FockBuilderImpl()

template<typename Tile , typename Policy , typename J_Builder = DirectFockBuilder<Tile, Policy>, typename K_Builder = DirectFockBuilder<Tile, Policy>, typename MA_Builder = ::mpqc::pbc::ma::MA_CFFLatticeSum< pbc::gaussian::AOFactory<Tile, Policy>>>
mpqc::lcao::pbc::FockBuilderImpl< Tile, Policy, J_Builder, K_Builder, MA_Builder >::~FockBuilderImpl ( )
inline

Member Function Documentation

◆ coulomb_builder()

template<typename Tile , typename Policy , typename J_Builder = DirectFockBuilder<Tile, Policy>, typename K_Builder = DirectFockBuilder<Tile, Policy>, typename MA_Builder = ::mpqc::pbc::ma::MA_CFFLatticeSum< pbc::gaussian::AOFactory<Tile, Policy>>>
const J_Builder* mpqc::lcao::pbc::FockBuilderImpl< Tile, Policy, J_Builder, K_Builder, MA_Builder >::coulomb_builder ( ) const
inline

◆ exchange_builder()

template<typename Tile , typename Policy , typename J_Builder = DirectFockBuilder<Tile, Policy>, typename K_Builder = DirectFockBuilder<Tile, Policy>, typename MA_Builder = ::mpqc::pbc::ma::MA_CFFLatticeSum< pbc::gaussian::AOFactory<Tile, Policy>>>
const K_Builder* mpqc::lcao::pbc::FockBuilderImpl< Tile, Policy, J_Builder, K_Builder, MA_Builder >::exchange_builder ( ) const
inline

◆ fock_lattice_range()

template<typename Tile , typename Policy , typename J_Builder = DirectFockBuilder<Tile, Policy>, typename K_Builder = DirectFockBuilder<Tile, Policy>, typename MA_Builder = ::mpqc::pbc::ma::MA_CFFLatticeSum< pbc::gaussian::AOFactory<Tile, Policy>>>
Vector3i mpqc::lcao::pbc::FockBuilderImpl< Tile, Policy, J_Builder, K_Builder, MA_Builder >::fock_lattice_range ( ) const
inlineoverridevirtual

This returns the lattice range of Fock representation F(μ_0, ν_R)

Returns
a vector of number of unit cells included in each positive direction

Implements mpqc::lcao::pbc::FockBuilder< Tile, Policy >.

◆ multipole_builder()

template<typename Tile , typename Policy , typename J_Builder = DirectFockBuilder<Tile, Policy>, typename K_Builder = DirectFockBuilder<Tile, Policy>, typename MA_Builder = ::mpqc::pbc::ma::MA_CFFLatticeSum< pbc::gaussian::AOFactory<Tile, Policy>>>
const MA_Builder* mpqc::lcao::pbc::FockBuilderImpl< Tile, Policy, J_Builder, K_Builder, MA_Builder >::multipole_builder ( ) const
inline

◆ operator()()

template<typename Tile , typename Policy , typename J_Builder = DirectFockBuilder<Tile, Policy>, typename K_Builder = DirectFockBuilder<Tile, Policy>, typename MA_Builder = ::mpqc::pbc::ma::MA_CFFLatticeSum< pbc::gaussian::AOFactory<Tile, Policy>>>
result_type mpqc::lcao::pbc::FockBuilderImpl< Tile, Policy, J_Builder, K_Builder, MA_Builder >::operator() ( array_type const &  D,
double  target_precision,
bool  is_density_diagonal 
)
inlineoverridevirtual

This computes the 2-e part of the closed-shell Fock matrix in periodic HF.

Parameters
Dthe (1-particle) density matrix D(μ_0, ν_Rd); this should be single-occupancy, i.e. integrate to half the number of electrons
target_precisionthe (absolute) target precision of the Fock matrix; this precision is not in general guaranteed, but the actual precision of the Fock matrix is usually proportional to this value.
is_density_diagonaltrue if D only has a D(μ_0, ν_0) block
Returns
Fock matrix + an optional energy correction

Implements mpqc::lcao::pbc::FockBuilder< Tile, Policy >.

◆ register_fock()

template<typename Tile , typename Policy , typename J_Builder = DirectFockBuilder<Tile, Policy>, typename K_Builder = DirectFockBuilder<Tile, Policy>, typename MA_Builder = ::mpqc::pbc::ma::MA_CFFLatticeSum< pbc::gaussian::AOFactory<Tile, Policy>>>
void mpqc::lcao::pbc::FockBuilderImpl< Tile, Policy, J_Builder, K_Builder, MA_Builder >::register_fock ( const array_type fock,
FormulaRegistry< array_type > &  registry 
) const
inlineoverridevirtual

Member Data Documentation

◆ have_ma

template<typename Tile , typename Policy , typename J_Builder = DirectFockBuilder<Tile, Policy>, typename K_Builder = DirectFockBuilder<Tile, Policy>, typename MA_Builder = ::mpqc::pbc::ma::MA_CFFLatticeSum< pbc::gaussian::AOFactory<Tile, Policy>>>
constexpr const bool mpqc::lcao::pbc::FockBuilderImpl< Tile, Policy, J_Builder, K_Builder, MA_Builder >::have_ma = !std::is_same_v<MA_Builder, void>
staticconstexpr

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