File: scf_utils.h

package info (click to toggle)
ergo 3.8.2-1.1
  • links: PTS, VCS
  • area: main
  • in suites: sid, trixie
  • size: 17,568 kB
  • sloc: cpp: 94,763; ansic: 17,785; sh: 10,701; makefile: 1,403; yacc: 127; lex: 116; awk: 23
file content (307 lines) | stat: -rw-r--r-- 10,494 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
/* Ergo, version 3.8.2, a program for linear scaling electronic structure
 * calculations.
 * Copyright (C) 2023 Elias Rudberg, Emanuel H. Rubensson, Pawel Salek,
 * and Anastasia Kruchinina.
 * 
 * This program is free software: you can redistribute it and/or modify
 * it under the terms of the GNU General Public License as published by
 * the Free Software Foundation, either version 3 of the License, or
 * (at your option) any later version.
 * 
 * This program is distributed in the hope that it will be useful,
 * but WITHOUT ANY WARRANTY; without even the implied warranty of
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 * GNU General Public License for more details.
 * 
 * You should have received a copy of the GNU General Public License
 * along with this program.  If not, see <http://www.gnu.org/licenses/>.
 * 
 * Primary academic reference:
 * Ergo: An open-source program for linear-scaling electronic structure
 * calculations,
 * Elias Rudberg, Emanuel H. Rubensson, Pawel Salek, and Anastasia
 * Kruchinina,
 * SoftwareX 7, 107 (2018),
 * <http://dx.doi.org/10.1016/j.softx.2018.03.005>
 * 
 * For further information about Ergo, see <http://www.ergoscf.org>.
 */

/** @file scf_utils.h

    @brief Various utilities used by self-consistent field (SCF)
    code. For example, interface routines converting between HML
    matrix format used in main SCF code and elementwise format used by
    integral code.

    @author: Elias Rudberg <em>responsible</em>
*/

#ifndef SCF_UTILS_HEADER
#define SCF_UTILS_HEADER

#include "molecule.h"
#include "basisinfo.h"
#include "integrals_2el.h"
#include "matrix_typedefs.h"
#include "grid_stream.h"
#include "SCF_statistics.h"


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,
			     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);

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);

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);

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 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);

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);

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);

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);

int 
determine_number_of_electrons_unrestricted(int noOfElectrons, 
					   int alpha_beta_diff, 
					   int* noOfElectrons_alpha, 
					   int* noOfElectrons_beta);

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 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);


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 = "");

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);

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);

#endif