38#include "./base/base_uses.f90"
43 LOGICAL,
PRIVATE,
PARAMETER :: debug_this_module = .false.
44 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'input_cp2k_hfx'
62 cpassert(.NOT.
ASSOCIATED(section))
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., &
69 NULLIFY (keyword, print_key, subsection)
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)
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.)
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.)
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)
104 description=
"Controls the printing basic info about hf method", &
109 CALL create_hf_pbc_section(subsection)
113 CALL create_hf_screening_section(subsection)
117 CALL create_hf_potential_section(subsection)
121 CALL create_hf_load_balance_section(subsection)
125 CALL create_hf_memory_section(subsection)
129 CALL create_hf_ri_section(subsection)
133 CALL create_hf_ace_section(subsection)
146 SUBROUTINE create_hf_load_balance_section(section)
152 cpassert(.NOT.
ASSOCIATED(section))
154 description=
"Parameters influencing the load balancing of the HF", &
155 n_keywords=1, n_subsections=0, repeats=.false., &
160 keyword, __location__, &
162 description=
"Number of bins per process used to group atom quartets.", &
169 keyword, __location__, &
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", &
180 keyword, __location__, &
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.)
192 description=
"Controls the printing of info about load balance", &
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.)
206 END SUBROUTINE create_hf_load_balance_section
215 SUBROUTINE create_hf_potential_section(section)
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., &
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)}$", &
247 "Truncated coulomb potential: if (r < R_c) 1/r else 0", &
248 "Truncated Mix coulomb and longrange potential, assumes/requires that the erf has fully decayed at R_c"), &
255 keyword, __location__, &
257 description=
"Parameter $\omega$ for short/longrange interaction", &
259 default_r_val=0.0_dp)
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)
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)
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)
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, &
297 keyword, __location__, &
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")
306 END SUBROUTINE create_hf_potential_section
315 SUBROUTINE create_hf_screening_section(section)
320 cpassert(.NOT.
ASSOCIATED(section))
322 description=
"Controls screening thresholds for Hartree-Fock exchange integrals.", &
323 n_keywords=1, n_subsections=0, repeats=.false., &
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)
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)
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.)
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.)
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)
380 END SUBROUTINE create_hf_screening_section
387 SUBROUTINE create_hf_pbc_section(section)
392 cpassert(.NOT.
ASSOCIATED(section))
394 description=
"Sets up periodic boundary condition parameters if requested ", &
395 n_keywords=1, n_subsections=0, repeats=.false., &
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", &
410 END SUBROUTINE create_hf_pbc_section
417 SUBROUTINE create_hf_memory_section(section)
422 cpassert(.NOT.
ASSOCIATED(section))
424 description=
"Sets up memory parameters for the storage of the ERI's if requested ", &
425 n_keywords=1, n_subsections=0, repeats=.false., &
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)
440 keyword, __location__, &
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", &
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", &
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", &
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.)
481 END SUBROUTINE create_hf_memory_section
487 SUBROUTINE create_hf_ri_section(section)
493 NULLIFY (keyword, print_key, subsection)
495 cpassert(.NOT.
ASSOCIATED(section))
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_")
501 CALL keyword_create(keyword, __location__, name=
"_SECTION_PARAMETERS_", &
502 description=
"controls the activation of RI", &
504 default_l_val=.false., &
505 lone_keyword_l_val=.true.)
510 description=
"Filter threshold for DBT tensor contraction.", &
511 variants=[
"KP_EPS_FILTER"], &
512 default_r_val=1.0e-09_dp)
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)
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)
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)
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, &
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, &
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)
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)
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}", &
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}", &
592 default_l_val=.true.)
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}", &
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, &
626 description=
"The type of RI operator. "// &
627 "Default is POTENTIAL_TYPE from INTERACTION_POTENTIAL. "// &
629 "Coulomb operator cannot be used in periodic systems.", &
630 usage=
"RI_METRIC {string}", &
632 variants=[
"KP_RI_METRIC"], &
634 enum_c_vals=
s2a(
"HFX",
"COULOMB",
"IDENTITY",
"TRUNCATED",
"SHORTRANGE"), &
635 enum_desc=
s2a(
"Same as HFX operator", &
636 "Standard Coulomb operator: 1/r", &
638 "Truncated Coulomb operator: 1/r if (r<R_c), 0 otherwise ", &
639 "Short range: erfc(omega*r)/r"), &
645 CALL keyword_create(keyword, __location__, name=
"2C_MATRIX_FUNCTIONS", &
646 description=
"Methods for matrix inverse and matrix square root.", &
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"), &
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)
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.)
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.)
682 description=
"Order of the iteration method for the calculation of "// &
683 "the sqrt of 2-center integral matrix.", &
689 description=
"Threshold used for lanczos estimates.", &
690 default_r_val=1.0e-3_dp)
695 description=
"Sets precision of the integral tensors.", &
696 variants=[
"KP_EPS_PGF_ORB"], &
697 default_r_val=1.0e-5_dp)
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)
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"), &
716 CALL keyword_create(keyword, __location__, name=
"MIN_BLOCK_SIZE", &
717 description=
"Minimum tensor block size.", &
722 CALL keyword_create(keyword, __location__, name=
"MAX_BLOCK_SIZE_MO", &
723 description=
"Maximum tensor block size for MOs.", &
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.", &
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.", &
745 description=
"Section of possible print options in the RI-HFX code.", &
746 n_keywords=0, n_subsections=1, repeats=.false.)
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.)
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.)
776 description=
"Controls the printing of DBCSR tensor log in RI HFX.", &
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", &
786 common_iter_levels=3)
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.)
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.)
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")
818 description=
"Controls the printing of RI 2-center integrals for the "// &
821 common_iter_levels=3)
828 END SUBROUTINE create_hf_ri_section
832 SUBROUTINE create_hf_ace_section(section)
837 cpassert(.NOT.
ASSOCIATED(section))
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.", &
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.)
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", &
878 END SUBROUTINE create_hf_ace_section
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
Defines the basic variable types.
integer, parameter, public dp
Utilities for string manipulations.