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: 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)", &
239 "Truncated coulomb potential: if (r < R_c) 1/r else 0", &
240 "Truncated Mix coulomb and longrange potential, assumes/requires that the erf has fully decayed at R_c"), &
247 keyword, __location__, &
249 description=
"Parameter for short/longrange interaction", &
251 default_r_val=0.0_dp)
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)
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)
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)
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, &
285 keyword, __location__, &
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")
294 END SUBROUTINE create_hf_potential_section
303 SUBROUTINE create_hf_screening_section(section)
308 cpassert(.NOT.
ASSOCIATED(section))
310 description=
"Sets up screening parameters if requested ", &
311 n_keywords=1, n_subsections=0, repeats=.false., &
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 "// &
321 usage=
"EPS_SCHWARZ 1.0E-6", &
322 default_r_val=1.0e-10_dp)
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)
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.)
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.)
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)
370 END SUBROUTINE create_hf_screening_section
377 SUBROUTINE create_hf_pbc_section(section)
382 cpassert(.NOT.
ASSOCIATED(section))
384 description=
"Sets up periodic boundary condition parameters if requested ", &
385 n_keywords=1, n_subsections=0, repeats=.false., &
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", &
400 END SUBROUTINE create_hf_pbc_section
407 SUBROUTINE create_hf_memory_section(section)
412 cpassert(.NOT.
ASSOCIATED(section))
414 description=
"Sets up memory parameters for the storage of the ERI's if requested ", &
415 n_keywords=1, n_subsections=0, repeats=.false., &
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)
430 keyword, __location__, &
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", &
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", &
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", &
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.)
471 END SUBROUTINE create_hf_memory_section
477 SUBROUTINE create_hf_ri_section(section)
483 NULLIFY (keyword, print_key, subsection)
485 cpassert(.NOT.
ASSOCIATED(section))
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_")
491 CALL keyword_create(keyword, __location__, name=
"_SECTION_PARAMETERS_", &
492 description=
"controls the activation of RI", &
494 default_l_val=.false., &
495 lone_keyword_l_val=.true.)
500 description=
"Filter threshold for DBT tensor contraction.", &
501 variants=(/
"KP_EPS_FILTER"/), &
502 default_r_val=1.0e-09_dp)
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)
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)
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)
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, &
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, &
547 description=
"The number of MPI subgroup that work in parallel during the SCF. "// &
548 "The default value is 1. Using N subgroups should speed up the "// &
549 "calculation by a factor ~N, at the cost of N times more memory usage.", &
550 variants=(/
"NGROUPS"/), &
551 usage=
"KP_NGROUPS {int}", &
557 CALL keyword_create(keyword, __location__, name=
"KP_USE_DELTA_P", &
558 description=
"This kweyword controls whether the KS matrix at each SCF cycle "// &
559 "is built by adding the contribution of the denisty difference (wrt to previous step) "// &
560 "to the KS matrix of the previous step. As the SCF converges, the density fluctuations "// &
561 "get smaller and sparsity increases, leading to faster SCF steps. Not always "// &
562 "numerically stable => turn off if SCF struggles to converge.", &
563 variants=
s2a(
"USE_DELTA_P",
"KP_USE_P_DIFF",
"USE_P_DIFF"), &
564 usage=
"KP_USE_DELTA_P {logical}", &
566 default_l_val=.true.)
570 CALL keyword_create(keyword, __location__, name=
"KP_STACK_SIZE", &
571 description=
"When doing contraction over periodic cells of the type: "// &
572 "T_mu^a,nu^b,P^c = (mu^a nu^b | Q^d) * (Q^d | P^c), with "// &
573 "a,b,c,d labeling cells, there are in principle Ncells "// &
574 "contractions taking place. Because a smaller number of "// &
575 "contractions involving larger tensors is more efficient, "// &
576 "the tensors can be stacked along the d direction. STCK_SIZE "// &
577 "controls the size of this stack. Larger stacks are more efficient, "// &
578 "but required more memory.", &
579 variants=(/
"STACK_SIZE"/), &
580 usage=
"KP_STACK_SIZE {int}", &
586 CALL keyword_create(keyword, __location__, name=
"KP_RI_BUMP_FACTOR", &
587 variants=
s2a(
"RI_BUMP",
"BUMP",
"BUMP_FACTOR"), &
588 description=
"In KP-RI-HFX, the extended RI basis set has a bump radius. "// &
589 "All basis elements within that radius contribute with full weight. "// &
590 "All basis elements beyond that radius have decaying weight, from "// &
591 "1 at the bump radius, to zero at the RI extension radius. The "// &
592 "bump radius is calculated as a fraction of the RI extension radius: "// &
593 "bump radius = KP_RI_NUMP_FACTOR * RI extension radius", &
594 default_r_val=0.85_dp, &
600 description=
"The type of RI operator. "// &
601 "Default is POTENTIAL_TYPE from INTERACTION_POTENTIAL. "// &
603 "Coulomb operator cannot be used in periodic systems.", &
604 usage=
"RI_METRIC {string}", &
606 variants=(/
"KP_RI_METRIC"/), &
608 enum_c_vals=
s2a(
"HFX",
"COULOMB",
"IDENTITY",
"TRUNCATED",
"SHORTRANGE"), &
609 enum_desc=
s2a(
"Same as HFX operator", &
610 "Standard Coulomb operator: 1/r", &
612 "Truncated Coulomb operator: 1/r if (r<R_c), 0 otherwise ", &
613 "Short range: erfc(omega*r)/r"), &
619 CALL keyword_create(keyword, __location__, name=
"2C_MATRIX_FUNCTIONS", &
620 description=
"Methods for matrix inverse and matrix square root.", &
622 enum_c_vals=
s2a(
"DIAG",
"CHOLESKY",
"ITER"), &
623 enum_desc=
s2a(
"Diagonalization with eigenvalue quenching: stable", &
624 "Cholesky: not stable in case of ill-conditioned RI basis", &
625 "Iterative algorithms: linear scaling "// &
626 "Hotelling's method for inverse and Newton-Schulz iteration for matrix square root"), &
632 description=
"Throw away linear combinations of RI basis functions with a small eigenvalue, "// &
633 "this is applied only if 2C_MATRIX_FUNCTIONS DIAG", &
634 default_r_val=1.0e-7_dp)
638 CALL keyword_create(keyword, __location__, name=
"CHECK_2C_MATRIX", &
639 description=
"Report accuracy for the inverse/sqrt of the 2-center integral matrix.", &
640 default_l_val=.false., lone_keyword_l_val=.true.)
645 keyword, __location__, &
646 name=
"CALC_COND_NUM", &
647 variants=(/
"CALC_CONDITION_NUMBER"/), &
648 description=
"Calculate the condition number of integral matrices.", &
649 usage=
"CALC_COND_NUM", &
650 default_l_val=.false., &
651 lone_keyword_l_val=.true.)
656 description=
"Order of the iteration method for the calculation of "// &
657 "the sqrt of 2-center integral matrix.", &
663 description=
"Threshold used for lanczos estimates.", &
664 default_r_val=1.0e-3_dp)
669 description=
"Sets precision of the integral tensors.", &
670 variants=(/
"KP_EPS_PGF_ORB"/), &
671 default_r_val=1.0e-5_dp)
675 CALL keyword_create(keyword, __location__, name=
"MAX_ITER_LANCZOS", &
676 description=
"Maximum number of lanczos iterations.", &
677 usage=
"MAX_ITER_LANCZOS ", default_i_val=500)
682 description=
"Flavor of RI: how to contract 3-center integrals", &
683 enum_c_vals=
s2a(
"MO",
"RHO"), &
684 enum_desc=
s2a(
"with MO coefficients",
"with density matrix"), &
690 CALL keyword_create(keyword, __location__, name=
"MIN_BLOCK_SIZE", &
691 description=
"Minimum tensor block size.", &
696 CALL keyword_create(keyword, __location__, name=
"MAX_BLOCK_SIZE_MO", &
697 description=
"Maximum tensor block size for MOs.", &
703 description=
"Memory reduction factor. This keyword controls the batching of tensor "// &
704 "contractions into smaller, more manageable chunks. The details vary "// &
705 "depending on the RI_FLAVOR.", &
710 CALL keyword_create(keyword, __location__, name=
"FLAVOR_SWITCH_MEMORY_CUT", &
711 description=
"Memory reduction factor to be applied upon RI_FLAVOR switching "// &
712 "from MO to RHO. The RHO flavor typically requires more memory, "// &
713 "and depending on the ressources available, a higher MEMORY_CUT.", &
719 description=
"Section of possible print options in the RI-HFX code.", &
720 n_keywords=0, n_subsections=1, repeats=.false.)
723 description=
"Controls the printing of DBCSR tensor log in RI HFX.", &
729 description=
"Controls the printing of the projection of the elecontric "// &
730 "density on the RI_HFX basis. n(r) = sum_s c_s*phi_RI_s(r), "// &
731 "where c_s = sum_pqr P_pq (pq|r) (r|s)^-1 and | is the RI_METRIC", &
733 common_iter_levels=3)
735 CALL keyword_create(keyword, __location__, name=
"MULTIPLY_BY_RI_2C_INTEGRALS", &
736 variants=
s2a(
"MULT_BY_RI",
"MULT_BY_S",
"MULT_BY_RI_INTS"), &
737 description=
"Whether the RI density coefficients to be printed should "// &
738 "be pre-multiplied by the RI_METRIC 2c-integrals: (r|s)*C_s.", &
739 default_l_val=.false., lone_keyword_l_val=.true.)
747 description=
"Controls the printing of RI 2-center integrals for the "// &
750 common_iter_levels=3)
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.