MPQC
3.0.0-alpha
|
The MBPT2 class implements several second-order perturbation theory methods. More...
#include <chemistry/qc/mbpt/mbpt.h>
Public Member Functions | |
MBPT2 (StateIn &) | |
MBPT2 (const Ref< KeyVal > &) | |
The KeyVal constructor. More... | |
void | save_data_state (StateOut &) |
Save the base classes (with save_data_state) and the members in the same order that the StateIn CTOR initializes them. More... | |
Ref< SCF > | ref () |
double | ref_energy () |
double | corr_energy () |
RefSCVector | ref_energy_gradient () |
RefSCVector | corr_energy_gradient () |
int | nelectron () |
Returns the number of electrons. | |
int | nfzcore () const |
int | nfzvirt () const |
RefSymmSCMatrix | density () |
Returns the SO density. | |
double | magnetic_moment () const |
Computes the S (or J) magnetic moment of the target state(s), in units of . More... | |
bool | analytic_gradient_implemented () const |
must overload this in a derived class if analytic gradient can be computed More... | |
int | value_implemented () const |
void | set_desired_value_accuracy (double) |
set the value accuracy | |
void | symmetry_changed () |
Call this if you have changed the molecular symmetry of the molecule contained by this MolecularEnergy. | |
void | obsolete () |
Marks all results as being out of date. More... | |
void | print (std::ostream &o=ExEnv::out0()) const |
Print information about the object. | |
Public Member Functions inherited from sc::Wavefunction | |
Wavefunction (StateIn &) | |
Wavefunction (const Ref< KeyVal > &) | |
The KeyVal constructor. More... | |
void | save_data_state (StateOut &) |
Save the base classes (with save_data_state) and the members in the same order that the StateIn CTOR initializes them. More... | |
double | density (const SCVector3 &) |
double | density_gradient (const SCVector3 &, double *) |
double | natural_orbital (const SCVector3 &r, int iorb) |
double | natural_orbital_density (const SCVector3 &r, int orb, double *orbval=0) |
double | orbital (const SCVector3 &r, int iorb, const RefSCMatrix &orbs) |
void | orbitals (const SCVector3 &r, const RefSCMatrix &orbs, RefSCVector &values) |
double | orbital_density (const SCVector3 &r, int iorb, const RefSCMatrix &orbs, double *orbval=0) |
double | total_charge () const |
Returns the total charge of the system. | |
virtual RefSymmSCMatrix | ao_density () |
Returns the AO density. | |
virtual RefSCMatrix | natural_orbitals () |
Returns the natural orbitals, in SO basis. | |
virtual RefDiagSCMatrix | natural_density () |
Returns the natural density (a diagonal matrix). | |
int | spin_polarized () |
Return 1 if the magnetic moment != 0. | |
int | dk () const |
Returns the level the of the Douglas-Kroll approximation. | |
virtual RefSymmSCMatrix | alpha_density () |
Return alpha electron densities in the SO basis. | |
virtual RefSymmSCMatrix | beta_density () |
Return beta electron densities in the SO basis. | |
virtual RefSymmSCMatrix | alpha_ao_density () |
Return alpha electron densities in the AO basis. | |
virtual RefSymmSCMatrix | beta_ao_density () |
Return beta electron densities in the AO basis. | |
virtual RefSCMatrix | nao (double *atom_charges=0) |
returns the ao to nao transformation matrix | |
virtual RefSymmSCMatrix | overlap () |
Returns the SO overlap matrix. | |
virtual RefSymmSCMatrix | core_hamiltonian_for_basis (const Ref< GaussianBasisSet > &bas, const Ref< GaussianBasisSet > &pbas=0) |
Returns the SO core Hamiltonian in the given basis and momentum basis. More... | |
virtual RefSymmSCMatrix | core_hamiltonian () |
Returns the SO core Hamiltonian. | |
RefSymmSCMatrix | core_hamiltonian_nr (const Ref< GaussianBasisSet > &bas) |
virtual double | nuclear_repulsion_energy () |
Returns the nuclear repulsion energy. More... | |
void | nuclear_repulsion_energy_gradient (double *g) |
Computes the nuclear repulsion gradient. More... | |
virtual void | nuclear_repulsion_energy_gradient (double **g) |
Computes the nuclear repulsion gradient. More... | |
RefSCDimension | ao_dimension () |
Atomic orbital dimension. | |
RefSCDimension | so_dimension () |
Symmetry adapted orbital dimension. | |
RefSCDimension | oso_dimension () |
Orthogonalized symmetry adapted orbital dimension. | |
Ref< SCMatrixKit > | basis_matrixkit () |
Matrix kit for AO, SO, orthogonalized SO, and MO dimensioned matrices. | |
Ref< Molecule > | molecule () const |
Returns the Molecule. | |
Ref< GaussianBasisSet > | basis () const |
Returns the basis set. | |
Ref< GaussianBasisSet > | momentum_basis () const |
Returns the basis used for p^2 in the DK correction. | |
Ref< GaussianBasisSet > | atom_basis () const |
Returns the basis set describing the nuclear charge distributions. | |
const double * | atom_basis_coef () const |
Returns the coefficients of the nuclear charge distribution basis functions. | |
Ref< Integral > | integral () |
Returns the integral evaluator. | |
void | symmetry_changed () |
Call this if you have changed the molecular symmetry of the molecule contained by this MolecularEnergy. | |
RefSCMatrix | so_to_orthog_so () |
Returns a matrix which does the default transform from SO's to orthogonal SO's. More... | |
RefSCMatrix | so_to_orthog_so_inverse () |
Returns the inverse of the transformation returned by so_to_orthog_so. | |
OverlapOrthog::OrthogMethod | orthog_method () const |
Returns the orthogonalization method. | |
virtual void | set_orthog_method (const OverlapOrthog::OrthogMethod &) |
(Re)Sets the orthogonalization method and makes this obsolete. More... | |
double | lindep_tol () const |
Returns the tolerance for linear dependencies. | |
void | set_lindep_tol (double) |
Re(Sets) the tolerance for linear dependencies. | |
void | obsolete () |
Marks all results as being out of date. More... | |
void | print (std::ostream &=ExEnv::out0()) const |
Print information about the object. | |
void | writeorbitals () |
output orbitals to some files to facilitate plotting, with the help of the WriteOrbital class. | |
Public Member Functions inherited from sc::MolecularEnergy | |
MolecularEnergy (const MolecularEnergy &) | |
MolecularEnergy (const Ref< KeyVal > &) | |
The KeyVal constructor. More... | |
MolecularEnergy (StateIn &) | |
void | set_checkpoint () |
Set up checkpointing. | |
void | set_checkpoint_file (const char *) |
void | set_checkpoint_freq (int freq) |
bool | if_to_checkpoint () const |
Check if need to checkpoint. | |
const char * | checkpoint_file () const |
int | checkpoint_freq () const |
MolecularEnergy & | operator= (const MolecularEnergy &) |
virtual double | energy () |
A wrapper around value();. | |
virtual RefSCDimension | moldim () const |
void | guess_hessian (RefSymmSCMatrix &) |
Compute a quick, approximate hessian. | |
RefSymmSCMatrix | inverse_hessian (RefSymmSCMatrix &) |
int | gradient_implemented () const |
Reports whether gradient is implemented either analytically or using MolecularGradient object. More... | |
int | hessian_implemented () const |
Reports whether hessian is implemented either analytically or using MolecularHessian object. More... | |
void | set_desired_gradient_accuracy (double acc) |
These functions overload their Function counterparts. More... | |
void | set_desired_hessian_accuracy (double acc) |
void | set_molhess (const Ref< MolecularHessian > &molhess) |
Use this function to provide MolecularHessian object that will be used to compute hessian. More... | |
const Ref< MolecularHessian > & | molhess () const |
RefSymmSCMatrix | hessian () |
Will throw if hessian_implemented() returns 0. | |
void | set_molgrad (const Ref< MolecularGradient > &molgrad) |
Use this function to provide MolecularGradient object that will be used to compute gradient. More... | |
const Ref< MolecularGradient > & | molgrad () const |
RefSCVector | gradient () |
Will throw if gradient_implemented() returns 0. | |
void | set_x (const RefSCVector &) |
Set and retrieve the coordinate values. | |
RefSCVector | get_cartesian_x () |
Return the cartesian coordinates. | |
RefSCVector | get_cartesian_gradient () |
Return the cartesian gradient. | |
RefSymmSCMatrix | get_cartesian_hessian () |
Return the cartesian hessian. | |
Ref< MolecularCoor > | molecularcoor () |
Ref< NonlinearTransform > | change_coordinates () |
An optimizer can call change coordinates periodically to give the function an opportunity to change its coordinate system. More... | |
virtual void | purge () |
This function purges any caches of data in MolecularEnergy. More... | |
const RefSCVector & | electric_field () const |
returns the electric field vector | |
void | print_natom_3 (const RefSCVector &, const char *t=0, std::ostream &o=ExEnv::out0()) const |
Nicely print n x 3 data that are stored in a vector. | |
void | print_natom_3 (double **, const char *t=0, std::ostream &o=ExEnv::out0()) const |
void | print_natom_3 (double *, const char *t=0, std::ostream &o=ExEnv::out0()) const |
Public Member Functions inherited from sc::Function | |
int | gradient_needed () const |
int | do_gradient (int) |
virtual double | actual_gradient_accuracy () const |
virtual double | desired_gradient_accuracy () const |
AccResultRefSCVector & | gradient_result () |
int | hessian_needed () const |
int | do_hessian (int) |
virtual double | actual_hessian_accuracy () const |
virtual double | desired_hessian_accuracy () const |
AccResultRefSymmSCMatrix & | hessian_result () |
virtual bool | desired_value_accuracy_set_to_default () const |
virtual bool | desired_gradient_accuracy_set_to_default () const |
virtual bool | desired_hessian_accuracy_set_to_default () const |
RefSCVector | get_x () const |
const RefSCVector & | get_x_no_copy () const |
void | print_desired_accuracy (std::ostream &=ExEnv::out0()) const |
similar to print(), but only prins desired accuracies | |
virtual bool | throw_if_tolerance_exceeded () const |
Overridden Compute member. | |
Function () | |
Function (StateIn &) | |
Function (const Function &) | |
Function (const Ref< KeyVal > &, double funcacc=DBL_EPSILON, double gradacc=DBL_EPSILON, double hessacc=DBL_EPSILON) | |
The keyval constructor reads the following keywords: More... | |
virtual | ~Function () |
Function & | operator= (const Function &) |
Ref< SCMatrixKit > | matrixkit () const |
Return the SCMatrixKit used to construct vectors and matrices. | |
RefSCDimension | dimension () const |
Return the SCDimension of the problem. | |
virtual double | value () |
Return the value of the function. | |
int | value_needed () const |
Returns nonzero if the current value is not up-to-date. | |
int | do_value (int) |
If passed a nonzero number, compute the value the next time compute() is called. More... | |
AccResultdouble & | value_result () |
virtual double | actual_value_accuracy () const |
Return the accuracy with which the value has been computed. | |
virtual double | desired_value_accuracy () const |
Return the accuracy with which the value is to be computed. | |
Public Member Functions inherited from sc::SavableState | |
SavableState & | operator= (const SavableState &) |
void | save_state (StateOut &) |
Save the state of the object as specified by the StateOut object. More... | |
void | save_object_state (StateOut &) |
This can be used for saving state when the exact type of the object is known for both the save and the restore. More... | |
virtual void | save_vbase_state (StateOut &) |
Save the virtual bases for the object. More... | |
Public Member Functions inherited from sc::DescribedClass | |
DescribedClass (const DescribedClass &) | |
DescribedClass & | operator= (const DescribedClass &) |
ClassDesc * | class_desc () const MPQC__NOEXCEPT |
This returns the unique pointer to the ClassDesc corresponding to the given type_info object. More... | |
const char * | class_name () const |
Return the name of the object's exact type. | |
int | class_version () const |
Return the version of the class. | |
Ref< DescribedClass > | ref () |
Return this object wrapped up in a Ref smart pointer. More... | |
Public Member Functions inherited from sc::RefCount | |
size_t | identifier () const |
Return the unique identifier for this object that can be compared for different objects of different types. More... | |
int | lock_ptr () const |
Lock this object. | |
int | unlock_ptr () const |
Unlock this object. | |
void | use_locks (bool inVal) |
start and stop using locks on this object | |
refcount_t | nreference () const |
Return the reference count. | |
refcount_t | reference () |
Increment the reference count and return the new count. | |
refcount_t | dereference () |
Decrement the reference count and return the new count. | |
int | managed () const |
void | unmanage () |
Turn off the reference counting mechanism for this object. More... | |
Protected Member Functions | |
void | init_variables () |
void | compute () |
Recompute at least the results that have compute true and are not already computed. More... | |
void | eigen (RefDiagSCMatrix &vals, RefSCMatrix &vecs, RefDiagSCMatrix &occs) |
void | compute_hsos_v1 () |
distsize_t | compute_v2_memory (int ni, int nfuncmax, int nbfme, int nshell, int ndocc, int nsocc, int nvir, int nproc) |
void | compute_hsos_v2 () |
void | compute_hsos_v2_lb () |
int | compute_cs_batchsize (size_t mem_static, int nocc_act) |
distsize_t | compute_cs_dynamic_memory (int ni, int nocc_act) |
int | make_cs_gmat (RefSymmSCMatrix &Gmat, double *DPmat) |
int | make_cs_gmat_new (RefSymmSCMatrix &Gmat, const RefSymmSCMatrix &DPmat) |
void | form_max_dens (double *DPmat, signed char *maxp) |
int | init_cs_gmat () |
void | done_cs_gmat () |
int | make_g_d_nor (RefSymmSCMatrix &Gmat, double *DPmat, const double *mgdbuff) |
void | cs_cphf (double **scf_vector, double *Laj, double *eigval, RefSCMatrix &P2aj) |
void | s2pdm_contrib (const double *intderbuf, double *PHF, double *P2AO, double **hf_ginter, double **ginter) |
void | hcore_cs_grad (double *PHF, double *PMP2, double **hf_ginter, double **ginter) |
void | overlap_cs_grad (double *WHF, double *WMP2, double **hf_ginter, double **ginter) |
void | compute_cs_grad () |
Protected Member Functions inherited from sc::Wavefunction | |
double | min_orthog_res () |
double | max_orthog_res () |
void | copy_orthog_info (const Ref< Wavefunction > &) |
Protected Member Functions inherited from sc::MolecularEnergy | |
void | failure (const char *) |
virtual void | set_energy (double) |
This is just a wrapper around set_value(). | |
virtual void | set_gradient (RefSCVector &) |
These are passed gradients and hessian in cartesian coordinates. More... | |
virtual void | set_hessian (RefSymmSCMatrix &) |
void | x_to_molecule () |
void | molecule_to_x () |
virtual bool | analytic_hessian_implemented () const |
must overload this in a derived class if analytic hessian can be computed More... | |
Protected Member Functions inherited from sc::Function | |
virtual void | set_value (double) |
virtual void | set_matrixkit (const Ref< SCMatrixKit > &) |
Set the SCMatrixKit that should be used to construct the requisite vectors and matrices. | |
virtual void | set_dimension (const RefSCDimension &) |
virtual void | set_actual_value_accuracy (double) |
virtual void | set_actual_gradient_accuracy (double) |
virtual void | set_actual_hessian_accuracy (double) |
RefSCVector & | get_x_reference () |
Get read/write access to the coordinates for modification. | |
void | do_change_coordinates (const Ref< NonlinearTransform > &) |
Change the coordinate system and apply the given transform to intermediates matrices and vectors. | |
Protected Member Functions inherited from sc::SavableState | |
SavableState (const SavableState &) | |
SavableState (StateIn &) | |
Each derived class StateIn CTOR handles the restore corresponding to calling save_object_state, save_vbase_state, and save_data_state listed above. More... | |
Protected Member Functions inherited from sc::RefCount | |
RefCount (const RefCount &) | |
RefCount & | operator= (const RefCount &) |
Static Protected Member Functions | |
static double | ref_to_mp2_acc () |
Protected Attributes | |
Ref< SCF > | reference_ |
Ref< MemoryGrp > | mem |
int | nfzc |
int | nfzv |
size_t | mem_alloc |
double | cphf_epsilon_ |
int | eliminate_in_gmat_ |
const double * | intbuf_ |
Ref< TwoBodyInt > | tbint_ |
Ref< TwoBodyInt > * | tbints_ |
Ref< TwoBodyDerivInt > * | tbintder_ |
int | nbasis |
int | noso |
Ref< MessageGrp > | msg_ |
int | nvir |
int | nocc |
int | nsocc |
Ref< ThreadGrp > | thr_ |
int | dynamic_ |
double | print_percent_ |
int | max_norb_ |
int * | symorb_irrep_ |
int * | symorb_num_ |
std::string | method_ |
std::string | algorithm_ |
int | do_d1_ |
int | do_d2_ |
int | nfuncmax |
double | hf_energy_ |
RefSCVector | hf_gradient_ |
double | restart_ecorr_ |
int | restart_orbital_v1_ |
int | restart_orbital_memgrp_ |
Protected Attributes inherited from sc::Wavefunction | |
ResultRefSCMatrix | natural_orbitals_ |
ResultRefDiagSCMatrix | natural_density_ |
Ref< GaussianBasisSet > | gbs_ |
int | debug_ |
Protected Attributes inherited from sc::MolecularEnergy | |
Ref< PointGroup > | initial_pg_ |
int | print_molecule_when_changed_ |
Protected Attributes inherited from sc::Function | |
Ref< SCMatrixKit > | matrixkit_ |
Used to construct new matrices. | |
RefSCVector | x_ |
The variables. | |
RefSCDimension | dim_ |
The dimension of x_. | |
AccResultdouble | value_ |
The value of the function at x_. | |
AccResultRefSCVector | gradient_ |
The gradient at x_. | |
AccResultRefSymmSCMatrix | hessian_ |
The hessian at x_. | |
bool | desired_value_accuracy_set_to_default_ |
bool | desired_gradient_accuracy_set_to_default_ |
bool | desired_hessian_accuracy_set_to_default_ |
bool | throw_if_tolerance_exceeded_ |
Additional Inherited Members | |
Static Public Member Functions inherited from sc::Wavefunction | |
static void | orbitals (const Ref< OrbitalSpace > &orbs, const std::vector< SCVector3 > &r, std::vector< double > &values) |
Static Public Member Functions inherited from sc::SavableState | |
static void | save_state (SavableState *s, StateOut &) |
static SavableState * | restore_state (StateIn &si) |
Restores objects saved with save_state. More... | |
static SavableState * | key_restore_state (StateIn &si, const char *keyword) |
Like restore_state, but keyword is used to override values while restoring. | |
static SavableState * | dir_restore_state (StateIn &si, const char *objectname, const char *keyword=0) |
The MBPT2 class implements several second-order perturbation theory methods.
The KeyVal constructor.
reference
This gives the reference wavefunction. It must be an object of type CLSCF for closed-shell molecules and HSOSSCF for open-shell molecules. The is no default.
nfzc
The number of frozen core orbitals. The default is 0. If no atoms have an atomic number greater than 30, then the number of orbitals to be frozen can be automatically determined by specifying nfzc = auto.
nfzv
The number of frozen virtual orbitals. The default is 0.
memory
The amount of memory, in bytes, that each processor may use.
method
This gives a string that must take on one of the values below. The default is mp for closed-shell systems and zapt for open-shell systems.
mp
Use Møller-Plesset perturbation theory. This is only valid for closed-shell systems. Energies and gradients can be computed with this method.
opt1
Use the OPT1 variant of open-shell perturbation theory. Only energies can be computed for open-shell systems.
opt2
Use the OPT2 variant of open-shell perturbation theory. Only energies can be computed for open-shell systems.
zapt
Use the ZAPT variant of open-shell perturbation theory. Only energies can be computed for open-shell systems.
algorithm
This gives a string that must take on one of the values given below. The default is memgrp for closed-shell systems. For open-shell systems v1 is used for a small number of processors and v2 is used otherwise.
memgrp
Use the distributed shared memory algorithm (which uses a MemoryGrp object). This is only valid for MP2 energies and gradients.
v1
Use algorithm V1. Only energies can be computed. The maximum number of processors that can be utilized is the number of virtual orbitals. This algorithm computes few integrals than the others, but has higher communication requirements.
v2
Use algorithm V2. Only energies can be computed. The maximum number of processors that can be utilized is the number of shells.
v2lb
Use a modified V2 algorithm that may compute more two electron integrals, but may get better load balance on the part of the calculation. Only energies can be computed. This is recommended only for computations involving large molecules (where the transformation is dominant) on very many processors (approaching the number of shells).
The v1 and v2 algorithms are discussed in Ida M. B. Nielsen and Edward T. Seidl, J. Comp. Chem. 16, 1301 (1995). The memgrp algorithm is discussed in Ida M. B. Nielsen, Chem. Phys. Lett. 255, 210 (1996).
memorygrp
A MemoryGrp object is used by the memgrp algorithm. If this is not given the program will try to find an appropriate default.
dynamic
This boolean keyword specifies whether dynamic load balancing is used. The default is false.
|
virtual |
must overload this in a derived class if analytic gradient can be computed
Reimplemented from sc::MolecularEnergy.
Reimplemented in sc::MBPT2_R12.
|
protectedvirtual |
Recompute at least the results that have compute true and are not already computed.
This should only be called by Result's members.
Implements sc::Compute.
Reimplemented in sc::MBPT2_R12.
|
virtual |
Computes the S (or J) magnetic moment of the target state(s), in units of .
Can be evaluated from density and overlap, as;
but derived Wavefunction may have this value as user input.
Reimplemented from sc::Wavefunction.
|
virtual |
Marks all results as being out of date.
Any subsequent access to results will cause Compute::compute() to be called.
Reimplemented from sc::Compute.
Reimplemented in sc::MBPT2_R12.
|
virtual |
Save the base classes (with save_data_state) and the members in the same order that the StateIn CTOR initializes them.
This must be implemented by the derived class if the class has data.
Reimplemented from sc::MolecularEnergy.
Reimplemented in sc::MBPT2_R12.
|
virtual |
Reimplemented from sc::Function.
Reimplemented in sc::MBPT2_R12.