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