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