Line data Source code
1 : !--------------------------------------------------------------------------------------------------!
2 : ! CP2K: A general program to perform molecular dynamics simulations !
3 : ! Copyright 2000-2026 CP2K developers group <https://cp2k.org> !
4 : ! !
5 : ! SPDX-License-Identifier: GPL-2.0-or-later !
6 : !--------------------------------------------------------------------------------------------------!
7 :
8 : ! **************************************************************************************************
9 : !> \brief function that builds the hartree fock exchange section of the input
10 : !> \par History
11 : !> 09.2007 created
12 : !> \author Manuel Guidon
13 : ! **************************************************************************************************
14 : MODULE input_cp2k_hfx
15 : USE bibliography, ONLY: Guidon2008, &
16 : Guidon2009
17 : USE cp_output_handling, ONLY: add_last_numeric, &
18 : cp_print_key_section_create, &
19 : high_print_level, &
20 : medium_print_level
21 : USE input_constants, ONLY: &
22 : do_potential_coulomb, do_potential_gaussian, do_potential_id, do_potential_long, &
23 : do_potential_mix_cl, do_potential_mix_cl_trunc, do_potential_mix_lg, do_potential_short, &
24 : do_potential_truncated, ehrenfest, gaussian, hfx_ri_do_2c_cholesky, hfx_ri_do_2c_diag, &
25 : hfx_ri_do_2c_iter
26 : USE input_keyword_types, ONLY: keyword_create, &
27 : keyword_release, &
28 : keyword_type
29 : USE input_section_types, ONLY: section_add_keyword, &
30 : section_add_subsection, &
31 : section_create, &
32 : section_release, &
33 : section_type
34 : USE input_val_types, ONLY: real_t
35 : USE kinds, ONLY: dp
36 : USE string_utilities, ONLY: s2a
37 : #include "./base/base_uses.f90"
38 :
39 : IMPLICIT NONE
40 : PRIVATE
41 :
42 : LOGICAL, PRIVATE, PARAMETER :: debug_this_module = .FALSE.
43 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'input_cp2k_hfx'
44 : INTEGER, PARAMETER, PUBLIC :: ri_mo = 1, ri_pmat = 2
45 :
46 : PUBLIC :: create_hfx_section
47 :
48 : CONTAINS
49 :
50 : ! **************************************************************************************************
51 : !> \brief creates the input section for the hf part
52 : !> \param section the section to create
53 : !> \author Manuel Guidon
54 : ! **************************************************************************************************
55 235944 : SUBROUTINE create_hfx_section(section)
56 : TYPE(section_type), POINTER :: section
57 :
58 : TYPE(keyword_type), POINTER :: keyword
59 : TYPE(section_type), POINTER :: print_key, subsection
60 :
61 235944 : CPASSERT(.NOT. ASSOCIATED(section))
62 : CALL section_create(section, __LOCATION__, name="HF", &
63 : description="Controls Hartree-Fock exchange for hybrid DFT, Hartree-Fock, "// &
64 : "and related post-Hartree-Fock workflows.", &
65 : n_keywords=5, n_subsections=2, repeats=.TRUE., &
66 707832 : citations=[Guidon2008, Guidon2009])
67 :
68 235944 : NULLIFY (keyword, print_key, subsection)
69 :
70 : CALL keyword_create(keyword, __LOCATION__, name="FRACTION", &
71 : description="Fraction of Hartree-Fock exchange to add to the total energy. "// &
72 : "1.0 implies standard Hartree-Fock if used with XC_FUNCTIONAL NONE. "// &
73 : "NOTE: In a mixed potential calculation this should be set to 1.0, otherwise "// &
74 : "all parts are multiplied with this factor. ", &
75 235944 : usage="FRACTION 1.0", default_r_val=1.0_dp)
76 235944 : CALL section_add_keyword(section, keyword)
77 235944 : CALL keyword_release(keyword)
78 :
79 : CALL keyword_create(keyword, __LOCATION__, name="TREAT_LSD_IN_CORE", &
80 : description="Determines how spin densities are taken into account. "// &
81 : "If true, the beta spin density is included via a second in core call. "// &
82 : "If false, alpha and beta spins are done in one shot ", &
83 235944 : usage="TREAT_LSD_IN_CORE TRUE", default_l_val=.FALSE.)
84 235944 : CALL section_add_keyword(section, keyword)
85 235944 : CALL keyword_release(keyword)
86 :
87 : CALL keyword_create(keyword, __LOCATION__, name="PW_HFX", &
88 : description="Compute the Hartree-Fock energy also in the plane wave basis. "// &
89 : "The value is ignored, and intended for debugging only.", &
90 235944 : usage="PW_HFX FALSE", default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
91 235944 : CALL section_add_keyword(section, keyword)
92 235944 : CALL keyword_release(keyword)
93 :
94 : CALL keyword_create(keyword, __LOCATION__, name="PW_HFX_BLOCKSIZE", &
95 : description="Improve the performance of pw_hfx at the cost of some additional memory "// &
96 : "by storing the realspace representation of PW_HFX_BLOCKSIZE states.", &
97 235944 : usage="PW_HFX_BLOCKSIZE 20", default_i_val=20)
98 235944 : CALL section_add_keyword(section, keyword)
99 235944 : CALL keyword_release(keyword)
100 :
101 235944 : NULLIFY (print_key)
102 : CALL cp_print_key_section_create(print_key, __LOCATION__, "HF_INFO", &
103 : description="Controls the printing basic info about hf method", &
104 235944 : print_level=medium_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
105 235944 : CALL section_add_subsection(section, print_key)
106 235944 : CALL section_release(print_key)
107 :
108 235944 : CALL create_hf_pbc_section(subsection)
109 235944 : CALL section_add_subsection(section, subsection)
110 235944 : CALL section_release(subsection)
111 :
112 235944 : CALL create_hf_screening_section(subsection)
113 235944 : CALL section_add_subsection(section, subsection)
114 235944 : CALL section_release(subsection)
115 :
116 235944 : CALL create_hf_potential_section(subsection)
117 235944 : CALL section_add_subsection(section, subsection)
118 235944 : CALL section_release(subsection)
119 :
120 235944 : CALL create_hf_load_balance_section(subsection)
121 235944 : CALL section_add_subsection(section, subsection)
122 235944 : CALL section_release(subsection)
123 :
124 235944 : CALL create_hf_memory_section(subsection)
125 235944 : CALL section_add_subsection(section, subsection)
126 235944 : CALL section_release(subsection)
127 :
128 235944 : CALL create_hf_ri_section(subsection)
129 235944 : CALL section_add_subsection(section, subsection)
130 235944 : CALL section_release(subsection)
131 :
132 235944 : END SUBROUTINE create_hfx_section
133 :
134 : ! **************************************************************************************************
135 : !> \brief !****f* input_cp2k_dft/create_hf_load_balance_section [1.0] *
136 : !>
137 : !> creates the input section for the hf potential part
138 : !> \param section the section to create
139 : !> \author Manuel Guidon
140 : ! **************************************************************************************************
141 235944 : SUBROUTINE create_hf_load_balance_section(section)
142 : TYPE(section_type), POINTER :: section
143 :
144 : TYPE(keyword_type), POINTER :: keyword
145 : TYPE(section_type), POINTER :: print_key
146 :
147 235944 : CPASSERT(.NOT. ASSOCIATED(section))
148 : CALL section_create(section, __LOCATION__, name="LOAD_BALANCE", &
149 : description="Parameters influencing the load balancing of the HF", &
150 : n_keywords=1, n_subsections=0, repeats=.FALSE., &
151 471888 : citations=[guidon2008])
152 :
153 235944 : NULLIFY (keyword)
154 : CALL keyword_create( &
155 : keyword, __LOCATION__, &
156 : name="NBINS", &
157 : description="Number of bins per process used to group atom quartets.", &
158 : usage="NBINS 32", &
159 235944 : default_i_val=64)
160 235944 : CALL section_add_keyword(section, keyword)
161 235944 : CALL keyword_release(keyword)
162 :
163 : CALL keyword_create( &
164 : keyword, __LOCATION__, &
165 : name="BLOCK_SIZE", &
166 : description="Determines the blocking used for the atomic quartet loops. "// &
167 : "A proper choice can speedup the calculation. The default (-1) is automatic.", &
168 : usage="BLOCK_SIZE 4", &
169 235944 : default_i_val=-1)
170 235944 : CALL section_add_keyword(section, keyword)
171 235944 : CALL keyword_release(keyword)
172 :
173 235944 : NULLIFY (keyword)
174 : CALL keyword_create( &
175 : keyword, __LOCATION__, &
176 : name="RANDOMIZE", &
177 : description="This flag controls the randomization of the bin assignment to processes. "// &
178 : "For highly ordered input structures with a bad load balance, setting "// &
179 : "this flag to TRUE might improve.", &
180 : usage="RANDOMIZE TRUE", &
181 235944 : default_l_val=.FALSE.)
182 235944 : CALL section_add_keyword(section, keyword)
183 235944 : CALL keyword_release(keyword)
184 :
185 235944 : NULLIFY (print_key)
186 : CALL cp_print_key_section_create(print_key, __LOCATION__, "PRINT", &
187 : description="Controls the printing of info about load balance", &
188 235944 : print_level=medium_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
189 235944 : CALL section_add_subsection(section, print_key)
190 :
191 235944 : CALL keyword_release(keyword)
192 : CALL keyword_create(keyword, __LOCATION__, &
193 : name="LOAD_BALANCE_INFO", &
194 : description="Activates the printing of load balance information ", &
195 : default_l_val=.FALSE., &
196 235944 : lone_keyword_l_val=.TRUE.)
197 235944 : CALL section_add_keyword(print_key, keyword)
198 235944 : CALL keyword_release(keyword)
199 235944 : CALL section_release(print_key)
200 :
201 235944 : END SUBROUTINE create_hf_load_balance_section
202 :
203 : ! **************************************************************************************************
204 : !> \brief !****f* input_cp2k_dft/create_hf_potential_section [1.0] *
205 : !>
206 : !> creates the input section for the hf potential part
207 : !> \param section the section to create
208 : !> \author Manuel Guidon
209 : ! **************************************************************************************************
210 235944 : SUBROUTINE create_hf_potential_section(section)
211 : TYPE(section_type), POINTER :: section
212 :
213 : TYPE(keyword_type), POINTER :: keyword
214 :
215 235944 : CPASSERT(.NOT. ASSOCIATED(section))
216 : CALL section_create(section, __LOCATION__, name="INTERACTION_POTENTIAL", &
217 : description="Defines the Coulomb, range-separated, mixed, or truncated interaction "// &
218 : "operator used for Hartree-Fock exchange.", &
219 : n_keywords=1, n_subsections=0, repeats=.FALSE., &
220 707832 : citations=[guidon2008, guidon2009])
221 :
222 235944 : NULLIFY (keyword)
223 : CALL keyword_create( &
224 : keyword, __LOCATION__, &
225 : name="POTENTIAL_TYPE", &
226 : description="Selects the interaction potential used for Hartree-Fock exchange. "// &
227 : "Periodic hybrid calculations commonly use a short-range, truncated, or mixed potential.", &
228 : usage="POTENTIAL_TYPE SHORTRANGE", &
229 : enum_c_vals=s2a("COULOMB", "SHORTRANGE", "LONGRANGE", "MIX_CL", "GAUSSIAN", &
230 : "MIX_LG", "IDENTITY", "TRUNCATED", "MIX_CL_TRUNC"), &
231 : enum_i_vals=[do_potential_coulomb, do_potential_short, do_potential_long, &
232 : do_potential_mix_cl, do_potential_gaussian, do_potential_mix_lg, &
233 : do_potential_id, do_potential_truncated, do_potential_mix_cl_trunc], &
234 : enum_desc=s2a("Coulomb potential: $\frac{1}{r}$", &
235 : "Shortrange potential: $\frac{\mathrm{erfc}(\omega \cdot r)}{r}$", &
236 : "Longrange potential: $\frac{\mathrm{erf}(\omega \cdot r)}{r}$", &
237 : "Mix coulomb and longrange potential: $\frac{1}{r} + \frac{\mathrm{erf}(\omega \cdot r)}{r}$", &
238 : "Damped Gaussian potential: $\exp{(-\omega^2 \cdot r^2)}$", &
239 : "Mix Gaussian and longrange potential: "// &
240 : "$\frac{\mathrm{erf}(\omega \cdot r)}{r} + \exp{(-\omega^2 \cdot r^2)}$", &
241 : "Overlap", &
242 : "Truncated coulomb potential: if (r < R_c) 1/r else 0", &
243 : "Truncated Mix coulomb and longrange potential, assumes/requires that the erf has fully decayed at R_c"), &
244 235944 : default_i_val=do_potential_coulomb)
245 235944 : CALL section_add_keyword(section, keyword)
246 235944 : CALL keyword_release(keyword)
247 :
248 235944 : NULLIFY (keyword)
249 : CALL keyword_create( &
250 : keyword, __LOCATION__, &
251 : name="OMEGA", &
252 : description="Parameter $\omega$ for short/longrange interaction", &
253 : usage="OMEGA 0.5", &
254 235944 : default_r_val=0.0_dp)
255 235944 : CALL section_add_keyword(section, keyword)
256 235944 : CALL keyword_release(keyword)
257 :
258 : CALL keyword_create(keyword, __LOCATION__, name="SCALE_COULOMB", &
259 : description="Scales Hartree-Fock contribution arising from a coulomb potential. "// &
260 : "Only valid when doing a mixed potential calculation", &
261 235944 : usage="SCALE_COULOMB 1.0", default_r_val=1.0_dp)
262 235944 : CALL section_add_keyword(section, keyword)
263 235944 : CALL keyword_release(keyword)
264 :
265 : CALL keyword_create(keyword, __LOCATION__, name="SCALE_LONGRANGE", &
266 : description="Scales Hartree-Fock contribution arising from a longrange potential. "// &
267 : "Only valid when doing a mixed potential calculation", &
268 235944 : usage="SCALE_LONGRANGE 1.0", default_r_val=1.0_dp)
269 235944 : CALL section_add_keyword(section, keyword)
270 235944 : CALL keyword_release(keyword)
271 :
272 : CALL keyword_create(keyword, __LOCATION__, name="SCALE_GAUSSIAN", &
273 : description="Scales Hartree-Fock contribution arising from a gaussian potential. "// &
274 : "Only valid when doing a mixed potential calculation", &
275 235944 : usage="SCALE_GAUSSIAN 1.0", default_r_val=1.0_dp)
276 235944 : CALL section_add_keyword(section, keyword)
277 235944 : CALL keyword_release(keyword)
278 :
279 : CALL keyword_create(keyword, __LOCATION__, name="CUTOFF_RADIUS", &
280 : description="Cutoff radius for the truncated $\frac{1}{r}$ potential or the short-range "// &
281 : "$\frac{\mathrm{erfc}(\omega \cdot r)}{r}$ potential. For truncated Coulomb in a "// &
282 : "periodic cell, choose a radius compatible with the cell dimensions. The default value "// &
283 : "for short-range potentials when this keyword is omitted is solved from "// &
284 : "$\frac{\mathrm{erfc}(\omega \cdot r)}{r} = \epsilon_{\mathrm{schwarz}}$ "// &
285 : "by Newton-Raphson method, with $\epsilon_{\mathrm{schwarz}}$ set by SCREENING/EPS_SCHWARZ", &
286 : usage="CUTOFF_RADIUS 10.0", type_of_var=real_t, & ! default_r_val=10.0_dp,&
287 235944 : unit_str="angstrom")
288 235944 : CALL section_add_keyword(section, keyword)
289 235944 : CALL keyword_release(keyword)
290 :
291 : CALL keyword_create( &
292 : keyword, __LOCATION__, &
293 : name="T_C_G_DATA", &
294 : description="Location of the file t_c_g.dat that contains the data for the "// &
295 : "evaluation of the truncated gamma function ", &
296 : usage="T_C_G_DATA /data/t_c_g.dat", &
297 235944 : default_c_val="t_c_g.dat")
298 235944 : CALL section_add_keyword(section, keyword)
299 235944 : CALL keyword_release(keyword)
300 :
301 235944 : END SUBROUTINE create_hf_potential_section
302 :
303 : !****f* input_cp2k_dft/create_hf_screening_section [1.0] *
304 :
305 : ! **************************************************************************************************
306 : !> \brief creates the input section for the hf screening part
307 : !> \param section the section to create
308 : !> \author Manuel Guidon
309 : ! **************************************************************************************************
310 235944 : SUBROUTINE create_hf_screening_section(section)
311 : TYPE(section_type), POINTER :: section
312 :
313 : TYPE(keyword_type), POINTER :: keyword
314 :
315 235944 : CPASSERT(.NOT. ASSOCIATED(section))
316 : CALL section_create(section, __LOCATION__, name="SCREENING", &
317 : description="Controls screening thresholds for Hartree-Fock exchange integrals.", &
318 : n_keywords=1, n_subsections=0, repeats=.FALSE., &
319 707832 : citations=[guidon2008, guidon2009])
320 :
321 235944 : NULLIFY (keyword)
322 : CALL keyword_create( &
323 : keyword, __LOCATION__, &
324 : name="EPS_SCHWARZ", &
325 : description="Schwarz inequality threshold for screening near-field electronic repulsion integrals. "// &
326 : "Tighter values reduce screening error but increase cost.", &
327 : usage="EPS_SCHWARZ 1.0E-6", &
328 235944 : default_r_val=1.0E-10_dp)
329 235944 : CALL section_add_keyword(section, keyword)
330 235944 : CALL keyword_release(keyword)
331 :
332 235944 : NULLIFY (keyword)
333 : CALL keyword_create( &
334 : keyword, __LOCATION__, &
335 : name="EPS_SCHWARZ_FORCES", &
336 : description="Schwarz threshold used for force-related electronic repulsion integrals. "// &
337 : "This is approximately the force accuracy and should normally be similar to EPS_SCF. "// &
338 : "Default value is 100*EPS_SCHWARZ.", &
339 : usage="EPS_SCHWARZ_FORCES 1.0E-5", &
340 235944 : default_r_val=1.0E-6_dp)
341 235944 : CALL section_add_keyword(section, keyword)
342 235944 : CALL keyword_release(keyword)
343 :
344 235944 : NULLIFY (keyword)
345 : CALL keyword_create( &
346 : keyword, __LOCATION__, &
347 : name="SCREEN_P_FORCES", &
348 : description="Screens the electronic repulsion integrals for the forces "// &
349 : "using the density matrix. Will be disabled for the "// &
350 : "response part of forces in MP2/RPA/TDDFT. "// &
351 : "This results in a significant speedup for large systems, "// &
352 : "but might require a somewhat tigher EPS_SCHWARZ_FORCES.", &
353 : usage="SCREEN_P_FORCES TRUE", &
354 235944 : default_l_val=.TRUE.)
355 235944 : CALL section_add_keyword(section, keyword)
356 235944 : CALL keyword_release(keyword)
357 :
358 235944 : NULLIFY (keyword)
359 : CALL keyword_create(keyword, __LOCATION__, name="SCREEN_ON_INITIAL_P", &
360 : description="Screen on an initial density matrix. For the first MD step"// &
361 : " this matrix must be provided by a Restart File.", &
362 235944 : usage="SCREEN_ON_INITIAL_P TRUE", default_l_val=.FALSE.)
363 235944 : CALL section_add_keyword(section, keyword)
364 235944 : CALL keyword_release(keyword)
365 :
366 235944 : NULLIFY (keyword)
367 : CALL keyword_create(keyword, __LOCATION__, name="P_SCREEN_CORRECTION_FACTOR", &
368 : description="Recalculates integrals on the fly if the actual density matrix is"// &
369 : " larger by a given factor than the initial one. If the factor is set"// &
370 : " to 0.0_dp, this feature is disabled.", &
371 235944 : usage="P_SCREEN_CORRECTION_FACTOR 0.0_dp", default_r_val=0.0_dp)
372 235944 : CALL section_add_keyword(section, keyword)
373 235944 : CALL keyword_release(keyword)
374 :
375 235944 : END SUBROUTINE create_hf_screening_section
376 :
377 : ! **************************************************************************************************
378 : !> \brief creates the input section for the hf-pbc part
379 : !> \param section the section to create
380 : !> \author Manuel Guidon
381 : ! **************************************************************************************************
382 235944 : SUBROUTINE create_hf_pbc_section(section)
383 : TYPE(section_type), POINTER :: section
384 :
385 : TYPE(keyword_type), POINTER :: keyword
386 :
387 235944 : CPASSERT(.NOT. ASSOCIATED(section))
388 : CALL section_create(section, __LOCATION__, name="PERIODIC", &
389 : description="Sets up periodic boundary condition parameters if requested ", &
390 : n_keywords=1, n_subsections=0, repeats=.FALSE., &
391 707832 : citations=[guidon2008, guidon2009])
392 235944 : NULLIFY (keyword)
393 : CALL keyword_create( &
394 : keyword, __LOCATION__, &
395 : name="NUMBER_OF_SHELLS", &
396 : description="Number of shells taken into account for periodicity. "// &
397 : "By default, cp2k tries to automatically evaluate this number. "// &
398 : "This algorithm might be to conservative, resulting in some overhead. "// &
399 : "You can try to adjust this number in order to make a calculation cheaper. ", &
400 : usage="NUMBER_OF_SHELLS 2", &
401 235944 : default_i_val=-1)
402 235944 : CALL section_add_keyword(section, keyword)
403 235944 : CALL keyword_release(keyword)
404 :
405 235944 : END SUBROUTINE create_hf_pbc_section
406 :
407 : ! **************************************************************************************************
408 : !> \brief creates the input section for the hf-memory part
409 : !> \param section the section to create
410 : !> \author Manuel Guidon
411 : ! **************************************************************************************************
412 235944 : SUBROUTINE create_hf_memory_section(section)
413 : TYPE(section_type), POINTER :: section
414 :
415 : TYPE(keyword_type), POINTER :: keyword
416 :
417 235944 : CPASSERT(.NOT. ASSOCIATED(section))
418 : CALL section_create(section, __LOCATION__, name="MEMORY", &
419 : description="Sets up memory parameters for the storage of the ERI's if requested ", &
420 : n_keywords=1, n_subsections=0, repeats=.FALSE., &
421 471888 : citations=[guidon2008])
422 235944 : NULLIFY (keyword)
423 : CALL keyword_create( &
424 : keyword, __LOCATION__, &
425 : name="EPS_STORAGE_SCALING", &
426 : variants=["EPS_STORAGE"], &
427 : description="Scaling factor to scale eps_schwarz. Storage threshold for compression "// &
428 : "will be EPS_SCHWARZ*EPS_STORAGE_SCALING.", &
429 : usage="EPS_STORAGE 1.0E-2", &
430 471888 : default_r_val=1.0E0_dp)
431 235944 : CALL section_add_keyword(section, keyword)
432 235944 : CALL keyword_release(keyword)
433 :
434 : CALL keyword_create( &
435 : keyword, __LOCATION__, &
436 : name="MAX_MEMORY", &
437 : description="Defines the maximum amount of memory [MiB] to be consumed by the full HFX module. "// &
438 : "All temporary buffers and helper arrays are subtracted from this number. "// &
439 : "What remains will be used for storage of integrals. NOTE: This number "// &
440 : "is assumed to represent the memory available to one MPI process. "// &
441 : "When running a threaded version, cp2k automatically takes care of "// &
442 : "distributing the memory among all the threads within a process.", &
443 : usage="MAX_MEMORY 256", &
444 235944 : default_i_val=512)
445 235944 : CALL section_add_keyword(section, keyword)
446 235944 : CALL keyword_release(keyword)
447 :
448 : CALL keyword_create( &
449 : keyword, __LOCATION__, &
450 : name="STORAGE_LOCATION", &
451 : description="Loaction where ERI's are stored if MAX_DISK_SPACE /=0 "// &
452 : "Expects a path to a directory. ", &
453 : usage="STORAGE_LOCATION /data/scratch", &
454 235944 : default_c_val=".")
455 235944 : CALL section_add_keyword(section, keyword)
456 235944 : CALL keyword_release(keyword)
457 :
458 : CALL keyword_create( &
459 : keyword, __LOCATION__, &
460 : name="MAX_DISK_SPACE", &
461 : description="Defines the maximum amount of disk space [MiB] used to store precomputed "// &
462 : "compressed four-center integrals. If 0, nothing is stored to disk", &
463 : usage="MAX_DISK_SPACE 256", &
464 235944 : default_i_val=0)
465 235944 : CALL section_add_keyword(section, keyword)
466 235944 : CALL keyword_release(keyword)
467 :
468 : CALL keyword_create(keyword, __LOCATION__, name="TREAT_FORCES_IN_CORE", &
469 : description="Determines whether the derivative ERI's should be stored to RAM or not. "// &
470 : "Only meaningful when performing Ehrenfest MD. "// &
471 : "Memory usage is defined via MAX_MEMORY, i.e. the memory is shared wit the energy ERI's.", &
472 235944 : usage="TREAT_FORCES_IN_CORE TRUE", default_l_val=.FALSE.)
473 235944 : CALL section_add_keyword(section, keyword)
474 235944 : CALL keyword_release(keyword)
475 :
476 235944 : END SUBROUTINE create_hf_memory_section
477 :
478 : ! **************************************************************************************************
479 : !> \brief ...
480 : !> \param section ...
481 : ! **************************************************************************************************
482 235944 : SUBROUTINE create_hf_ri_section(section)
483 : TYPE(section_type), POINTER :: section
484 :
485 : TYPE(keyword_type), POINTER :: keyword
486 : TYPE(section_type), POINTER :: print_key, subsection
487 :
488 235944 : NULLIFY (keyword, print_key, subsection)
489 :
490 235944 : CPASSERT(.NOT. ASSOCIATED(section))
491 : CALL section_create(section, __LOCATION__, name="RI", &
492 : description="Parameters for RI methods in HFX, including RI-HFXk with "// &
493 : "k-point sampling. All keywords relevant to RI-HFXk have an "// &
494 235944 : "alias starting with KP_")
495 :
496 : CALL keyword_create(keyword, __LOCATION__, name="_SECTION_PARAMETERS_", &
497 : description="controls the activation of RI", &
498 : usage="&RI T", &
499 : default_l_val=.FALSE., &
500 235944 : lone_keyword_l_val=.TRUE.)
501 235944 : CALL section_add_keyword(section, keyword)
502 235944 : CALL keyword_release(keyword)
503 :
504 : CALL keyword_create(keyword, __LOCATION__, name="EPS_FILTER", &
505 : description="Filter threshold for DBT tensor contraction.", &
506 : variants=["KP_EPS_FILTER"], &
507 471888 : default_r_val=1.0E-09_dp)
508 235944 : CALL section_add_keyword(section, keyword)
509 235944 : CALL keyword_release(keyword)
510 :
511 : CALL keyword_create(keyword, __LOCATION__, name="EPS_FILTER_2C", &
512 : description="Filter threshold for 2c integrals. Default should be kept.", &
513 235944 : default_r_val=1.0E-12_dp)
514 235944 : CALL section_add_keyword(section, keyword)
515 235944 : CALL keyword_release(keyword)
516 :
517 : CALL keyword_create(keyword, __LOCATION__, &
518 : name="EPS_STORAGE_SCALING", &
519 : description="Scaling factor to scale EPS_FILTER for storage of 3-center integrals. Storage threshold "// &
520 : "will be EPS_FILTER*EPS_STORAGE_SCALING.", &
521 235944 : default_r_val=0.01_dp)
522 235944 : CALL section_add_keyword(section, keyword)
523 235944 : CALL keyword_release(keyword)
524 :
525 : CALL keyword_create(keyword, __LOCATION__, name="EPS_FILTER_MO", &
526 : description="Filter threshold for contraction of 3-center integrals with MOs. "// &
527 : "Default should be kept.", &
528 235944 : default_r_val=1.0E-12_dp)
529 235944 : CALL section_add_keyword(section, keyword)
530 235944 : CALL keyword_release(keyword)
531 :
532 : CALL keyword_create(keyword, __LOCATION__, name="OMEGA", &
533 : description="The range parameter for the short range operator (in 1/a0). "// &
534 : "Default is OMEGA from INTERACTION_POTENTIAL. ", &
535 : variants=["KP_OMEGA"], &
536 : default_r_val=0.0_dp, &
537 471888 : repeats=.FALSE.)
538 235944 : CALL section_add_keyword(section, keyword)
539 235944 : CALL keyword_release(keyword)
540 :
541 : CALL keyword_create(keyword, __LOCATION__, name="CUTOFF_RADIUS", &
542 : description="The cutoff radius (in Angstroms) for the truncated Coulomb operator. "// &
543 : "Default is CUTOFF_RADIUS from INTERACTION_POTENTIAL. ", &
544 : variants=["KP_CUTOFF_RADIUS"], &
545 : default_r_val=0.0_dp, &
546 : repeats=.FALSE., &
547 471888 : unit_str="angstrom")
548 235944 : CALL section_add_keyword(section, keyword)
549 235944 : CALL keyword_release(keyword)
550 :
551 : CALL keyword_create(keyword, __LOCATION__, name="SCALE_COULOMB", &
552 : description="Scales Hartree-Fock contribution arising from a coulomb potential. "// &
553 : "Only valid when doing a mixed potential calculation. "// &
554 : "Default is SCALE_COULOMB from INTERACTION_POTENTIAL", &
555 235944 : usage="SCALE_COULOMB 1.0", default_r_val=1.0_dp)
556 235944 : CALL section_add_keyword(section, keyword)
557 235944 : CALL keyword_release(keyword)
558 :
559 : CALL keyword_create(keyword, __LOCATION__, name="SCALE_LONGRANGE", &
560 : description="Scales Hartree-Fock contribution arising from a longrange potential. "// &
561 : "Only valid when doing a mixed potential calculation. "// &
562 : "Default if SCALE_LONGRANGE from INTERACTION_POTENTIAL", &
563 235944 : usage="SCALE_LONGRANGE 1.0", default_r_val=1.0_dp)
564 235944 : CALL section_add_keyword(section, keyword)
565 235944 : CALL keyword_release(keyword)
566 :
567 : CALL keyword_create(keyword, __LOCATION__, name="KP_NGROUPS", &
568 : description="The number of MPI subgroup that work in parallel during the SCF. "// &
569 : "The default value is 1. Using N subgroups should speed up the "// &
570 : "calculation by a factor ~N, at the cost of N times more memory usage.", &
571 : variants=["NGROUPS"], &
572 : usage="KP_NGROUPS {int}", &
573 : repeats=.FALSE., &
574 471888 : default_i_val=1)
575 235944 : CALL section_add_keyword(section, keyword)
576 235944 : CALL keyword_release(keyword)
577 :
578 : CALL keyword_create(keyword, __LOCATION__, name="KP_USE_DELTA_P", &
579 : description="This kweyword controls whether the KS matrix at each SCF cycle "// &
580 : "is built by adding the contribution of the denisty difference (wrt to previous step) "// &
581 : "to the KS matrix of the previous step. As the SCF converges, the density fluctuations "// &
582 : "get smaller and sparsity increases, leading to faster SCF steps. Not always "// &
583 : "numerically stable => turn off if SCF struggles to converge.", &
584 : variants=s2a("USE_DELTA_P", "KP_USE_P_DIFF", "USE_P_DIFF"), &
585 : usage="KP_USE_DELTA_P {logical}", &
586 : repeats=.FALSE., &
587 235944 : default_l_val=.TRUE.)
588 235944 : CALL section_add_keyword(section, keyword)
589 235944 : CALL keyword_release(keyword)
590 :
591 : CALL keyword_create(keyword, __LOCATION__, name="KP_STACK_SIZE", &
592 : description="When doing contraction over periodic cells of the type: "// &
593 : "T_mu^a,nu^b,P^c = (mu^a nu^b | Q^d) * (Q^d | P^c), with "// &
594 : "a,b,c,d labeling cells, there are in principle Ncells "// &
595 : "contractions taking place. Because a smaller number of "// &
596 : "contractions involving larger tensors is more efficient, "// &
597 : "the tensors can be stacked along the d direction. STCK_SIZE "// &
598 : "controls the size of this stack. Larger stacks are more efficient, "// &
599 : "but required more memory.", &
600 : variants=["STACK_SIZE"], &
601 : usage="KP_STACK_SIZE {int}", &
602 : repeats=.FALSE., &
603 471888 : default_i_val=16)
604 235944 : CALL section_add_keyword(section, keyword)
605 235944 : CALL keyword_release(keyword)
606 :
607 : CALL keyword_create(keyword, __LOCATION__, name="KP_RI_BUMP_FACTOR", &
608 : variants=s2a("RI_BUMP", "BUMP", "BUMP_FACTOR"), &
609 : description="In KP-RI-HFX, the extended RI basis set has a bump radius. "// &
610 : "All basis elements within that radius contribute with full weight. "// &
611 : "All basis elements beyond that radius have decaying weight, from "// &
612 : "1 at the bump radius, to zero at the RI extension radius. The "// &
613 : "bump radius is calculated as a fraction of the RI extension radius: "// &
614 : "bump radius = KP_RI_NUMP_FACTOR * RI extension radius", &
615 : default_r_val=0.85_dp, &
616 235944 : repeats=.FALSE.)
617 235944 : CALL section_add_keyword(section, keyword)
618 235944 : CALL keyword_release(keyword)
619 :
620 : CALL keyword_create(keyword, __LOCATION__, name="RI_METRIC", &
621 : description="The type of RI operator. "// &
622 : "Default is POTENTIAL_TYPE from INTERACTION_POTENTIAL. "// &
623 : "The standard "// &
624 : "Coulomb operator cannot be used in periodic systems.", &
625 : usage="RI_METRIC {string}", &
626 : repeats=.FALSE., &
627 : variants=["KP_RI_METRIC"], &
628 : default_i_val=0, &
629 : enum_c_vals=s2a("HFX", "COULOMB", "IDENTITY", "TRUNCATED", "SHORTRANGE"), &
630 : enum_desc=s2a("Same as HFX operator", &
631 : "Standard Coulomb operator: 1/r", &
632 : "Overlap", &
633 : "Truncated Coulomb operator: 1/r if (r<R_c), 0 otherwise ", &
634 : "Short range: erfc(omega*r)/r"), &
635 : enum_i_vals=[0, do_potential_coulomb, do_potential_id, do_potential_truncated, &
636 471888 : do_potential_short])
637 235944 : CALL section_add_keyword(section, keyword)
638 235944 : CALL keyword_release(keyword)
639 :
640 : CALL keyword_create(keyword, __LOCATION__, name="2C_MATRIX_FUNCTIONS", &
641 : description="Methods for matrix inverse and matrix square root.", &
642 : default_i_val=hfx_ri_do_2c_cholesky, &
643 : enum_c_vals=s2a("DIAG", "CHOLESKY", "ITER"), &
644 : enum_desc=s2a("Diagonalization with eigenvalue quenching: stable", &
645 : "Cholesky: not stable in case of ill-conditioned RI basis", &
646 : "Iterative algorithms: linear scaling "// &
647 : "Hotelling's method for inverse and Newton-Schulz iteration for matrix square root"), &
648 235944 : enum_i_vals=[hfx_ri_do_2c_diag, hfx_ri_do_2c_cholesky, hfx_ri_do_2c_iter])
649 235944 : CALL section_add_keyword(section, keyword)
650 235944 : CALL keyword_release(keyword)
651 :
652 : CALL keyword_create(keyword, __LOCATION__, name="EPS_EIGVAL", &
653 : description="Throw away linear combinations of RI basis functions with a small eigenvalue, "// &
654 : "this is applied only if 2C_MATRIX_FUNCTIONS DIAG", &
655 235944 : default_r_val=1.0e-7_dp)
656 235944 : CALL section_add_keyword(section, keyword)
657 235944 : CALL keyword_release(keyword)
658 :
659 : CALL keyword_create(keyword, __LOCATION__, name="CHECK_2C_MATRIX", &
660 : description="Report accuracy for the inverse/sqrt of the 2-center integral matrix.", &
661 235944 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
662 235944 : CALL section_add_keyword(section, keyword)
663 235944 : CALL keyword_release(keyword)
664 :
665 : CALL keyword_create( &
666 : keyword, __LOCATION__, &
667 : name="CALC_COND_NUM", &
668 : variants=["CALC_CONDITION_NUMBER"], &
669 : description="Calculate the condition number of integral matrices.", &
670 : usage="CALC_COND_NUM", &
671 : default_l_val=.FALSE., &
672 471888 : lone_keyword_l_val=.TRUE.)
673 235944 : CALL section_add_keyword(section, keyword)
674 235944 : CALL keyword_release(keyword)
675 :
676 : CALL keyword_create(keyword, __LOCATION__, name="SQRT_ORDER", &
677 : description="Order of the iteration method for the calculation of "// &
678 : "the sqrt of 2-center integral matrix.", &
679 235944 : default_i_val=3)
680 235944 : CALL section_add_keyword(section, keyword)
681 235944 : CALL keyword_release(keyword)
682 :
683 : CALL keyword_create(keyword, __LOCATION__, name="EPS_LANCZOS", &
684 : description="Threshold used for lanczos estimates.", &
685 235944 : default_r_val=1.0E-3_dp)
686 235944 : CALL section_add_keyword(section, keyword)
687 235944 : CALL keyword_release(keyword)
688 :
689 : CALL keyword_create(keyword, __LOCATION__, name="EPS_PGF_ORB", &
690 : description="Sets precision of the integral tensors.", &
691 : variants=["KP_EPS_PGF_ORB"], &
692 471888 : default_r_val=1.0E-5_dp)
693 235944 : CALL section_add_keyword(section, keyword)
694 235944 : CALL keyword_release(keyword)
695 :
696 : CALL keyword_create(keyword, __LOCATION__, name="MAX_ITER_LANCZOS", &
697 : description="Maximum number of lanczos iterations.", &
698 235944 : usage="MAX_ITER_LANCZOS ", default_i_val=500)
699 235944 : CALL section_add_keyword(section, keyword)
700 235944 : CALL keyword_release(keyword)
701 :
702 : CALL keyword_create(keyword, __LOCATION__, name="RI_FLAVOR", &
703 : description="Flavor of RI: how to contract 3-center integrals", &
704 : enum_c_vals=s2a("MO", "RHO"), &
705 : enum_desc=s2a("with MO coefficients", "with density matrix"), &
706 : enum_i_vals=[ri_mo, ri_pmat], &
707 235944 : default_i_val=ri_pmat)
708 235944 : CALL section_add_keyword(section, keyword)
709 235944 : CALL keyword_release(keyword)
710 :
711 : CALL keyword_create(keyword, __LOCATION__, name="MIN_BLOCK_SIZE", &
712 : description="Minimum tensor block size.", &
713 235944 : default_i_val=4)
714 235944 : CALL section_add_keyword(section, keyword)
715 235944 : CALL keyword_release(keyword)
716 :
717 : CALL keyword_create(keyword, __LOCATION__, name="MAX_BLOCK_SIZE_MO", &
718 : description="Maximum tensor block size for MOs.", &
719 235944 : default_i_val=64)
720 235944 : CALL section_add_keyword(section, keyword)
721 235944 : CALL keyword_release(keyword)
722 :
723 : CALL keyword_create(keyword, __LOCATION__, name="MEMORY_CUT", &
724 : description="Memory reduction factor. This keyword controls the batching of tensor "// &
725 : "contractions into smaller, more manageable chunks. The details vary "// &
726 : "depending on the RI_FLAVOR.", &
727 235944 : default_i_val=3)
728 235944 : CALL section_add_keyword(section, keyword)
729 235944 : CALL keyword_release(keyword)
730 :
731 : CALL keyword_create(keyword, __LOCATION__, name="FLAVOR_SWITCH_MEMORY_CUT", &
732 : description="Memory reduction factor to be applied upon RI_FLAVOR switching "// &
733 : "from MO to RHO. The RHO flavor typically requires more memory, "// &
734 : "and depending on the ressources available, a higher MEMORY_CUT.", &
735 235944 : default_i_val=3)
736 235944 : CALL section_add_keyword(section, keyword)
737 235944 : CALL keyword_release(keyword)
738 :
739 : CALL section_create(subsection, __LOCATION__, name="PRINT", &
740 : description="Section of possible print options in the RI-HFX code.", &
741 235944 : n_keywords=0, n_subsections=1, repeats=.FALSE.)
742 :
743 : CALL keyword_create(keyword, __LOCATION__, name="KP_RI_PROGRESS_BAR", &
744 : variants=s2a("PROGRESS_BAR", "PROGRESS", "KP_PROGRESS", "KP_PROGRESS_BAR"), &
745 : description="Whether a progress bar for individual SCF steps should be printed. "// &
746 : "In RI-HFXk, an expensive triple loop runs over periodic images and "// &
747 : "atomic pairs. This printing option tracks the progress of this loop "// &
748 : "in real time. Note that some work also takes place before the loop "// &
749 : "starts, and the time spent doing it depends on the value of "// &
750 : "KP_STACK_SIZE (larger = faster, but more memory used).", &
751 235944 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
752 235944 : CALL section_add_keyword(subsection, keyword)
753 235944 : CALL keyword_release(keyword)
754 :
755 : !TODO: add a incentive to restart from GGA calculation? Test that it improves things
756 : CALL keyword_create(keyword, __LOCATION__, name="KP_RI_MEMORY_ESTIMATE", &
757 : variants=s2a("MEMORY_ESTIMATE"), &
758 : description="Calculate and print a rough upper bound estimate of the memory "// &
759 : "required to run a RI-HFXk ENERGY calculation. Note that a fair "// &
760 : "amount of computing must take place before this estimate can be "// &
761 : "produced. If the calculation runs out of memory beforehand, "// &
762 : "use more resources, or change the value of memory related keywords "// &
763 : "(e.g. KP_NGROUPS, interaction potential parameters, KP_STACK_SIZE). "// &
764 : "The estimate is more accurate when restarting from a (GGA) wavefunction. "// &
765 : "Calculations involving forces will require more memory than this estimate.", &
766 235944 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
767 235944 : CALL section_add_keyword(subsection, keyword)
768 235944 : CALL keyword_release(keyword)
769 :
770 : CALL cp_print_key_section_create(print_key, __LOCATION__, "RI_INFO", &
771 : description="Controls the printing of DBCSR tensor log in RI HFX.", &
772 235944 : print_level=high_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
773 235944 : CALL section_add_subsection(subsection, print_key)
774 235944 : CALL section_release(print_key)
775 :
776 : CALL cp_print_key_section_create(print_key, __LOCATION__, "RI_DENSITY_COEFFS", &
777 : description="Controls the printing of the projection of the elecontric "// &
778 : "density on the RI_HFX basis. n(r) = sum_s c_s*phi_RI_s(r), "// &
779 : "where c_s = sum_pqr P_pq (pq|r) (r|s)^-1 and | is the RI_METRIC", &
780 : print_level=high_print_level, filename="RI_DENSITY_COEFFS", &
781 235944 : common_iter_levels=3)
782 :
783 : CALL keyword_create(keyword, __LOCATION__, name="MULTIPLY_BY_RI_2C_INTEGRALS", &
784 : variants=s2a("MULT_BY_RI", "MULT_BY_S", "MULT_BY_RI_INTS"), &
785 : description="Whether the RI density coefficients to be printed should "// &
786 : "be pre-multiplied by the RI_METRIC 2c-integrals: (r|s)*C_s. "// &
787 : "Not compatible with the SKIP_RI_METRIC keyword.", &
788 235944 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
789 235944 : CALL section_add_keyword(print_key, keyword)
790 235944 : CALL keyword_release(keyword)
791 :
792 : CALL keyword_create(keyword, __LOCATION__, name="SKIP_RI_METRIC", &
793 : variants=s2a("SKIP_INVERSE", "SKIP_2C_INTS", "SKIP_2C_INTEGRALS"), &
794 : description="Skip the calculation, inversion, and contraction of the 2-center RI "// &
795 : "metric integrals. The printed coefficients are only the contraction of the "// &
796 : "density matrix with the 3-center integrals, i.e. c_r = sum_pq P_pq (pq|r) "// &
797 : "Allows for memory savings when printing the RI density coefficients.", &
798 235944 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
799 235944 : CALL section_add_keyword(print_key, keyword)
800 235944 : CALL keyword_release(keyword)
801 :
802 : CALL keyword_create(keyword, __LOCATION__, name="FILE_FORMAT", &
803 : description="Format of file containing density fitting coefficients: "// &
804 : "BASIC(default)-original format; EXTENDED-format with basis set info.", &
805 235944 : default_c_val="BASIC")
806 235944 : CALL section_add_keyword(print_key, keyword)
807 235944 : CALL keyword_release(keyword)
808 :
809 235944 : CALL section_add_subsection(subsection, print_key)
810 235944 : CALL section_release(print_key)
811 :
812 : CALL cp_print_key_section_create(print_key, __LOCATION__, "RI_METRIC_2C_INTS", &
813 : description="Controls the printing of RI 2-center integrals for the "// &
814 : "HFX potential.", &
815 : print_level=high_print_level, filename="RI_2C_INTS", &
816 235944 : common_iter_levels=3)
817 235944 : CALL section_add_subsection(subsection, print_key)
818 235944 : CALL section_release(print_key)
819 :
820 235944 : CALL section_add_subsection(section, subsection)
821 235944 : CALL section_release(subsection)
822 :
823 235944 : END SUBROUTINE create_hf_ri_section
824 :
825 : END MODULE input_cp2k_hfx
|