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