61#include "./base/base_uses.f90"
66 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'input_cp2k_qs'
83 cpassert(.NOT.
ASSOCIATED(section))
85 description=
"parameters needed to set up the Quickstep framework", &
86 n_keywords=1, n_subsections=0, repeats=.false.)
88 NULLIFY (keyword, subsection)
92 description=
"Try setting all EPS_xxx to values leading to an energy correct up to EPS_DEFAULT", &
93 usage=
"EPS_DEFAULT real", default_r_val=1.0e-10_dp)
97 CALL keyword_create(keyword, __location__, name=
"EPS_CORE_CHARGE", &
98 description=
"Precision for mapping the core charges.Overrides EPS_DEFAULT/100.0 value", &
99 usage=
"EPS_CORE_CHARGE real", type_of_var=
real_t)
104 keyword, __location__, name=
"EPS_GVG_RSPACE", &
105 variants=[
"EPS_GVG"], &
106 description=
"Sets precision of the realspace KS matrix element integration. Overrides SQRT(EPS_DEFAULT) value", &
107 usage=
"EPS_GVG_RSPACE real", type_of_var=
real_t)
112 description=
"Sets precision of the overlap matrix elements. Overrides SQRT(EPS_DEFAULT) value", &
113 usage=
"EPS_PGF_ORB real", type_of_var=
real_t)
118 keyword, __location__, name=
"EPS_KG_ORB", &
119 description=
"Sets precision used in coloring the subsets for the Kim-Gordon method. Overrides SQRT(EPS_DEFAULT) value", &
120 usage=
"EPS_KG_ORB 1.0E-8", &
126 description=
"Adjusts the precision for the local part of the pseudo potential. ", &
127 usage=
"EPS_PPL real", type_of_var=
real_t, default_r_val=1.0e-2_dp)
132 keyword, __location__, name=
"EPS_PPNL", &
133 description=
"Sets precision of the non-local part of the pseudo potential. Overrides sqrt(EPS_DEFAULT) value", &
134 usage=
"EPS_PPNL real", type_of_var=
real_t)
139 description=
"Sets precision of the GAPW projection. Overrides EPS_DEFAULT value", &
140 usage=
"EPS_CPC real", type_of_var=
real_t)
145 description=
"Sets precision of the density mapping on the grids.Overrides EPS_DEFAULT value", &
146 usage=
"EPS_RHO real", type_of_var=
real_t)
150 CALL keyword_create(keyword, __location__, name=
"EPS_RHO_RSPACE", &
151 description=
"Sets precision of the density mapping in rspace.Overrides EPS_DEFAULT value."// &
152 " Overrides EPS_RHO value", &
153 usage=
"EPS_RHO_RSPACE real", type_of_var=
real_t)
157 CALL keyword_create(keyword, __location__, name=
"EPS_RHO_GSPACE", &
158 description=
"Sets precision of the density mapping in gspace.Overrides EPS_DEFAULT value."// &
159 " Overrides EPS_RHO value", &
160 usage=
"EPS_RHO_GSPACE real", type_of_var=
real_t)
164 CALL keyword_create(keyword, __location__, name=
"EPS_FILTER_MATRIX", &
165 description=
"Sets the threshold for filtering matrix elements.", &
166 usage=
"EPS_FILTER_MATRIX 1.0E-6", type_of_var=
real_t, default_r_val=0.0e0_dp)
171 variants=[
"EPS_FIT"], &
172 description=
"GAPW: tolerance controlling the split of Gaussian basis functions "// &
173 "into hard atom-centered and soft grid-expanded parts. Smaller values include "// &
174 "harder Gaussians in the soft density and can require a larger MGRID cutoff.", &
175 usage=
"EPSFIT real", default_r_val=1.0e-4_dp)
180 variants=[
"EPS_ISO"], &
181 description=
"GAPW: precision to determine an isolated projector", &
182 usage=
"EPSISO real", default_r_val=1.0e-12_dp)
187 variants=[
"EPS_SVD"], &
188 description=
"GAPW: tolerance used in the singular value decomposition of the "// &
189 "projector matrix. Smaller values can improve numerical accuracy at increased cost.", &
190 usage=
"EPS_SVD real", default_r_val=1.0e-8_dp)
195 variants=
s2a(
"EPSVRHO0",
"EPS_VRHO0"), &
196 description=
"GAPW: tolerance used to determine the range of the "// &
197 "V(rho0-rho0_soft) compensation contribution.", &
198 usage=
"EPSRHO0 real", default_r_val=1.0e-6_dp)
203 variants=
s2a(
"ALPHA0_H",
"ALPHA0"), &
204 description=
"GAPW: Exponent for hard compensation charge", &
205 usage=
"ALPHA0_HARD real", default_r_val=0.0_dp)
210 keyword, __location__, name=
"FORCE_PAW", &
211 description=
"Use the GAPW scheme also for atoms with soft basis sets, i.e. "// &
212 "the local densities are computed even if hard and soft should be equal. "// &
213 "If this keyword is not set to true, those atoms with soft basis sets are treated by a GPW scheme, i.e. "// &
214 "the corresponding density contribution goes on the global grid and is expanded in PW. "// &
215 "This option nullifies the effect of the GPW_TYPE in the atomic KIND", &
217 default_l_val=.false., lone_keyword_l_val=.true.)
221 CALL keyword_create(keyword, __location__, name=
"MAX_RAD_LOCAL", &
222 description=
"GAPW : maximum radius of gaussian functions"// &
223 " included in the generation of projectors", &
224 usage=
"MAX_RAD_LOCAL real", default_r_val=25.0_dp)
228 CALL keyword_create(keyword, __location__, name=
"GAPW_1C_BASIS", &
229 description=
"Specifies how to construct the GAPW one center basis set. "// &
230 "Default is to use the primitives from the orbital basis.", &
231 usage=
"GAPW_1C_BASIS MEDIUM", &
232 enum_c_vals=
s2a(
"ORB",
"EXT_SMALL",
"EXT_MEDIUM",
"EXT_LARGE",
"EXT_VERY_LARGE"), &
233 enum_desc=
s2a(
"Use orbital basis set.", &
234 "Extension using Small number of primitive Gaussians.", &
235 "Extension using Medium number of primitive Gaussians.", &
236 "Extension using Large number of primitive Gaussians.", &
237 "Extension using Very Large number of primitive Gaussians."), &
245 keyword, __location__, name=
"GAPW_ACCURATE_XCINT", &
246 description=
"Use the accurate GAPW/GAPW_XC XC integration scheme for one-center hard/soft "// &
247 "density differences. This opt-in path covers regular GAPW/GAPW_XC XC energy, potential, "// &
248 "forces, mGGA/tau terms, NLCC, Fine-XC grids, local XC energy-density transfer, analytical "// &
249 "stress, TDDFPT/response forces, ADMM-GAPW force paths, an ADMM-GAPW diagonal stress "// &
250 "debug path, nonlocal vdW smoke coverage, representative k-point, XAS_TDP, and RTBSE "// &
251 "smoke cases, and KG GAPW/GAPW_XC EMBED, EMBED_RI, ATOMIC, and NONE cases. Regular-grid "// &
252 "local energy and stress cube print keys keep their existing soft-grid semantics and are "// &
253 "not changed by this keyword. The default remains off while this coverage is being prepared "// &
254 "for a future default change.", &
255 usage=
"GAPW_ACCURATE_XCINT", &
256 default_l_val=.false., lone_keyword_l_val=.true.)
261 keyword, __location__, name=
"ALPHA_WEIGHTS", &
262 description=
"Gaussian exponent reference (rc=1.2 Bohr) for accurate integration in GAPW.", &
263 usage=
"ALPHA_WEIGHTS 10.0", default_r_val=6.0_dp)
267 CALL keyword_create(keyword, __location__, name=
"MIN_PAIR_LIST_RADIUS", &
268 description=
"Set the minimum value [Bohr] for the overlap pair list radius."// &
269 " Default is 0.0 Bohr, negative values are changed to the cell size."// &
270 " This allows to control the sparsity of the KS matrix for HFX calculations.", &
271 usage=
"MIN_PAIR_LIST_RADIUS real", default_r_val=0.0_dp)
277 description=
"Perform a linear scaling SCF", &
278 usage=
"LS_SCF", lone_keyword_l_val=.true., &
279 default_l_val=.false.)
284 description=
"Perform ALMO SCF", &
285 usage=
"ALMO_SCF", lone_keyword_l_val=.true., &
286 default_l_val=.false.)
291 description=
"Perform transport calculations (coupling CP2K and OMEN)", &
292 usage=
"TRANSPORT", lone_keyword_l_val=.true., &
293 default_l_val=.false.)
298 description=
"Use a Kim-Gordon-like scheme.", &
299 usage=
"KG_METHOD", lone_keyword_l_val=.true., &
304 CALL keyword_create(keyword, __location__, name=
"REF_EMBED_SUBSYS", &
305 description=
"A total, reference, system in DFT embedding. ", &
306 usage=
"REF_EMBED_SUBSYS FALSE", &
307 default_l_val=.false., lone_keyword_l_val=.true.)
311 CALL keyword_create(keyword, __location__, name=
"CLUSTER_EMBED_SUBSYS", &
312 description=
"A cluster treated with DFT in DFT embedding. ", &
313 usage=
"CLUSTER_EMBED_SUBSYS FALSE", &
314 default_l_val=.false., lone_keyword_l_val=.true.)
318 CALL keyword_create(keyword, __location__, name=
"HIGH_LEVEL_EMBED_SUBSYS", &
319 description=
"A cluster treated with a high-level method in DFT embedding. ", &
320 usage=
"HIGH_LEVEL_EMBED_SUBSYS FALSE", &
321 default_l_val=.false., lone_keyword_l_val=.true.)
325 CALL keyword_create(keyword, __location__, name=
"DFET_EMBEDDED", &
326 description=
"Calculation with DFT-embedding potential. ", &
327 usage=
"DFET_EMBEDDED FALSE", &
328 default_l_val=.false., lone_keyword_l_val=.true.)
332 CALL keyword_create(keyword, __location__, name=
"DMFET_EMBEDDED", &
333 description=
"Calculation with DM embedding potential. ", &
334 usage=
"DMFET_EMBEDDED FALSE", &
335 default_l_val=.false., lone_keyword_l_val=.true.)
341 description=
"Order of Gaussian type expansion of Slater orbital basis sets.", &
342 usage=
"STO_NG", default_i_val=6)
347 variants=[
"LMAXRHO1"], &
348 description=
"GAPW : max L number for expansion of the atomic densities in spherical gaussians", &
349 usage=
"LMAXN1 integer", &
355 variants=[
"LMAXRHO0"], &
356 description=
"GAPW : max L number for the expansion compensation densities in spherical gaussians", &
357 usage=
"LMAXN0 integer", &
363 description=
"GAPW : integer added to the max L of the basis set, used to determine the "// &
364 "maximum value of L for the compensation charge density.", &
365 usage=
"LADDN0 integer", &
372 description=
"GAPW: algorithm to construct the atomic radial grids", &
373 usage=
"QUADRATURE GC_SIMPLE", &
374 enum_c_vals=
s2a(
"GC_SIMPLE",
"GC_TRANSFORMED",
"GC_LOG"), &
376 enum_desc=
s2a(
"Gauss-Chebyshev quadrature", &
377 "Transformed Gauss-Chebyshev quadrature", &
378 "Logarithmic transformed Gauss-Chebyshev quadrature"), &
384 description=
"What kind of PW_GRID should be employed", &
385 usage=
"PW_GRID NS-FULLSPACE", &
386 enum_c_vals=
s2a(
"SPHERICAL",
"NS-FULLSPACE",
"NS-HALFSPACE"), &
387 enum_desc=
s2a(
"- not tested",
" tested",
" - not tested"), &
393 CALL keyword_create(keyword, __location__, name=
"PW_GRID_LAYOUT", &
394 description=
"Force a particular real-space layout for the plane waves grids. "// &
395 "Numbers ≤ 0 mean that this dimension is free, incorrect layouts will be ignored. "// &
396 "The default (/-1,-1/) causes CP2K to select a good value, "// &
397 "i.e. plane distributed for large grids, more general distribution for small grids.", &
398 usage=
"PW_GRID_LAYOUT 4 16", &
399 repeats=.false., n_var=2, &
400 default_i_vals=[-1, -1])
404 CALL keyword_create(keyword, __location__, name=
"PW_GRID_BLOCKED", &
405 description=
"Can be used to set the distribution in g-space for the pw grids and their FFT.", &
406 usage=
"PW_GRID_BLOCKED FREE", &
407 enum_c_vals=
s2a(
"FREE",
"TRUE",
"FALSE"), &
408 enum_desc=
s2a(
"CP2K will select an appropriate value",
"blocked",
"not blocked"), &
415 keyword, __location__, name=
"EXTRAPOLATION", &
416 variants=
s2a(
"INTERPOLATION",
"WF_INTERPOLATION"), &
417 description=
"Extrapolation strategy for the wavefunction during e.g. MD. "// &
418 "Not all options are available for all simulation methods. "// &
419 "PS and ASPC are recommended, see also EXTRAPOLATION_ORDER.", &
421 usage=
"EXTRAPOLATION PS", &
422 enum_c_vals=
s2a(
"USE_GUESS",
"USE_PREV_P",
"USE_PREV_RHO_R",
"LINEAR_WF", &
423 "LINEAR_P",
"LINEAR_PS",
"USE_PREV_WF",
"PS",
"FROZEN",
"ASPC", &
424 "GEXT_PROJ",
"GEXT_PROJ_QTR"), &
426 "Use the method specified with SCF_GUESS, i.e. no extrapolation", &
427 "Use the previous density matrix", &
428 "Use the previous density in real space", &
429 "Linear extrapolation of the wavefunction (not available for k-points)", &
430 "Linear extrapolation of the density matrix", &
431 "Linear extrapolation of the density matrix times the overlap matrix (not available for k-points)", &
432 "Use the previous wavefunction", &
433 "Higher order extrapolation of the density matrix times the overlap matrix", &
434 "Frozen ... (not available for k-points)", &
435 "Always stable predictor corrector, similar to PS, but going for MD stability instead of initial guess accuracy.", &
436 "GExt extrapolation for the density matrix times the overlap matrix.", &
437 "Quasi time reversible GExt extrapolation for the density matrix times the overlap matrix."), &
455 CALL keyword_create(keyword, __location__, name=
"EXTRAPOLATION_ORDER", &
456 description=
"Order for the PS, ASPC extrapolation (typically 2-4) or "// &
457 "order for the GEXT_PROJ, GEXT_PROJ_QTR extrapolation (typically 4-10). "// &
458 "Higher order might bring more accuracy, but comes, "// &
459 "for large systems, also at some cost. "// &
460 "In some cases, a high order extrapolation is not stable,"// &
461 " and the order needs to be reduced.", &
462 usage=
"EXTRAPOLATION_ORDER {integer}", default_i_val=3)
467 description=
"Specifies the electronic structure method that should be employed", &
468 usage=
"METHOD GAPW", &
469 enum_c_vals=
s2a(
"GAPW",
"GAPW_XC",
"GPW",
"LRIGPW",
"RIGPW", &
470 "MNDO",
"MNDOD",
"AM1",
"PM3",
"PM6",
"PM6-FM",
"PDG",
"RM1",
"PNNL",
"DFTB",
"xTB",
"OFGPW"), &
471 enum_desc=
s2a(
"Gaussian and augmented plane waves method", &
472 "Gaussian and augmented plane waves method only for XC", &
473 "Gaussian and plane waves method", &
474 "Local resolution of identity method", &
475 "Resolution of identity method for HXC terms", &
476 "MNDO semiempirical",
"MNDO-d semiempirical",
"AM1 semiempirical", &
477 "PM3 semiempirical",
"PM6 semiempirical",
"PM6-FM semiempirical",
"PDG semiempirical", &
478 "RM1 semiempirical", &
479 "PNNL semiempirical", &
480 "DFTB Density Functional based Tight-Binding", &
481 "GFN-xTB Extended Tight-Binding", &
482 "OFGPW Orbital-free GPW method"), &
495 description=
"Specifies the method used to calculate the local pseudopotential contribution.", &
496 usage=
"CORE_PPL ANALYTIC", &
497 enum_c_vals=
s2a(
"ANALYTIC",
"GRID"), &
498 enum_desc=
s2a(
"Analytic integration of integrals", &
499 "Numerical integration on real space grid. Lumped together with core charge"), &
505 CALL keyword_create(keyword, __location__, name=
"EMBED_RESTART_FILE_NAME", &
506 description=
"Root of the file name where to read the embedding "// &
507 "potential guess.", &
508 usage=
"EMBED_RESTART_FILE_NAME <FILENAME>", &
513 CALL keyword_create(keyword, __location__, name=
"EMBED_CUBE_FILE_NAME", &
514 description=
"Root of the file name where to read the embedding "// &
515 "potential (guess) as a cube.", &
516 usage=
"EMBED_CUBE_FILE_NAME <FILENAME>", &
521 CALL keyword_create(keyword, __location__, name=
"EMBED_SPIN_CUBE_FILE_NAME", &
522 description=
"Root of the file name where to read the spin part "// &
523 "of the embedding potential (guess) as a cube.", &
524 usage=
"EMBED_SPIN_CUBE_FILE_NAME <FILENAME>", &
545 CALL create_mulliken_section(subsection)
557 CALL create_s2_restraint_section(subsection)
584 SUBROUTINE create_mulliken_section(section)
590 cpassert(.NOT.
ASSOCIATED(section))
591 CALL section_create(section, __location__, name=
"MULLIKEN_RESTRAINT", &
592 description=
"Use mulliken charges in a restraint (check code for details)", &
593 n_keywords=7, n_subsections=0, repeats=.false.)
596 description=
"force constant of the restraint", &
597 usage=
"STRENGTH {real} ", default_r_val=0.1_dp)
602 description=
"target value of the restraint", &
603 usage=
"TARGET {real} ", default_r_val=1._dp)
608 description=
"Specifies the list of atoms that is summed in the restraint", &
609 usage=
"ATOMS {integer} {integer} .. {integer}", &
614 END SUBROUTINE create_mulliken_section
623 CHARACTER(len=*),
INTENT(in) :: section_name
628 NULLIFY (keyword, print_key)
629 cpassert(.NOT.
ASSOCIATED(section))
630 CALL section_create(section, __location__, name=trim(adjustl(section_name)), &
631 description=
"Use DDAPC charges in a restraint (check code for details)", &
632 n_keywords=7, n_subsections=0, repeats=.true.)
634 CALL keyword_create(keyword, __location__, name=
"TYPE_OF_DENSITY", &
635 description=
"Specifies the type of density used for the fitting", &
636 usage=
"TYPE_OF_DENSITY (FULL|SPIN)", &
637 enum_c_vals=
s2a(
"FULL",
"SPIN"), &
639 enum_desc=
s2a(
"Full density",
"Spin density"), &
645 description=
"force constant of the restraint", &
646 usage=
"STRENGTH {real} ", default_r_val=0.1_dp)
651 description=
"target value of the restraint", &
652 usage=
"TARGET {real} ", default_r_val=1._dp)
657 description=
"Specifies the list of atoms that is summed in the restraint", &
658 usage=
"ATOMS {integer} {integer} .. {integer}", &
664 description=
"Defines the the coefficient of the atom in the atom list (default is one) ", &
665 usage=
"COEFF 1.0 -1.0", &
666 type_of_var=
real_t, n_var=-1)
670 CALL keyword_create(keyword, __location__, name=
"FUNCTIONAL_FORM", &
671 description=
"Specifies the functional form of the term added", &
672 usage=
"FUNCTIONAL_FORM RESTRAINT", &
673 enum_c_vals=
s2a(
"RESTRAINT",
"CONSTRAINT"), &
675 enum_desc=
s2a(
"Harmonic potential: s*(q-t)**2",
"Constraint form: s*(q-t)"), &
681 description=
"Controls the printing basic info about the method", &
692 SUBROUTINE create_s2_restraint_section(section)
698 cpassert(.NOT.
ASSOCIATED(section))
701 description=
"Use S2 in a re/constraint (OT only)", &
702 n_keywords=7, n_subsections=0, repeats=.false.)
705 description=
"force constant of the restraint", &
706 usage=
"STRENGTH {real} ", default_r_val=0.1_dp)
711 description=
"target value of the restraint", &
712 usage=
"TARGET {real} ", default_r_val=1._dp)
716 CALL keyword_create(keyword, __location__, name=
"FUNCTIONAL_FORM", &
717 description=
"Specifies the functional form of the term added", &
718 usage=
"FUNCTIONAL_FORM RESTRAINT", &
719 enum_c_vals=
s2a(
"RESTRAINT",
"CONSTRAINT"), &
721 enum_desc=
s2a(
"Harmonic potential: s*(q-t)**2",
"Constraint form: s*(q-t)"), &
726 END SUBROUTINE create_s2_restraint_section
739 cpassert(.NOT.
ASSOCIATED(section))
741 description=
"This section specifies optional parameters for LRIGPW.", &
742 n_keywords=3, n_subsections=0, repeats=.false., citations=[
golze2017b])
746 CALL keyword_create(keyword, __location__, name=
"LRI_OVERLAP_MATRIX", &
747 description=
"Specifies whether to calculate the inverse or the "// &
748 "pseudoinverse of the overlap matrix of the auxiliary "// &
749 "basis set. Calculating the pseudoinverse is necessary "// &
750 "for very large auxiliary basis sets, but more expensive. "// &
751 "Using the pseudoinverse, consistent forces are not "// &
753 usage=
"LRI_OVERLAP_MATRIX INVERSE", &
754 enum_c_vals=
s2a(
"INVERSE",
"PSEUDO_INVERSE_SVD",
"PSEUDO_INVERSE_DIAG", &
756 enum_desc=
s2a(
"Calculate inverse of the overlap matrix.", &
757 "Calculate the pseuodinverse of the overlap matrix "// &
758 "using singular value decomposition.", &
759 "Calculate the pseudoinverse of the overlap matrix "// &
760 "by prior diagonalization.", &
761 "Choose automatically for each pair whether to "// &
762 "calculate the inverse or pseudoinverse based on the "// &
763 "condition number of the overlap matrix for each pair. "// &
764 "Calculating the pseudoinverse is much more expensive."), &
771 CALL keyword_create(keyword, __location__, name=
"MAX_CONDITION_NUM", &
772 description=
"If AUTOSELECT is chosen for LRI_OVERLAP_MATRIX, this "// &
773 "keyword specifies that the pseudoinverse is calculated "// &
774 "only if the LOG of the condition number of the lri "// &
775 "overlap matrix is larger than the given value.", &
776 usage=
"MAX_CONDITION_NUM 20.0", default_r_val=20.0_dp)
781 description=
"Threshold for ABA and ABB integrals in LRI. "// &
782 "This is used for screening in the KS and "// &
783 "force calculations (tensor contractions).", &
784 usage=
"EPS_O3_INT 1.e-10", default_r_val=1.0e-14_dp)
788 CALL keyword_create(keyword, __location__, name=
"DEBUG_LRI_INTEGRALS", &
789 description=
"Debug the integrals needed for LRIGPW.", &
790 usage=
"DEBUG_LRI_INTEGRALS TRUE", &
791 default_l_val=.false., lone_keyword_l_val=.true.)
795 CALL keyword_create(keyword, __location__, name=
"EXACT_1C_TERMS", &
796 description=
"Don't use LRI for one center densities.", &
797 usage=
"EXACT_1C_TERMS TRUE", &
798 default_l_val=.false., lone_keyword_l_val=.true.)
803 description=
"Use LRI/RI for local pseudopotential.", &
804 usage=
"PPL_RI TRUE", &
805 default_l_val=.false., lone_keyword_l_val=.true.)
810 description=
"Print statistical information on the RI calculation.", &
811 usage=
"RI_STATISTIC TRUE", &
812 default_l_val=.false., lone_keyword_l_val=.true.)
816 CALL keyword_create(keyword, __location__, name=
"DISTANT_PAIR_APPROXIMATION", &
817 description=
"Calculate distant pairs using an independent atom approximation.", &
818 usage=
"DISTANT_PAIR_APPROXIMATION TRUE", &
819 default_l_val=.false., lone_keyword_l_val=.true.)
823 CALL keyword_create(keyword, __location__, name=
"DISTANT_PAIR_METHOD", &
824 description=
"Method used to separate pair density for distant pairs. "// &
825 "Options: EW (equal weights); AW (atomic weights); SW (set weights); "// &
826 "LW (shell function weights)", &
827 usage=
"DISTANT_PAIR_METHOD {method}", &
832 CALL keyword_create(keyword, __location__, name=
"DISTANT_PAIR_RADII", &
833 description=
"Inner and outer radii used in distant "// &
834 "pair separation. Smooth interpolation between inner and outer "// &
836 usage=
"DISTANT_PAIR_RADII r_inner {real} r_outer {real} ", &
837 n_var=2, default_r_vals=[8._dp, 12._dp], unit_str=
'bohr', &
842 CALL keyword_create(keyword, __location__, name=
"SHG_LRI_INTEGRALS", &
843 description=
"Uses the SHG (solid harmonic Gaussian) integral "// &
844 "scheme instead of Obara-Saika", &
845 usage=
"SHG_LRI_INTEGRALS TRUE", &
846 default_l_val=.true., lone_keyword_l_val=.true., &
852 description=
"Approximation to be used for the inverse of the "// &
853 "RI overlap matrix. INVF, INVS: exact inverse, apply directly "// &
854 "for solver (F:full matrix, S:sparsematrix). AINV approximate inverse, use with PCG. "// &
855 "NONE: no approximation used with CG solver.", &
856 usage=
"RI_SINV NONE", default_c_val=
"INVF")
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public stewart2007
integer, save, public vandevondele2006
integer, save, public golze2017b
integer, save, public vandevondele2005a
integer, save, public kuhne2007
integer, save, public golze2017a
integer, save, public lippert1999
integer, save, public dewar1977
integer, save, public vanvoorhis2015
integer, save, public repasky2002
integer, save, public iannuzzi2006
integer, save, public rocha2006
integer, save, public lippert1997
integer, save, public andermatt2016
integer, save, public thiel1992
integer, save, public schenter2008
integer, save, public brelaz1979
integer, save, public krack2000
integer, save, public dewar1985
integer, save, public stewart1989
integer, save, public kolafa2004
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer, parameter, public low_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
This module defines the grid data type and some basic operations on it.
integer, parameter, public do_pw_grid_blocked_false
integer, parameter, public do_pw_grid_blocked_true
integer, parameter, public do_pw_grid_blocked_free
Utilities for string manipulations.