MPQC 3.0.0-alpha
Loading...
Searching...
No Matches
df.h
1//
2// df.h
3//
4// Copyright (C) 2009 Edward Valeev
5//
6// Author: Edward Valeev <evaleev@vt.edu>
7// Maintainer: EV
8//
9// This file is part of the SC Toolkit.
10//
11// The SC Toolkit is free software; you can redistribute it and/or modify
12// it under the terms of the GNU Library General Public License as published by
13// the Free Software Foundation; either version 2, or (at your option)
14// any later version.
15//
16// The SC Toolkit is distributed in the hope that it will be useful,
17// but WITHOUT ANY WARRANTY; without even the implied warranty of
18// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19// GNU Library General Public License for more details.
20//
21// You should have received a copy of the GNU Library General Public License
22// along with the SC Toolkit; see the file COPYING.LIB. If not, write to
23// the Free Software Foundation, 675 Mass Ave, Cambridge, MA 02139, USA.
24//
25// The U.S. Government is granted a limited license as per AL 91-7.
26//
27
28
29#ifndef _mpqc_src_lib_chemistry_qc_df_df_h
30#define _mpqc_src_lib_chemistry_qc_df_df_h
31
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>
38
39namespace sc {
40
41 class OrbitalSpace;
42 class DistArray4;
43
61 class DensityFitting: virtual public SavableState {
62 public:
63 enum SolveMethod {
64 SolveMethod_InverseCholesky = 0,
65 SolveMethod_Cholesky = 1,
66 SolveMethod_RefinedCholesky = 2,
67 SolveMethod_InverseBunchKaufman = 3,
68 SolveMethod_BunchKaufman = 4,
69 SolveMethod_RefinedBunchKaufman = 5
70 };
71
73
89 const std::string& kernel_key,
90 SolveMethod solver,
91 const Ref<OrbitalSpace>& space1,
92 const Ref<OrbitalSpace>& space2,
93 const Ref<GaussianBasisSet>& fitting_basis,
94 bool ri = false);
97
98 const Ref<MOIntsRuntime>& runtime() const { return runtime_; }
99 const Ref<Integral>& integral() const {
100 return runtime()->factory()->integral();
101 }
102 const Ref<OrbitalSpace>& space1() const { return space1_; }
103 const Ref<OrbitalSpace>& space2() const { return space2_; }
104 const Ref<GaussianBasisSet>& fbasis() const { return fbasis_; }
105 const std::string& kernel_key() const { return kernel_key_; }
106 SolveMethod solver() const { return solver_; }
108 const Ref<DistArray4>& C() const {
109 return C_;
110 }
113
114 virtual void compute();
115
123 static int definite_kernel(TwoBodyOper::type kernel_type, const Ref<IntParams>& kernel_params);
124 static int definite_kernel(TwoBodyOper::type kernel_type, const std::string& kernel_params_key);
125
126 private:
127
128 Ref<MOIntsRuntime> runtime_;
129 Ref<GaussianBasisSet> fbasis_;
130 Ref<OrbitalSpace> space1_;
131 Ref<OrbitalSpace> space2_;
132 std::string kernel_key_;
133 SolveMethod solver_;
134 bool ri_;
135
136 bool evaluated_;
137
138 static ClassDesc class_desc_;
139
140 protected:
141
143 RefSymmSCMatrix kernel_;
144 Ref<DistArray4> cC_;
145
147 const RefSymmSCMatrix& kernel() const {
148 return kernel_;
149 }
156 return cC_;
157 }
158
159 };
160
163 public:
165
169 const std::string& kernel_key,
170 SolveMethod solver,
171 const Ref<OrbitalSpace>& space1,
172 const Ref<OrbitalSpace>& space2,
173 const Ref<GaussianBasisSet>& fitting_basis,
174 const Ref<DistArray4>& mo1_ao2_df);
177
178 private:
179
180 Ref<DistArray4> mo1_ao2_df_;
181
182 // overloads DensityFitting::compute()
183 void compute();
184
185 static ClassDesc class_desc_;
186
187 };
188
191 public:
197 const std::string& kernel_key,
198 SolveMethod solver,
199 const Ref<OrbitalSpace>& space1,
200 const Ref<OrbitalSpace>& space2,
201 const Ref<GaussianBasisSet>& fitting_basis,
202 const Ref<DistArray4>& df21);
205
206 private:
207 Ref<DistArray4> df21_;
208
209 // overloads DensityFitting::compute()
210 void compute();
211
212 static ClassDesc class_desc_;
213
214 };
215
216} // end of namespace sc
217
218#endif // end of header guard
219
220// Local Variables:
221// mode: c++
222// c-file-style: "CLJ-CONDENSED"
223// End:
This class is used to contain information about classes.
Definition class.h:147
Decomposition by density fitting with respect to some kernel.
Definition df.h:61
static int definite_kernel(TwoBodyOper::type kernel_type, const Ref< IntParams > &kernel_params)
Test whether a density-fitting kernel has definite sign.
const Ref< DistArray4 > & conjugateC() const
returns the conjugate expansion matrix defined as .
Definition df.h:155
void save_data_state(StateOut &)
Save the base classes (with save_data_state) and the members in the same order that the StateIn CTOR ...
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.
RefSCDimension product_dimension() const
produces RefSCDimension that corresponds to the product space
const Ref< DistArray4 > & C() const
returns the fitting coeffcients
Definition df.h:108
const RefSymmSCMatrix & kernel() const
returns the kernel matrix in the fitting basis, i.e. .
Definition df.h:147
Computes density fitting for |ij) density from fitting of |ji) DensityFitting.
Definition df.h:190
void save_data_state(StateOut &)
Save the base classes (with save_data_state) and the members in the same order that the StateIn CTOR ...
PermutedDensityFitting(const Ref< MOIntsRuntime > &rtime, const std::string &kernel_key, SolveMethod solver, const Ref< OrbitalSpace > &space1, const Ref< OrbitalSpace > &space2, const Ref< GaussianBasisSet > &fitting_basis, const Ref< DistArray4 > &df21)
compute density fitting for |mo1 mo2) from |mo1 ao2)
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
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
Base class for objects that can save/restore state.
Definition state.h:45
Restores fundamental and user-defined types from images created with StateOut.
Definition statein.h:79
Serializes fundamental and user-defined types.
Definition stateout.h:71
Computes density fitting for |ij) density from fitting of |iq) DensityFitting where q is the AO space...
Definition df.h:162
void save_data_state(StateOut &)
Save the base classes (with save_data_state) and the members in the same order that the StateIn CTOR ...
TransformedDensityFitting(const Ref< MOIntsRuntime > &rtime, const std::string &kernel_key, SolveMethod solver, const Ref< OrbitalSpace > &space1, const Ref< OrbitalSpace > &space2, const Ref< GaussianBasisSet > &fitting_basis, const Ref< DistArray4 > &mo1_ao2_df)
compute density fitting for |mo1 mo2) from |mo1 ao2)
TwoBodyMOIntsRuntimeUnion23 packages 2-center and 3-center runtimes; it also keeps track of 2-center ...
Definition tbint_runtime.h:523
Contains all MPQC code up to version 3.
Definition mpqcin.h:14
type
types of known two-body operators
Definition operator.h:318

Generated at Wed Sep 25 2024 02:45:29 for MPQC 3.0.0-alpha using the documentation package Doxygen 1.12.0.