29 #ifndef _mpqc_src_lib_chemistry_qc_df_df_h
30 #define _mpqc_src_lib_chemistry_qc_df_df_h
32 #include <util/state/state.h>
33 #include <util/misc/compute.h>
34 #include <util/group/memory.h>
35 #include <math/scmat/result.h>
36 #include <chemistry/qc/basis/integral.h>
37 #include <chemistry/qc/lcao/tbint_runtime.h>
64 SolveMethod_InverseCholesky = 0,
65 SolveMethod_Cholesky = 1,
66 SolveMethod_RefinedCholesky = 2,
67 SolveMethod_InverseBunchKaufman = 3,
68 SolveMethod_BunchKaufman = 4,
69 SolveMethod_RefinedBunchKaufman = 5
89 const std::string& kernel_key,
100 return runtime()->factory()->integral();
105 const std::string& kernel_key()
const {
return kernel_key_; }
106 SolveMethod solver()
const {
return solver_; }
114 virtual void compute();
132 std::string kernel_key_;
169 const std::string& kernel_key,
197 const std::string& kernel_key,
218 #endif // end of header guard
TwoBodyMOIntsRuntimeUnion23 packages 2-center and 3-center runtimes; it also keeps track of 2-center ...
Definition: tbint_runtime.h:523
const Ref< DistArray4 > & conjugateC() const
returns the conjugate expansion matrix defined as .
Definition: df.h:155
const RefSymmSCMatrix & kernel() const
returns the kernel matrix in the fitting basis, i.e. .
Definition: df.h:147
The RefSymmSCMatrix class is a smart pointer to an SCSymmSCMatrix specialization.
Definition: matrix.h:265
A template class that maintains references counts.
Definition: ref.h:361
Decomposition by density fitting with respect to some kernel.
Definition: df.h:61
PermutedDensityFitting(const Ref< DensityFitting > &df21)
compute density fitting for |mo1 mo2) from |mo1 ao2)
The RefSCDimension class is a smart pointer to an SCDimension specialization.
Definition: dim.h:152
This class is used to contain information about classes.
Definition: class.h:147
DensityFitting(const Ref< MOIntsRuntime > &rtime, const std::string &kernel_key, SolveMethod solver, const Ref< OrbitalSpace > &space1, const Ref< OrbitalSpace > &space2, const Ref< GaussianBasisSet > &fitting_basis, bool ri=false)
Determines density fitting of product (space1 space2) in terms of fitting basis.
const Ref< DistArray4 > & C() const
returns the fitting coeffcients
Definition: df.h:108
void save_data_state(StateOut &)
Save the base classes (with save_data_state) and the members in the same order that the StateIn CTOR ...
Definition: stateout.h:71
RefSCDimension product_dimension() const
produces RefSCDimension that corresponds to the product space
type
types of known two-body operators
Definition: operator.h:318
Base class for objects that can save/restore state.
Definition: state.h:45
Computes density fitting for |ij) density from fitting of |ji) DensityFitting.
Definition: df.h:190
static int definite_kernel(TwoBodyOper::type kernel_type, const Ref< IntParams > &kernel_params)
Test whether a density-fitting kernel has definite sign.
Contains all MPQC code up to version 3.
Definition: mpqcin.h:14
void save_data_state(StateOut &)
Save the base classes (with save_data_state) and the members in the same order that the StateIn CTOR ...
Generated at Sun Jan 26 2020 23:23:58 for MPQC
3.0.0-alpha using the documentation package Doxygen
1.8.16.