MPQC
3.0.0-alpha
|
The SCF class is the base for all classes that use a self-consistent field procedure to solve an effective one body problem. More...
#include <chemistry/qc/scf/scf.h>
Public Member Functions | |
SCF (StateIn &) | |
SCF (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... | |
RefSCMatrix | oso_eigenvectors () |
Returns the orthogonal MO (columns) to orthogonal-SO (rows) transformation matrix. | |
RefDiagSCMatrix | eigenvalues () |
Returns the MO basis eigenvalues. | |
int | spin_unrestricted () |
Return 1 if the alpha orbitals are not equal to the beta orbitals. | |
virtual int | n_fock_matrices () const =0 |
virtual RefSymmSCMatrix | fock (int)=0 |
virtual RefSymmSCMatrix | effective_fock ()=0 |
virtual Ref< DensityFittingInfo > | dfinfo () const |
return the DensityFittingInfo object used to implement compute() this is important to be able to reconstruct the Fock matrix | |
virtual double | one_body_energy () |
virtual void | two_body_energy (double &ec, double &ex) |
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 | purge () |
This function purges any caches of data in MolecularEnergy. More... | |
void | print (std::ostream &o=ExEnv::out0()) const |
Print information about the object. | |
Public Member Functions inherited from sc::OneBodyWavefunction | |
OneBodyWavefunction (StateIn &) | |
OneBodyWavefunction (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... | |
int | nelectron () |
Returns the number of electrons. | |
void | set_desired_value_accuracy (double eps) |
Overload of Function::set_desired_value_accuracy(). More... | |
RefSCMatrix | mo_to_orthog_so () |
/ Returns the SO to MO transformation matrix. More... | |
RefSCMatrix | eigenvectors () |
Deprecated. More... | |
virtual double | occupation (int irrep, int vectornum)=0 |
Returns the occupation. More... | |
double | occupation (int vectornum) |
Returns the occupation. More... | |
virtual double | alpha_occupation (int irrep, int vectornum) |
Returns the alpha occupation. More... | |
virtual double | beta_occupation (int irrep, int vectornum) |
Returns the beta occupation. More... | |
double | alpha_occupation (int vectornum) |
Returns the alpha occupation. More... | |
double | beta_occupation (int vectornum) |
Returns the beta occupation. More... | |
virtual RefSCMatrix | oso_alpha_eigenvectors () |
virtual RefSCMatrix | oso_beta_eigenvectors () |
virtual RefSCMatrix | alpha_eigenvectors () |
virtual RefSCMatrix | beta_eigenvectors () |
virtual RefDiagSCMatrix | alpha_eigenvalues () |
virtual RefDiagSCMatrix | beta_eigenvalues () |
virtual RefDiagSCMatrix | projected_eigenvalues (const Ref< OneBodyWavefunction > &guess_wfn, int alp=1) |
Imports the eigenvalues of guess_wfn . More... | |
virtual RefSCMatrix | projected_eigenvectors (const Ref< OneBodyWavefunction > &guess_wfn, int alp=1) |
Projects the density (not eigenvalues) of guess_wfn into the current basis set. More... | |
virtual RefSCMatrix | hcore_guess () |
Return a guess vector. More... | |
virtual RefSCMatrix | hcore_guess (RefDiagSCMatrix &val) |
Return a guess vector and the eigenvalues. More... | |
void | symmetry_changed () |
Call this if you have changed the molecular symmetry of the molecule contained by this MolecularEnergy. | |
double | orbital (const SCVector3 &r, int iorb) |
returns the value of MO iorb at point r. To compute several MOs at several points use orbitals() instead | |
void | orbitals (const std::vector< SCVector3 > &r, std::vector< double > &values, unsigned int first, unsigned int last, bool energy_ordered=false) |
computes values of MOs in range [first,last] at points r and store them in an array values More... | |
double | orbital_density (const SCVector3 &r, int iorb, double *orbval=0) |
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 double | magnetic_moment () const |
Computes the S (or J) magnetic moment of the target state(s), in units of . More... | |
virtual RefSymmSCMatrix | density ()=0 |
Returns the SO density. | |
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... | |
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 |
virtual int | value_implemented () 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 Types | |
enum | Access { Read, Write, Accum } |
Protected Member Functions | |
virtual void | init_threads () |
virtual void | done_threads () |
virtual void | compute () |
Recompute at least the results that have compute true and are not already computed. More... | |
virtual double | compute_vector (double &, double enuclear) |
virtual Ref< SCExtrapError > | extrap_error () |
virtual void | compute_gradient (const RefSCVector &) |
virtual void | compute_hessian (const RefSymmSCMatrix &) |
virtual void | savestate_iter (int) |
virtual void | savestate_to_file (const std::string &filename) |
signed char * | init_pmax (double *) |
RefSymmSCMatrix | get_local_data (const RefSymmSCMatrix &, double *&, Access) |
virtual void | initial_vector () |
virtual void | obsolete_vector () |
Obsolete scf vector so that next call to initial_vector() will cause recomputation. | |
void | init_mem (int) |
void | so_density (const RefSymmSCMatrix &d, double occ, int alp=1) |
int * | read_occ (const Ref< KeyVal > &, const char *name, int nirrep) |
virtual void | set_occupations (const RefDiagSCMatrix &)=0 |
virtual void | init_vector ()=0 |
virtual void | done_vector ()=0 |
virtual double | new_density ()=0 |
virtual void | reset_density ()=0 |
virtual double | scf_energy ()=0 |
virtual Ref< SCExtrapData > | initial_extrap_data () |
virtual Ref< SCExtrapData > | extrap_data ()=0 |
virtual void | ao_fock (double accuracy)=0 |
virtual void | init_gradient ()=0 |
virtual void | done_gradient ()=0 |
virtual RefSymmSCMatrix | lagrangian ()=0 |
virtual RefSymmSCMatrix | gradient_density ()=0 |
virtual void | two_body_deriv (double *)=0 |
virtual void | init_hessian ()=0 |
virtual void | done_hessian ()=0 |
Protected Member Functions inherited from sc::OneBodyWavefunction | |
void | init_sym_info () |
int | form_occupations (int *&newocc, const int *oldocc) |
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_gradient_implemented () const |
must overload this in a derived class if analytic gradient can be computed More... | |
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 | guess_acc_ratio () |
how much lower is the desired accuracy of the guess? | |
static void | iter_print (int iter, double energy, double delta, double walltime, std::ostream &os=ExEnv::out0()) |
prints iteration log | |
Protected Attributes | |
int | compute_guess_ |
RefDiagSCMatrix | current_evals_ |
int | keep_guess_wfn_ |
Ref< OneBodyWavefunction > | guess_wfn_ |
int | always_use_guess_wfn_ |
Ref< SelfConsistentExtrapolation > | extrap_ |
Ref< AccumH > | accumdih_ |
Ref< AccumH > | accumddh_ |
int | maxiter_ |
int | miniter_ |
int | dens_reset_freq_ |
int | reset_occ_ |
size_t | storage_ |
int | print_all_evals_ |
int | print_occ_evals_ |
bool | fake_scf_convergence_after_fock_build_ |
int | fake_scf_convergence_after_n_iter_ |
double | level_shift_ |
Ref< MessageGrp > | scf_grp_ |
Ref< ThreadGrp > | threadgrp_ |
int | local_ |
int | local_dens_ |
Ref< TwoBodyInt > * | tbis_ |
std::string | previous_savestate_file_ |
RefSCMatrix | oso_scf_vector_ |
RefSCMatrix | oso_scf_vector_beta_ |
RefSymmSCMatrix | hcore_ |
Protected Attributes inherited from sc::OneBodyWavefunction | |
ResultRefSymmSCMatrix | density_ |
AccResultRefSCMatrix | oso_eigenvectors_ |
AccResultRefDiagSCMatrix | eigenvalues_ |
int | nirrep_ |
int * | nvecperirrep_ |
double * | occupations_ |
double * | alpha_occupations_ |
double * | beta_occupations_ |
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 SCF class is the base for all classes that use a self-consistent field procedure to solve an effective one body problem.
The KeyVal constructor.
maxiter
This integer specifies the maximum number of SCF iterations. The default is 40.
density_reset_frequency
This integer specifies how often, in term of SCF iterations, will be reset to . The default is 10.
reset_occupations
Reassign the occupations after each iteration based on the eigenvalues. This only has an effect for molecules with higher than symmetry. The default is false.
level_shift
Specificies a shift for the diagonal Fock matrix elements. Doubly occupied orbitals are shifted by this amount and singly occupied orbitals are shifted by half this amount. This can improve convergence in difficult cases. The units are hartrees and the default is 0.
extrap
This specifies an object of type SelfConsistentExtrapolation. The default is a DIIS object.
memory
The amount of memory that each processor may use. The default is 0 (minimal memory use).
local_density
If this is true, a local copy of the density and matrix will be made on all nodes, even if a distributed matrix specialization is used. The default is true.
guess_wavefunction
This specifies the initial guess for the solution to the SCF equations. This can be either a OneBodyWavefunction object or the name of file that contains the saved state of a OneBodyWavefunction object. By default, SuperpositionOfAtomicDensities object will be used; if that fails, core hamiltonian guess will be used instead.
keep_guess_wavefunction
The guess wavefunction is normally discarded after it is projected. Setting this boolean variable to true will cause the guess to be kept. This is useful when doing frequencies of symmetric molecules by finite displacements, because the wavefunction is lost whenever the molecule is displaced into lower symmetry.
always_use_guess_wavefunction
If the orbitals must be recomputed after they have already been computed once, then the old orbitals are used as the initial guess by default. However, if this option is true, then the guess wavefunction will be used, if available. If a guess wavefunction is not available, then a core Hamiltonian guess will be used. If this option is set to true, then keep_guess_wavefunction should also be set to true.
print_evals
Takes a boolean value. If true, print all eigenvalues after the SCF procedure converges. Takes a boolean value. The default is false.
print_occ_evals
Takes a boolean value. If true, print the occupied eigenvalues after the SCF procedure converges. The default is false.
accumdih
Optional. Takes an AccumH derivative. This provides additional contributions to the energy and the Fock matrix that are summed in once for the entire SCF procedure. AccumH's that are independent of the density should be given here.
accumddh
Optional. Takes an AccumH derivative. This provides additional contributions to the energy and the Fock matrix that are summed in each iteration SCF procedure. AccumH's that depend on the density must be given here.
iter_log
Optional. An SCFIterationLogger object to log data on a per-iteration basis. See the SCFIterationLogger class for more information.
|
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.
|
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.
|
virtual |
This function purges any caches of data in MolecularEnergy.
It is useful with MolecularEnergy objects that keep state when obsolete() is called (for example, it makes sense for SCF to keep its old eigenvector and reuse it as a guess when geometry changes). The default implementation does nothing and must be overloaded in classes which need it
Reimplemented from sc::MolecularEnergy.
|
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::UnrestrictedSCF, sc::TCSCF, sc::UHF, and sc::TCHF.