ergo
|
#include "molecule.h"
#include "basisinfo.h"
#include "integrals_2el.h"
#include "matrix_typedefs.h"
#include "densityfitting.h"
#include "grid_stream.h"
#include "SCF_statistics.h"
Go to the source code of this file.
Functions | |
void | output_sparsity (int n, const normalMatrix &M, const char *matrixName) |
void | output_sparsity_symm (int n, const symmMatrix &M, const char *matrixName) |
void | output_sparsity_triang (int n, const triangMatrix &M, const char *matrixName) |
int | compute_h_core_matrix_sparse (const IntegralInfo &integralInfo, const Molecule &molecule, const Molecule &extraCharges, ergo_real electric_field_x, ergo_real electric_field_y, ergo_real electric_field_z, const BasisInfoStruct &basisInfo, symmMatrix &H_core_Matrix_sparse, ergo_real threshold_integrals_1el, int noOfThreadsForV, mat::SizesAndBlocks const &matrix_size_block_info, std::vector< int > const &permutationHML, int const create_dipole_mtx=0, std::vector< int > const *const inversePermutationHML=0, std::string const *const calculation_identifier=0, std::string const *const method_and_basis_set=0) |
int | compute_h_core_matrix_simple_dense (const IntegralInfo &integralInfo, const Molecule &molecule, const BasisInfoStruct &basisInfo, symmMatrix &H_core_Matrix_sparse, ergo_real threshold_integrals_1el, int noOfThreadsForV, mat::SizesAndBlocks const &matrix_size_block_info, std::vector< int > const &permutationHML) |
int | get_gradient_for_given_mol_and_dens (const IntegralInfo &integralInfo, const Molecule &molecule, const BasisInfoStruct &basisInfo, const symmMatrix &D, ergo_real threshold_integrals_1el, mat::SizesAndBlocks const &matrix_size_block_info, std::vector< int > const &permutationHML, ergo_real *result_gradient_list) |
int | save_symmetric_matrix (symmMatrix &A, const BasisInfoStruct &basisInfo, const char *name, std::vector< int > const &inversePermutationHML) |
Saves specified symmetic matrix to a file of specified name. More... | |
int | add_disturbance_to_matrix (int n, symmMatrix &A, ergo_real disturbance, int specificElementCount, const int *elementIndexVector, std::vector< int > const &permutationHML) |
int | get_simple_starting_guess_sparse (int n, int noOfElectrons, symmMatrix &densityMatrix) |
int | write_diag_elements_to_file (int n, const symmMatrix &M, const char *fileName, std::vector< int > const &permutationHML) |
int | get_diag_matrix_from_file (int n, symmMatrix &M, const char *fileName, std::vector< int > const &permutationHML) |
int | write_full_matrix (int n, const symmMatrix &M, const char *fileName, std::vector< int > const &inversePermutationHML) |
int | write_basis_func_coord_file (const BasisInfoStruct &basisInfo) |
int | write_2el_integral_m_file (const BasisInfoStruct &basisInfo, const IntegralInfo &integralInfo) |
int | get_2e_matrix_and_energy_sparse (const BasisInfoStruct &basisInfo, const BasisInfoStruct &basisInfoDensFit, const Molecule &molecule, const IntegralInfo &integralInfo, symmMatrix &twoelMatrix_sparse, symmMatrix &densityMatrix_sparse, const JK::Params &J_K_params, const JK::ExchWeights &CAM_params, const Dft::GridParams &gridParams, int do_xc, ergo_real *energy_2el, int noOfElectrons, DensfitData *df_data, mat::SizesAndBlocks const &matrix_size_block_info, std::vector< int > const &permutationHML, std::vector< int > const &inversePermutationHML, int get_J_K_Fxc_matrices, symmMatrix &J_matrix, symmMatrix &K_matrix, symmMatrix &Fxc_matrix, SCF_statistics &stats) |
General routine for computing the two-electron part of the Fock/KS matrix. More... | |
int | get_2e_matrices_and_energy_sparse_unrestricted (const BasisInfoStruct &basisInfo, const BasisInfoStruct &basisInfoDensFit, const Molecule &molecule, const IntegralInfo &integralInfo, const JK::ExchWeights &CAM_params, symmMatrix &twoelMatrix_sparse_alpha, symmMatrix &twoelMatrix_sparse_beta, symmMatrix &densityMatrix_sparse_alpha, symmMatrix &densityMatrix_sparse_beta, const JK::Params &J_K_params, const Dft::GridParams &gridParams, int do_xc, ergo_real *energy_2el, int noOfElectrons, DensfitData *df_data, mat::SizesAndBlocks const &matrix_size_block_info, std::vector< int > const &permutationHML, std::vector< int > const &inversePermutationHML) |
int | get_2e_matrices_and_energy_restricted_open (const BasisInfoStruct &basisInfo, const BasisInfoStruct &basisInfoDensFit, const Molecule &molecule, const IntegralInfo &integralInfo, const JK::ExchWeights &CAM_params, symmMatrix &twoelMatrix_Fc, symmMatrix &twoelMatrix_Fo, symmMatrix &densityMatrix_sparse_alpha, symmMatrix &densityMatrix_sparse_beta, const JK::Params &J_K_params, const Dft::GridParams &gridParams, int do_xc, ergo_real *energy_2el, int noOfElectrons, DensfitData *df_data, mat::SizesAndBlocks const &matrix_size_block_info, std::vector< int > const &permutationHML, std::vector< int > const &inversePermutationHML) |
Computes G_c and G_o. More... | |
int | compute_FDSminusSDF_sparse (int n, symmMatrix &F_symm, symmMatrix &D_symm, symmMatrix &S_symm, normalMatrix &result, ergo_real sparse_threshold) |
int | determine_number_of_electrons_unrestricted (int noOfElectrons, int alpha_beta_diff, int *noOfElectrons_alpha, int *noOfElectrons_beta) |
void | get_hf_weight_and_cam_params (int use_dft, ergo_real *exch_param_alpha, ergo_real *exch_param_beta, ergo_real *exch_param_mu) |
void | get_dipole_moment (const symmMatrix &densityMatrix, const BasisInfoStruct &basisInfo, mat::SizesAndBlocks const &matrix_size_block_info, std::vector< int > const &permutationHML, const Molecule &molecule) |
void | do_mulliken_atomic_charges (const symmMatrix &densityMatrix, const symmMatrix &S_symm, const BasisInfoStruct &basisInfo, mat::SizesAndBlocks const &matrix_size_block_info, std::vector< int > const &permutationHML, std::vector< int > const &inversePermutationHML, const Molecule &molecule) |
void | do_mulliken_spin_densities (const symmMatrix &spinDensityMatrix, const symmMatrix &S_symm, const BasisInfoStruct &basisInfo, mat::SizesAndBlocks const &matrix_size_block_info, std::vector< int > const &permutationHML, std::vector< int > const &inversePermutationHML, const Molecule &molecule) |
void | do_density_images (const BasisInfoStruct &basisInfo, const Molecule &molecule, const ergo_real *densityMatrixFull_tot, const ergo_real *densityMatrixFull_spin, double output_density_images_boxwidth) |
void | do_acc_scan_J (const symmMatrix &D, const IntegralInfo &integralInfo, const BasisInfoStruct &basisInfo, triangMatrix &invCholFactor, bool doInvCholFactorTransformation, const JK::Params &J_K_params, mat::SizesAndBlocks const &matrix_size_block_info, std::vector< int > const &permutationHML, int nSteps, ergo_real startThresh, ergo_real stepFactor) |
void | do_acc_scan_K (symmMatrix &D, const IntegralInfo &integralInfo, const BasisInfoStruct &basisInfo, triangMatrix &invCholFactor, bool doInvCholFactorTransformation, const JK::ExchWeights &CAM_params, const JK::Params &J_K_params, mat::SizesAndBlocks const &matrix_size_block_info, std::vector< int > const &permutationHML, std::vector< int > const &inversePermutationHML, int nSteps, ergo_real startThresh, ergo_real stepFactor) |
void | do_acc_scan_Vxc (symmMatrix &D, const IntegralInfo &integralInfo, const BasisInfoStruct &basisInfo, const Molecule &molecule, const Dft::GridParams &gridParams, int noOfElectrons, triangMatrix &invCholFactor, bool doInvCholFactorTransformation, mat::SizesAndBlocks const &matrix_size_block_info, std::vector< int > const &permutationHML, std::vector< int > const &inversePermutationHML, int nSteps, ergo_real startThresh, ergo_real stepFactor) |
int add_disturbance_to_matrix | ( | int | n, |
symmMatrix & | A, | ||
ergo_real | disturbance, | ||
int | specificElementCount, | ||
const int * | elementIndexVector, | ||
std::vector< int > const & | permutationHML | ||
) |
int compute_FDSminusSDF_sparse | ( | int | n, |
symmMatrix & | F_symm, | ||
symmMatrix & | D_symm, | ||
symmMatrix & | S_symm, | ||
normalMatrix & | result, | ||
ergo_real | sparse_threshold | ||
) |
References mat::MatrixBase< Treal, Tmatrix >::clear(), do_output(), mat::MatrixGeneral< Treal, Tmatrix >::eucl_thresh(), mat::SingletonForTimings::instance(), LOG_AREA_DENSFROMF, LOG_AREA_SCF, LOG_CAT_INFO, output_current_memory_usage(), output_sparsity(), Util::TimeMeter::print(), mat::FileWritable::readFromFile(), mat::SingletonForTimings::reset(), mat::FileWritable::writeAndReadAll(), and mat::FileWritable::writeToFile().
Referenced by SCF_restricted::get_FDSminusSDF(), and SCF_unrestricted::get_FDSminusSDF().
int compute_h_core_matrix_simple_dense | ( | const IntegralInfo & | integralInfo, |
const Molecule & | molecule, | ||
const BasisInfoStruct & | basisInfo, | ||
symmMatrix & | H_core_Matrix_sparse, | ||
ergo_real | threshold_integrals_1el, | ||
int | noOfThreadsForV, | ||
mat::SizesAndBlocks const & | matrix_size_block_info, | ||
std::vector< int > const & | permutationHML | ||
) |
int compute_h_core_matrix_sparse | ( | const IntegralInfo & | integralInfo, |
const Molecule & | molecule, | ||
const Molecule & | extraCharges, | ||
ergo_real | electric_field_x, | ||
ergo_real | electric_field_y, | ||
ergo_real | electric_field_z, | ||
const BasisInfoStruct & | basisInfo, | ||
symmMatrix & | H_core_Matrix_sparse, | ||
ergo_real | threshold_integrals_1el, | ||
int | noOfThreadsForV, | ||
mat::SizesAndBlocks const & | matrix_size_block_info, | ||
std::vector< int > const & | permutationHML, | ||
int const | create_dipole_mtx = 0 , |
||
std::vector< int > const *const | inversePermutationHML = 0 , |
||
std::string const *const | calculation_identifier = 0 , |
||
std::string const *const | method_and_basis_set = 0 |
||
) |
References Atom::charge, mat::MatrixBase< Treal, Tmatrix >::clear(), compute_operator_matrix_sparse_symm(), compute_T_sparse(), compute_V_sparse(), Atom::coords, do_output(), Molecule::getAtom(), Molecule::getNoOfAtoms(), LOG_AREA_SCF, LOG_CAT_ERROR, LOG_CAT_INFO, max, BasisInfoStruct::noOfBasisFuncs, Util::TimeMeter::print(), mat::FileWritable::readFromFile(), mat::MatrixBase< Treal, Tmatrix >::resetSizesAndBlocks(), UNIT_one_Angstrom, write_matrix_in_matrix_market_format(), and mat::FileWritable::writeToFile().
Referenced by do_tdhf_dynamics(), and SCF_general::SCF_general().
int determine_number_of_electrons_unrestricted | ( | int | noOfElectrons, |
int | alpha_beta_diff, | ||
int * | noOfElectrons_alpha, | ||
int * | noOfElectrons_beta | ||
) |
References do_output(), LOG_AREA_SCF, LOG_CAT_ERROR, and LOG_CAT_INFO.
Referenced by SCF_unrestricted::SCF_unrestricted().
void do_acc_scan_J | ( | const symmMatrix & | D, |
const IntegralInfo & | integralInfo, | ||
const BasisInfoStruct & | basisInfo, | ||
triangMatrix & | invCholFactor, | ||
bool | doInvCholFactorTransformation, | ||
const JK::Params & | J_K_params, | ||
mat::SizesAndBlocks const & | matrix_size_block_info, | ||
std::vector< int > const & | permutationHML, | ||
int | nSteps, | ||
ergo_real | startThresh, | ||
ergo_real | stepFactor | ||
) |
References do_scan_and_report(), mat::FileWritable::readFromFile(), and mat::FileWritable::writeToFile().
Referenced by SCF_restricted::get_2e_part_and_energy(), and Ergo::registerInputVariables().
void do_acc_scan_K | ( | symmMatrix & | D, |
const IntegralInfo & | integralInfo, | ||
const BasisInfoStruct & | basisInfo, | ||
triangMatrix & | invCholFactor, | ||
bool | doInvCholFactorTransformation, | ||
const JK::ExchWeights & | CAM_params, | ||
const JK::Params & | J_K_params, | ||
mat::SizesAndBlocks const & | matrix_size_block_info, | ||
std::vector< int > const & | permutationHML, | ||
std::vector< int > const & | inversePermutationHML, | ||
int | nSteps, | ||
ergo_real | startThresh, | ||
ergo_real | stepFactor | ||
) |
References do_scan_and_report(), mat::FileWritable::readFromFile(), and mat::FileWritable::writeToFile().
Referenced by SCF_restricted::get_2e_part_and_energy(), and Ergo::registerInputVariables().
void do_acc_scan_Vxc | ( | symmMatrix & | D, |
const IntegralInfo & | integralInfo, | ||
const BasisInfoStruct & | basisInfo, | ||
const Molecule & | molecule, | ||
const Dft::GridParams & | gridParams, | ||
int | noOfElectrons, | ||
triangMatrix & | invCholFactor, | ||
bool | doInvCholFactorTransformation, | ||
mat::SizesAndBlocks const & | matrix_size_block_info, | ||
std::vector< int > const & | permutationHML, | ||
std::vector< int > const & | inversePermutationHML, | ||
int | nSteps, | ||
ergo_real | startThresh, | ||
ergo_real | stepFactor | ||
) |
References do_scan_and_report(), mat::FileWritable::readFromFile(), and mat::FileWritable::writeToFile().
Referenced by SCF_restricted::get_2e_part_and_energy(), and Ergo::registerInputVariables().
void do_density_images | ( | const BasisInfoStruct & | basisInfo, |
const Molecule & | molecule, | ||
const ergo_real * | densityMatrixFull_tot, | ||
const ergo_real * | densityMatrixFull_spin, | ||
double | output_density_images_boxwidth | ||
) |
References Atom::coords, do_output(), get_density(), get_no_of_primitives_for_density(), Molecule::getAtom(), Molecule::getNoOfAtoms(), integrate_density_in_box_2(), LOG_AREA_SCF, LOG_CAT_ERROR, LOG_CAT_INFO, Util::TimeMeter::print(), and write_gcube_file_header().
Referenced by SCF_restricted::output_density_images(), and SCF_unrestricted::output_density_images().
void do_mulliken_atomic_charges | ( | const symmMatrix & | densityMatrix, |
const symmMatrix & | S_symm, | ||
const BasisInfoStruct & | basisInfo, | ||
mat::SizesAndBlocks const & | matrix_size_block_info, | ||
std::vector< int > const & | permutationHML, | ||
std::vector< int > const & | inversePermutationHML, | ||
const Molecule & | molecule | ||
) |
References Atom::charge, do_output(), get_mulliken_charges(), Molecule::getAtom(), Molecule::getNoOfAtoms(), LOG_AREA_SCF, and LOG_CAT_INFO.
Referenced by SCF_restricted::do_mulliken_pop_stuff(), and SCF_unrestricted::do_mulliken_pop_stuff().
void do_mulliken_spin_densities | ( | const symmMatrix & | spinDensityMatrix, |
const symmMatrix & | S_symm, | ||
const BasisInfoStruct & | basisInfo, | ||
mat::SizesAndBlocks const & | matrix_size_block_info, | ||
std::vector< int > const & | permutationHML, | ||
std::vector< int > const & | inversePermutationHML, | ||
const Molecule & | molecule | ||
) |
References do_output(), get_mulliken_charges(), Molecule::getNoOfAtoms(), LOG_AREA_SCF, and LOG_CAT_INFO.
Referenced by SCF_unrestricted::do_mulliken_pop_stuff().
int get_2e_matrices_and_energy_restricted_open | ( | const BasisInfoStruct & | basisInfo, |
const BasisInfoStruct & | basisInfoDensFit, | ||
const Molecule & | molecule, | ||
const IntegralInfo & | integralInfo, | ||
const JK::ExchWeights & | CAM_params, | ||
symmMatrix & | twoelMatrix_Fc, | ||
symmMatrix & | twoelMatrix_Fo, | ||
symmMatrix & | densityMatrix_sparse_alpha, | ||
symmMatrix & | densityMatrix_sparse_beta, | ||
const JK::Params & | J_K_params, | ||
const Dft::GridParams & | gridParams, | ||
int | do_xc, | ||
ergo_real * | energy_2el, | ||
int | noOfElectrons, | ||
DensfitData * | df_data, | ||
mat::SizesAndBlocks const & | matrix_size_block_info, | ||
std::vector< int > const & | permutationHML, | ||
std::vector< int > const & | inversePermutationHML | ||
) |
Computes G_c and G_o.
G_c is defined as J_a + J_b + K_a + K_b. G_o is defined as J_a + J_b + K_a.
References JK::ExchWeights::alpha, mat::MatrixSymmetric< Treal, Tmatrix >::assignFromFull(), mat::MatrixBase< Treal, Tmatrix >::clear(), compute_J_by_boxes_sparse(), compute_K_by_boxes_sparse(), JK::ExchWeights::computeRangeSeparatedExchange, dft_get_uxc_mt(), do_output(), mat::MatrixSymmetric< Treal, Tmatrix >::fullMatrix(), LOG_AREA_SCF, LOG_CAT_ERROR, LOG_CAT_INFO, BasisInfoStruct::noOfBasisFuncs, output_current_memory_usage(), mat::MatrixBase< Treal, Tmatrix >::resetSizesAndBlocks(), mat::MatrixSymmetric< Treal, Tmatrix >::trace_ab(), JK::Params::use_densfit_for_J, and JK::Params::use_naive_fockmatrix_construction.
Referenced by SCF_unrestricted::get_2e_part_and_energy().
int get_2e_matrices_and_energy_sparse_unrestricted | ( | const BasisInfoStruct & | basisInfo, |
const BasisInfoStruct & | basisInfoDensFit, | ||
const Molecule & | molecule, | ||
const IntegralInfo & | integralInfo, | ||
const JK::ExchWeights & | CAM_params, | ||
symmMatrix & | twoelMatrix_sparse_alpha, | ||
symmMatrix & | twoelMatrix_sparse_beta, | ||
symmMatrix & | densityMatrix_sparse_alpha, | ||
symmMatrix & | densityMatrix_sparse_beta, | ||
const JK::Params & | J_K_params, | ||
const Dft::GridParams & | gridParams, | ||
int | do_xc, | ||
ergo_real * | energy_2el, | ||
int | noOfElectrons, | ||
DensfitData * | df_data, | ||
mat::SizesAndBlocks const & | matrix_size_block_info, | ||
std::vector< int > const & | permutationHML, | ||
std::vector< int > const & | inversePermutationHML | ||
) |
References JK::ExchWeights::alpha, mat::MatrixSymmetric< Treal, Tmatrix >::assignFromFull(), compute_JK_single_box(), JK::ExchWeights::computeRangeSeparatedExchange, dft_get_uxc_mt(), do_output(), mat::MatrixSymmetric< Treal, Tmatrix >::fullMatrix(), get_2e_matrices_and_energy_simple_HF_sparse_unrestricted(), get_2e_matrices_and_energy_simple_sparse_unrestricted(), get_trace(), LOG_AREA_SCF, LOG_CAT_ERROR, LOG_CAT_INFO, BasisInfoStruct::noOfBasisFuncs, output_current_memory_usage(), JK::Params::threshold_J, JK::Params::threshold_K, JK::Params::use_densfit_for_J, JK::Params::use_differential_density, JK::Params::use_fmm, and JK::Params::use_naive_fockmatrix_construction.
Referenced by SCF_unrestricted::get_2e_part_and_energy().
int get_2e_matrix_and_energy_sparse | ( | const BasisInfoStruct & | basisInfo, |
const BasisInfoStruct & | basisInfoDensFit, | ||
const Molecule & | molecule, | ||
const IntegralInfo & | integralInfo, | ||
symmMatrix & | twoelMatrix_sparse, | ||
symmMatrix & | densityMatrix_sparse, | ||
const JK::Params & | J_K_params, | ||
const JK::ExchWeights & | CAM_params, | ||
const Dft::GridParams & | gridParams, | ||
int | do_xc, | ||
ergo_real * | energy_2el, | ||
int | noOfElectrons, | ||
DensfitData * | df_data, | ||
mat::SizesAndBlocks const & | matrix_size_block_info, | ||
std::vector< int > const & | permutationHML, | ||
std::vector< int > const & | inversePermutationHML, | ||
int | get_J_K_Fxc_matrices, | ||
symmMatrix & | J_matrix, | ||
symmMatrix & | K_matrix, | ||
symmMatrix & | Fxc_matrix, | ||
SCF_statistics & | stats | ||
) |
General routine for computing the two-electron part of the Fock/KS matrix.
Both input and output matrices are in sparse form, but full matrix form may be used in intermediate steps. If FMM is not used, or if CAM is used, or if density fitting is used, full matrix format is applied.
basisInfo | the used basis set. |
basisInfoDensFit | the auxiliary basis set (may be NULL). |
molecule | position of atoms (used for eg. XC grid). |
integralInfo | - the integral evaluation object. |
twoelMatrix_sparse | - the evaluation result. |
densityMatrix_sparse | - the density for which 2e matrix is to be evaluated. |
J_K_params | the settings of the integral evaluation. |
do_xc | whether xc contribution to 2e matrix and energy are to be added. 1 means that the traditional full matrix code should be called, 2 means that the sparse variant is to be used. |
energy_2el | 2el energy contribution |
noOfElectrons | number of electrons... |
CAM_params | a structure containing parameters needed when functionals like CAMB3LYP are used. |
gridParams | a structure containing parameters for the grid. |
df_data | parameters related to density fitting. |
matrix_size_block_info | block sizes etc for hierarchic matrix library. |
get_J_K_Fxc_matrices | flag saying if matrices should be saved for statistics/testing purposes. If that feature is active, matrices are saved in parameters |
J_matrix | |
K_matrix | |
Fxc_matrix | . |
stats | a structure holding SCF statistics. |
permutationHML | - the permutation of basis functions, needed for transformation between the dense and sparse formats. |
inversePermutationHML | - the inverse permutation of basis functions, needed for transformation between the dense and sparse formats. |
References add_square_matrices(), JK::ExchWeights::alpha, mat::MatrixSymmetric< Treal, Tmatrix >::assignFromFull(), compute_2e_matrix_coulomb(), compute_2e_matrix_exchange(), compute_2e_matrix_list(), compute_2e_matrix_list_difden(), compute_2e_matrix_list_explicit(), JK::ExchWeights::computeRangeSeparatedExchange, dft_get_xc_mt(), do_output(), mat::MatrixSymmetric< Treal, Tmatrix >::fullMatrix(), get_2e_matrix_and_energy_simple_HF_sparse(), get_2e_matrix_and_energy_simple_sparse(), get_trace(), LOG_AREA_SCF, LOG_CAT_ERROR, LOG_CAT_INFO, BasisInfoStruct::noOfBasisFuncs, output_current_memory_usage(), JK::Params::threshold_J, JK::Params::threshold_K, JK::Params::use_densfit_for_J, JK::Params::use_differential_density, JK::Params::use_fmm, and JK::Params::use_naive_fockmatrix_construction.
Referenced by SCF_restricted::get_2e_part_and_energy().
int get_diag_matrix_from_file | ( | int | n, |
symmMatrix & | M, | ||
const char * | fileName, | ||
std::vector< int > const & | permutationHML | ||
) |
void get_dipole_moment | ( | const symmMatrix & | densityMatrix, |
const BasisInfoStruct & | basisInfo, | ||
mat::SizesAndBlocks const & | matrix_size_block_info, | ||
std::vector< int > const & | permutationHML, | ||
const Molecule & | molecule | ||
) |
References compute_dipole_moment_onecoord(), do_output(), LOG_AREA_SCF, and LOG_CAT_INFO.
Referenced by SCF_restricted::compute_dipole_moment(), and SCF_unrestricted::compute_dipole_moment().
int get_gradient_for_given_mol_and_dens | ( | const IntegralInfo & | integralInfo, |
const Molecule & | molecule, | ||
const BasisInfoStruct & | basisInfo, | ||
const symmMatrix & | D, | ||
ergo_real | threshold_integrals_1el, | ||
mat::SizesAndBlocks const & | matrix_size_block_info, | ||
std::vector< int > const & | permutationHML, | ||
ergo_real * | result_gradient_list | ||
) |
References compute_gradient_of_nucl_and_trDV().
Referenced by SCF_restricted::compute_gradient_fixeddens().
void get_hf_weight_and_cam_params | ( | int | use_dft, |
ergo_real * | exch_param_alpha, | ||
ergo_real * | exch_param_beta, | ||
ergo_real * | exch_param_mu | ||
) |
References fun_get_cam_param(), and fun_get_hf_weight.
Referenced by es_run(), SCF_general::SCF_general(), and ErgoE2Evaluator::transform().
int get_simple_starting_guess_sparse | ( | int | n, |
int | noOfElectrons, | ||
symmMatrix & | densityMatrix | ||
) |
void output_sparsity | ( | int | n, |
const normalMatrix & | M, | ||
const char * | matrixName | ||
) |
References output_sparsity_template().
Referenced by compute_FDSminusSDF_sparse(), SCF_restricted::get_FDSminusSDF(), and SCF_unrestricted::get_FDSminusSDF().
void output_sparsity_symm | ( | int | n, |
const symmMatrix & | M, | ||
const char * | matrixName | ||
) |
References output_sparsity_template().
Referenced by SCF_restricted::get_new_density_matrix(), SCF_restricted::get_starting_guess_density(), SCF_unrestricted::get_starting_guess_density(), SCF_restricted::output_sparsity_S_F_D(), SCF_unrestricted::output_sparsity_S_F_D(), and SCF_general::SCF_general().
void output_sparsity_triang | ( | int | n, |
const triangMatrix & | M, | ||
const char * | matrixName | ||
) |
References output_sparsity_template().
Referenced by SCF_general::SCF_general().
int save_symmetric_matrix | ( | symmMatrix & | A, |
const BasisInfoStruct & | basisInfo, | ||
const char * | fileName, | ||
std::vector< int > const & | inversePermutationHML | ||
) |
Saves specified symmetic matrix to a file of specified name.
A | the matrix to save. The matrix must be saved to a backing store already. |
basisInfo | the basis set description. |
fileName | The file that will contain the matrix. It will be overwritten without further questions. |
inversePermutationHML | permutation information needed when using hierarchic matrices. |
References ddf_writeShellListAndDensityMatricesToFile(), do_output(), mat::MatrixGeneral< Treal, Tmatrix >::fullMatrix(), LOG_AREA_SCF, LOG_CAT_INFO, BasisInfoStruct::noOfBasisFuncs, mat::FileWritable::readFromFile(), and mat::FileWritable::writeToFile().
Referenced by SCF_restricted::save_final_potential(), and SCF_general::SCF_general().
int write_2el_integral_m_file | ( | const BasisInfoStruct & | basisInfo, |
const IntegralInfo & | integralInfo | ||
) |
References do_2e_integral(), and BasisInfoStruct::noOfBasisFuncs.
Referenced by SCF_general::SCF_general().
int write_basis_func_coord_file | ( | const BasisInfoStruct & | basisInfo | ) |
References BasisInfoStruct::basisFuncList, BasisFuncStruct_::centerCoords, and BasisInfoStruct::noOfBasisFuncs.
Referenced by SCF_general::SCF_general().
int write_diag_elements_to_file | ( | int | n, |
const symmMatrix & | M, | ||
const char * | fileName, | ||
std::vector< int > const & | permutationHML | ||
) |
int write_full_matrix | ( | int | n, |
const symmMatrix & | M, | ||
const char * | fileName, | ||
std::vector< int > const & | inversePermutationHML | ||
) |
References mat::MatrixSymmetric< Treal, Tmatrix >::fullMatrix().
Referenced by SCF_restricted::save_full_matrices_for_matlab().