ergo
scf_utils.h
Go to the documentation of this file.
1 /* Ergo, version 3.7, a program for linear scaling electronic structure
2  * calculations.
3  * Copyright (C) 2018 Elias Rudberg, Emanuel H. Rubensson, Pawel Salek,
4  * and Anastasia Kruchinina.
5  *
6  * This program is free software: you can redistribute it and/or modify
7  * it under the terms of the GNU General Public License as published by
8  * the Free Software Foundation, either version 3 of the License, or
9  * (at your option) any later version.
10  *
11  * This program is distributed in the hope that it will be useful,
12  * but WITHOUT ANY WARRANTY; without even the implied warranty of
13  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  * GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License
17  * along with this program. If not, see <http://www.gnu.org/licenses/>.
18  *
19  * Primary academic reference:
20  * Ergo: An open-source program for linear-scaling electronic structure
21  * calculations,
22  * Elias Rudberg, Emanuel H. Rubensson, Pawel Salek, and Anastasia
23  * Kruchinina,
24  * SoftwareX 7, 107 (2018),
25  * <http://dx.doi.org/10.1016/j.softx.2018.03.005>
26  *
27  * For further information about Ergo, see <http://www.ergoscf.org>.
28  */
29 
40 #ifndef SCF_UTILS_HEADER
41 #define SCF_UTILS_HEADER
42 
43 #include "molecule.h"
44 #include "basisinfo.h"
45 #include "integrals_2el.h"
46 #include "matrix_typedefs.h"
47 #include "grid_stream.h"
48 #include "SCF_statistics.h"
49 
50 
51 void output_sparsity(int n, const normalMatrix & M, const char* matrixName);
52 void output_sparsity_symm(int n, const symmMatrix & M, const char* matrixName);
53 void output_sparsity_triang(int n, const triangMatrix & M, const char* matrixName);
54 
55 int
56 compute_h_core_matrix_sparse(const IntegralInfo& integralInfo,
57  const Molecule& molecule,
58  const Molecule& extraCharges,
59  ergo_real electric_field_x,
60  ergo_real electric_field_y,
61  ergo_real electric_field_z,
62  const BasisInfoStruct& basisInfo,
63  symmMatrix & H_core_Matrix_sparse,
64  ergo_real threshold_integrals_1el,
65  int noOfThreadsForV,
66  ergo_real boxSizeForVT,
67  ergo_real & result_nuclearRepulsionEnergy,
68  mat::SizesAndBlocks const & matrix_size_block_info,
69  std::vector<int> const & permutationHML,
70  int const create_dipole_mtx = 0,
71  std::vector<int> const * const inversePermutationHML = 0,
72  std::string const * const calculation_identifier = 0,
73  std::string const * const method_and_basis_set = 0);
74 
75 int
77  const Molecule& molecule,
78  const BasisInfoStruct& basisInfo,
79  symmMatrix & H_core_Matrix_sparse,
80  ergo_real threshold_integrals_1el,
81  int noOfThreadsForV,
82  mat::SizesAndBlocks const & matrix_size_block_info,
83  std::vector<int> const & permutationHML,
84  ergo_real & result_nuclearRepulsionEnergy);
85 
86 int
88  const Molecule& molecule,
89  const BasisInfoStruct& basisInfo,
90  const symmMatrix & D,
91  ergo_real threshold_integrals_1el,
92  mat::SizesAndBlocks const & matrix_size_block_info,
93  std::vector<int> const & permutationHML,
94  ergo_real* result_gradient_list);
95 
97  const BasisInfoStruct & basisInfo,
98  const char *name,
99  std::vector<int> const & inversePermutationHML);
100 
101 int
103  symmMatrix & A,
104  ergo_real disturbance,
105  int specificElementCount,
106  const int* elementIndexVector,
107  std::vector<int> const & permutationHML);
108 
109 int
111  int noOfElectrons,
112  symmMatrix & densityMatrix);
113 
114 int
116  const symmMatrix & M,
117  const char* fileName,
118  std::vector<int> const & permutationHML);
119 
120 int
122  symmMatrix & M,
123  const char* fileName,
124  std::vector<int> const & permutationHML);
125 
126 int
127 write_full_matrix(int n,
128  const symmMatrix & M,
129  const char* fileName,
130  std::vector<int> const & inversePermutationHML);
131 
132 int
134 
135 int
136 write_2el_integral_m_file(const BasisInfoStruct & basisInfo, const IntegralInfo & integralInfo);
137 
138 int
140  const Molecule& molecule,
141  const IntegralInfo& integralInfo,
142  symmMatrix & twoelMatrix_sparse,
143  symmMatrix & densityMatrix_sparse,
144  const JK::Params& J_K_params,
145  const JK::ExchWeights & CAM_params,
146  const Dft::GridParams& gridParams,
147  int do_xc,
148  ergo_real* energy_2el,
149  int noOfElectrons,
150  mat::SizesAndBlocks const & matrix_size_block_info,
151  std::vector<int> const & permutationHML,
152  std::vector<int> const & inversePermutationHML,
153  int get_J_K_Fxc_matrices,
154  symmMatrix & J_matrix,
155  symmMatrix & K_matrix,
156  symmMatrix & Fxc_matrix,
157  SCF_statistics & stats);
158 
159 int
161  const Molecule& molecule,
162  const IntegralInfo& integralInfo,
163  const JK::ExchWeights & CAM_params,
164  symmMatrix & twoelMatrix_sparse_alpha,
165  symmMatrix & twoelMatrix_sparse_beta,
166  symmMatrix & densityMatrix_sparse_alpha,
167  symmMatrix & densityMatrix_sparse_beta,
168  const JK::Params& J_K_params,
169  const Dft::GridParams& gridParams,
170  int do_xc,
171  ergo_real* energy_2el,
172  int noOfElectrons,
173  mat::SizesAndBlocks const & matrix_size_block_info,
174  std::vector<int> const & permutationHML,
175  std::vector<int> const & inversePermutationHML);
176 
177 int
179  const Molecule& molecule,
180  const IntegralInfo& integralInfo,
181  const JK::ExchWeights & CAM_params,
182  symmMatrix & twoelMatrix_Fc,
183  symmMatrix & twoelMatrix_Fo,
184  symmMatrix & densityMatrix_sparse_alpha,
185  symmMatrix & densityMatrix_sparse_beta,
186  const JK::Params& J_K_params,
187  const Dft::GridParams& gridParams,
188  int do_xc,
189  ergo_real* energy_2el,
190  int noOfElectrons,
191  mat::SizesAndBlocks const & matrix_size_block_info,
192  std::vector<int> const & permutationHML,
193  std::vector<int> const & inversePermutationHML);
194 
195 int
197  symmMatrix & F_symm,
198  symmMatrix & D_symm,
199  symmMatrix & S_symm,
200  normalMatrix & result,
201  ergo_real sparse_threshold);
202 
203 int
205  int alpha_beta_diff,
206  int* noOfElectrons_alpha,
207  int* noOfElectrons_beta);
208 
209 void
210 get_hf_weight_and_cam_params(int use_dft,
211  ergo_real* exch_param_alpha,
212  ergo_real* exch_param_beta,
213  ergo_real* exch_param_mu);
214 
215 int
217  int alpha_beta_diff,
218  int* noOfElectrons_alpha,
219  int* noOfElectrons_beta);
220 
221 void
222 do_mulliken_atomic_charges(const symmMatrix & densityMatrix,
223  const symmMatrix & S_symm,
224  const BasisInfoStruct & basisInfo,
225  mat::SizesAndBlocks const & matrix_size_block_info,
226  std::vector<int> const & permutationHML,
227  std::vector<int> const & inversePermutationHML,
228  const Molecule& molecule);
229 
230 void
231 do_mulliken_spin_densities(const symmMatrix & spinDensityMatrix,
232  const symmMatrix & S_symm,
233  const BasisInfoStruct & basisInfo,
234  mat::SizesAndBlocks const & matrix_size_block_info,
235  std::vector<int> const & permutationHML,
236  std::vector<int> const & inversePermutationHML,
237  const Molecule& molecule);
238 
239 void get_exp_value_pos_operator(const BasisInfoStruct & basisInfo,
240  const Molecule& molecule,
241  const symmMatrix & densityMatrix,
242  mat::SizesAndBlocks const & matrix_size_block_info,
243  std::vector<int> const & permutationHML,
244  std::vector<ergo_real> &mean,
245  std::vector<ergo_real> &std);
246 
247 
248 void
249 do_density_images(const BasisInfoStruct & basisInfo,
250  const Molecule& molecule,
251  const ergo_real* densityMatrixFull_tot,
252  const ergo_real* densityMatrixFull_spin,
253  double output_density_images_boxwidth,
254  const std::string &filename_id = "");
255 
256 void
257 do_acc_scan_J(const symmMatrix & D,
258  const IntegralInfo & integralInfo,
259  const BasisInfoStruct & basisInfo,
260  triangMatrix & invCholFactor,
261  bool doInvCholFactorTransformation,
262  const JK::Params & J_K_params,
263  mat::SizesAndBlocks const & matrix_size_block_info,
264  std::vector<int> const & permutationHML,
265  int nSteps,
266  ergo_real startThresh,
267  ergo_real stepFactor);
268 
269 void
271  const IntegralInfo & integralInfo,
272  const BasisInfoStruct & basisInfo,
273  triangMatrix & invCholFactor,
274  bool doInvCholFactorTransformation,
275  const JK::ExchWeights & CAM_params,
276  const JK::Params & J_K_params,
277  mat::SizesAndBlocks const & matrix_size_block_info,
278  std::vector<int> const & permutationHML,
279  std::vector<int> const & inversePermutationHML,
280  int nSteps,
281  ergo_real startThresh,
282  ergo_real stepFactor);
283 
284 void
286  const IntegralInfo & integralInfo,
287  const BasisInfoStruct & basisInfo,
288  const Molecule & molecule,
289  const Dft::GridParams & gridParams,
290  int noOfElectrons,
291  triangMatrix & invCholFactor,
292  bool doInvCholFactorTransformation,
293  mat::SizesAndBlocks const & matrix_size_block_info,
294  std::vector<int> const & permutationHML,
295  std::vector<int> const & inversePermutationHML,
296  int nSteps,
297  ergo_real startThresh,
298  ergo_real stepFactor);
299 
300 void
302  const std::string & calculation_identifier,
303  const std::string & method_and_basis_set,
304  const std::vector<int> & inversePermutationHML,
305  const BasisInfoStruct & basisInfo);
306 
307 #endif
#define A
void output_sparsity_symm(int n, const symmMatrix &M, const char *matrixName)
Definition: scf_utils.cc:376
double ergo_real
Definition: realtype.h:69
int get_simple_starting_guess_sparse(int n, int noOfElectrons, symmMatrix &densityMatrix)
Definition: scf_utils.cc:818
A structure describing the grid settings.
Definition: grid_params.h:59
Streaming grid generator.
int get_2e_matrices_and_energy_restricted_open(const BasisInfoStruct &basisInfo, 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, mat::SizesAndBlocks const &matrix_size_block_info, std::vector< int > const &permutationHML, std::vector< int > const &inversePermutationHML)
Computes G_c and G_o.
Definition: scf_utils.cc:1929
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)
Definition: scf_utils.cc:2325
Definition: integrals_2el.h:44
Definition: SCF_statistics.h:57
MatrixSymmetric< real, matri > symmMatrix
Definition: test_LanczosSeveralLargestEig.cc:69
int get_diag_matrix_from_file(int n, symmMatrix &M, const char *fileName, std::vector< int > const &permutationHML)
Definition: scf_utils.cc:840
Code for setting up basis functions starting from shells.
Class representing a molecule as a set of atoms with assiciated coordinates and charges of the atomic...
void output_sparsity_triang(int n, const triangMatrix &M, const char *matrixName)
Definition: scf_utils.cc:381
Describes dimensions of matrix and its blocks on all levels.
Definition: SizesAndBlocks.h:45
Representation of a molecule as a set of nuclei and total charge.
Definition: molecule.h:87
void create_mtx_files_with_different_orderings(const symmMatrix &A, const std::string &calculation_identifier, const std::string &method_and_basis_set, const std::vector< int > &inversePermutationHML, const BasisInfoStruct &basisInfo)
Definition: scf_utils.cc:2639
void output_sparsity(int n, const normalMatrix &M, const char *matrixName)
Definition: scf_utils.cc:371
Definition: integral_info.h:147
MatrixGeneral< real, matri > normalMatrix
Definition: test_LanczosSeveralLargestEig.cc:71
int write_diag_elements_to_file(int n, const symmMatrix &M, const char *fileName, std::vector< int > const &permutationHML)
Definition: scf_utils.cc:917
MatrixTriangular< real, matri > triangMatrix
Definition: test_LanczosSeveralLargestEig.cc:70
Contains coefficients needed for quick integral evaluation.
Definition: integral_info.h:93
int get_2e_matrix_and_energy_sparse(const BasisInfoStruct &basisInfo, 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, 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.
Definition: scf_utils.cc:1531
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, const std::string &filename_id="")
Definition: scf_utils.cc:2410
int write_basis_func_coord_file(const BasisInfoStruct &basisInfo)
Definition: scf_utils.cc:976
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, ergo_real &result_nuclearRepulsionEnergy)
Definition: scf_utils.cc:409
Definition: basisinfo.h:112
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, ergo_real boxSizeForVT, ergo_real &result_nuclearRepulsionEnergy, 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)
Definition: scf_utils.cc:457
int get_2e_matrices_and_energy_sparse_unrestricted(const BasisInfoStruct &basisInfo, 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, mat::SizesAndBlocks const &matrix_size_block_info, std::vector< int > const &permutationHML, std::vector< int > const &inversePermutationHML)
Definition: scf_utils.cc:1721
Header file with typedefs for matrix and vector types.
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)
Definition: scf_utils.cc:320
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.
Definition: scf_utils.cc:735
int determine_number_of_electrons_unrestricted(int noOfElectrons, int alpha_beta_diff, int *noOfElectrons_alpha, int *noOfElectrons_beta)
Definition: scf_utils.cc:2201
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)
Definition: scf_utils.cc:2611
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)
Definition: scf_utils.cc:158
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)
Definition: scf_utils.cc:235
void get_exp_value_pos_operator(const BasisInfoStruct &basisInfo, const Molecule &molecule, const symmMatrix &densityMatrix, mat::SizesAndBlocks const &matrix_size_block_info, std::vector< int > const &permutationHML, std::vector< ergo_real > &mean, std::vector< ergo_real > &std)
Definition: scf_utils.cc:2364
int add_disturbance_to_matrix(int n, symmMatrix &A, ergo_real disturbance, int specificElementCount, const int *elementIndexVector, std::vector< int > const &permutationHML)
Definition: scf_utils.cc:761
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)
Definition: scf_utils.cc:2296
int compute_FDSminusSDF_sparse(int n, symmMatrix &F_symm, symmMatrix &D_symm, symmMatrix &S_symm, normalMatrix &result, ergo_real sparse_threshold)
Definition: scf_utils.cc:2084
Parameters related to integral evaluation.
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)
Definition: scf_utils.cc:701
int write_2el_integral_m_file(const BasisInfoStruct &basisInfo, const IntegralInfo &integralInfo)
Definition: scf_utils.cc:995
Class for keeping timings and other statistics related to self-consistent field (SCF) procedure...
int write_full_matrix(int n, const symmMatrix &M, const char *fileName, std::vector< int > const &inversePermutationHML)
Definition: scf_utils.cc:959