37#include "./base/base_uses.f90"
42 LOGICAL,
PRIVATE,
PARAMETER :: debug_this_module = .false.
43 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'input_cp2k_hfx'
61 cpassert(.NOT.
ASSOCIATED(section))
63 description=
"Sets up the Hartree-Fock parameters if requested ", &
64 n_keywords=5, n_subsections=2, repeats=.true., &
67 NULLIFY (keyword, print_key, subsection)
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)
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.)
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.)
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)
102 description=
"Controls the printing basic info about hf method", &
107 CALL create_hf_pbc_section(subsection)
111 CALL create_hf_screening_section(subsection)
115 CALL create_hf_potential_section(subsection)
119 CALL create_hf_load_balance_section(subsection)
123 CALL create_hf_memory_section(subsection)
127 CALL create_hf_ri_section(subsection)
140 SUBROUTINE create_hf_load_balance_section(section)
146 cpassert(.NOT.
ASSOCIATED(section))
148 description=
"Parameters influencing the load balancing of the HF", &
149 n_keywords=1, n_subsections=0, repeats=.false., &
154 keyword, __location__, &
156 description=
"Number of bins per process used to group atom quartets.", &
163 keyword, __location__, &
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", &
174 keyword, __location__, &
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.)
186 description=
"Controls the printing of info about load balance", &
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.)
200 END SUBROUTINE create_hf_load_balance_section
209 SUBROUTINE create_hf_potential_section(section)
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., &
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)}$", &
240 "Truncated coulomb potential: if (r < R_c) 1/r else 0", &
241 "Truncated Mix coulomb and longrange potential, assumes/requires that the erf has fully decayed at R_c"), &
248 keyword, __location__, &
250 description=
"Parameter $\omega$ for short/longrange interaction", &
252 default_r_val=0.0_dp)
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)
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)
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)
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, &
289 keyword, __location__, &
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")
298 END SUBROUTINE create_hf_potential_section
307 SUBROUTINE create_hf_screening_section(section)
312 cpassert(.NOT.
ASSOCIATED(section))
314 description=
"Sets up screening parameters if requested ", &
315 n_keywords=1, n_subsections=0, repeats=.false., &
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 "// &
325 usage=
"EPS_SCHWARZ 1.0E-6", &
326 default_r_val=1.0e-10_dp)
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)
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.)
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.)
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)
374 END SUBROUTINE create_hf_screening_section
381 SUBROUTINE create_hf_pbc_section(section)
386 cpassert(.NOT.
ASSOCIATED(section))
388 description=
"Sets up periodic boundary condition parameters if requested ", &
389 n_keywords=1, n_subsections=0, repeats=.false., &
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", &
404 END SUBROUTINE create_hf_pbc_section
411 SUBROUTINE create_hf_memory_section(section)
416 cpassert(.NOT.
ASSOCIATED(section))
418 description=
"Sets up memory parameters for the storage of the ERI's if requested ", &
419 n_keywords=1, n_subsections=0, repeats=.false., &
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)
434 keyword, __location__, &
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", &
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", &
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", &
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.)
475 END SUBROUTINE create_hf_memory_section
481 SUBROUTINE create_hf_ri_section(section)
487 NULLIFY (keyword, print_key, subsection)
489 cpassert(.NOT.
ASSOCIATED(section))
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_")
495 CALL keyword_create(keyword, __location__, name=
"_SECTION_PARAMETERS_", &
496 description=
"controls the activation of RI", &
498 default_l_val=.false., &
499 lone_keyword_l_val=.true.)
504 description=
"Filter threshold for DBT tensor contraction.", &
505 variants=[
"KP_EPS_FILTER"], &
506 default_r_val=1.0e-09_dp)
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)
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)
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)
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, &
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, &
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)
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)
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}", &
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}", &
586 default_l_val=.true.)
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}", &
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, &
620 description=
"The type of RI operator. "// &
621 "Default is POTENTIAL_TYPE from INTERACTION_POTENTIAL. "// &
623 "Coulomb operator cannot be used in periodic systems.", &
624 usage=
"RI_METRIC {string}", &
626 variants=[
"KP_RI_METRIC"], &
628 enum_c_vals=
s2a(
"HFX",
"COULOMB",
"IDENTITY",
"TRUNCATED",
"SHORTRANGE"), &
629 enum_desc=
s2a(
"Same as HFX operator", &
630 "Standard Coulomb operator: 1/r", &
632 "Truncated Coulomb operator: 1/r if (r<R_c), 0 otherwise ", &
633 "Short range: erfc(omega*r)/r"), &
639 CALL keyword_create(keyword, __location__, name=
"2C_MATRIX_FUNCTIONS", &
640 description=
"Methods for matrix inverse and matrix square root.", &
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"), &
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)
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.)
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.)
676 description=
"Order of the iteration method for the calculation of "// &
677 "the sqrt of 2-center integral matrix.", &
683 description=
"Threshold used for lanczos estimates.", &
684 default_r_val=1.0e-3_dp)
689 description=
"Sets precision of the integral tensors.", &
690 variants=[
"KP_EPS_PGF_ORB"], &
691 default_r_val=1.0e-5_dp)
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)
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"), &
710 CALL keyword_create(keyword, __location__, name=
"MIN_BLOCK_SIZE", &
711 description=
"Minimum tensor block size.", &
716 CALL keyword_create(keyword, __location__, name=
"MAX_BLOCK_SIZE_MO", &
717 description=
"Maximum tensor block size for MOs.", &
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.", &
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.", &
739 description=
"Section of possible print options in the RI-HFX code.", &
740 n_keywords=0, n_subsections=1, repeats=.false.)
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.)
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.)
770 description=
"Controls the printing of DBCSR tensor log in RI HFX.", &
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", &
780 common_iter_levels=3)
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.)
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.)
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")
812 description=
"Controls the printing of RI 2-center integrals for the "// &
815 common_iter_levels=3)
822 END SUBROUTINE create_hf_ri_section
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
Defines the basic variable types.
integer, parameter, public dp
Utilities for string manipulations.