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: precision to give the extension of a hard gaussian ", &
173 usage=
"EPSFIT real", default_r_val=1.0e-4_dp)
178 variants=[
"EPS_ISO"], &
179 description=
"GAPW: precision to determine an isolated projector", &
180 usage=
"EPSISO real", default_r_val=1.0e-12_dp)
185 variants=[
"EPS_SVD"], &
186 description=
"GAPW: tolerance used in the singular value decomposition of the projector matrix", &
187 usage=
"EPS_SVD real", default_r_val=1.0e-8_dp)
192 variants=
s2a(
"EPSVRHO0",
"EPS_VRHO0"), &
193 description=
"GAPW : precision to determine the range of V(rho0-rho0soft)", &
194 usage=
"EPSRHO0 real", default_r_val=1.0e-6_dp)
199 variants=
s2a(
"ALPHA0_H",
"ALPHA0"), &
200 description=
"GAPW: Exponent for hard compensation charge", &
201 usage=
"ALPHA0_HARD real", default_r_val=0.0_dp)
206 keyword, __location__, name=
"FORCE_PAW", &
207 description=
"Use the GAPW scheme also for atoms with soft basis sets, i.e. "// &
208 "the local densities are computed even if hard and soft should be equal. "// &
209 "If this keyword is not set to true, those atoms with soft basis sets are treated by a GPW scheme, i.e. "// &
210 "the corresponding density contribution goes on the global grid and is expanded in PW. "// &
211 "This option nullifies the effect of the GPW_TYPE in the atomic KIND", &
213 default_l_val=.false., lone_keyword_l_val=.true.)
217 CALL keyword_create(keyword, __location__, name=
"MAX_RAD_LOCAL", &
218 description=
"GAPW : maximum radius of gaussian functions"// &
219 " included in the generation of projectors", &
220 usage=
"MAX_RAD_LOCAL real", default_r_val=25.0_dp)
224 CALL keyword_create(keyword, __location__, name=
"GAPW_1C_BASIS", &
225 description=
"Specifies how to construct the GAPW one center basis set. "// &
226 "Default is to use the primitives from the orbital basis.", &
227 usage=
"GAPW_1C_BASIS MEDIUM", &
228 enum_c_vals=
s2a(
"ORB",
"EXT_SMALL",
"EXT_MEDIUM",
"EXT_LARGE",
"EXT_VERY_LARGE"), &
229 enum_desc=
s2a(
"Use orbital basis set.", &
230 "Extension using Small number of primitive Gaussians.", &
231 "Extension using Medium number of primitive Gaussians.", &
232 "Extension using Large number of primitive Gaussians.", &
233 "Extension using Very Large number of primitive Gaussians."), &
241 keyword, __location__, name=
"GAPW_ACCURATE_XCINT", &
242 description=
"Use a more accurate integration scheme for the soft XC energy in GAPW.", &
243 usage=
"GAPW_ACCURATE_XCINT", &
244 default_l_val=.false., lone_keyword_l_val=.true.)
249 keyword, __location__, name=
"ALPHA_WEIGHTS", &
250 description=
"Gaussian exponent reference (rc=1.2 Bohr) for accurate integration in GAPW.", &
251 usage=
"ALPHA_WEIGHTS 10.0", default_r_val=6.0_dp)
255 CALL keyword_create(keyword, __location__, name=
"MIN_PAIR_LIST_RADIUS", &
256 description=
"Set the minimum value [Bohr] for the overlap pair list radius."// &
257 " Default is 0.0 Bohr, negative values are changed to the cell size."// &
258 " This allows to control the sparsity of the KS matrix for HFX calculations.", &
259 usage=
"MIN_PAIR_LIST_RADIUS real", default_r_val=0.0_dp)
265 description=
"Perform a linear scaling SCF", &
266 usage=
"LS_SCF", lone_keyword_l_val=.true., &
267 default_l_val=.false.)
272 description=
"Perform ALMO SCF", &
273 usage=
"ALMO_SCF", lone_keyword_l_val=.true., &
274 default_l_val=.false.)
279 description=
"Perform transport calculations (coupling CP2K and OMEN)", &
280 usage=
"TRANSPORT", lone_keyword_l_val=.true., &
281 default_l_val=.false.)
286 description=
"Use a Kim-Gordon-like scheme.", &
287 usage=
"KG_METHOD", lone_keyword_l_val=.true., &
292 CALL keyword_create(keyword, __location__, name=
"REF_EMBED_SUBSYS", &
293 description=
"A total, reference, system in DFT embedding. ", &
294 usage=
"REF_EMBED_SUBSYS FALSE", &
295 default_l_val=.false., lone_keyword_l_val=.true.)
299 CALL keyword_create(keyword, __location__, name=
"CLUSTER_EMBED_SUBSYS", &
300 description=
"A cluster treated with DFT in DFT embedding. ", &
301 usage=
"CLUSTER_EMBED_SUBSYS FALSE", &
302 default_l_val=.false., lone_keyword_l_val=.true.)
306 CALL keyword_create(keyword, __location__, name=
"HIGH_LEVEL_EMBED_SUBSYS", &
307 description=
"A cluster treated with a high-level method in DFT embedding. ", &
308 usage=
"HIGH_LEVEL_EMBED_SUBSYS FALSE", &
309 default_l_val=.false., lone_keyword_l_val=.true.)
313 CALL keyword_create(keyword, __location__, name=
"DFET_EMBEDDED", &
314 description=
"Calculation with DFT-embedding potential. ", &
315 usage=
"DFET_EMBEDDED FALSE", &
316 default_l_val=.false., lone_keyword_l_val=.true.)
320 CALL keyword_create(keyword, __location__, name=
"DMFET_EMBEDDED", &
321 description=
"Calculation with DM embedding potential. ", &
322 usage=
"DMFET_EMBEDDED FALSE", &
323 default_l_val=.false., lone_keyword_l_val=.true.)
329 description=
"Order of Gaussian type expansion of Slater orbital basis sets.", &
330 usage=
"STO_NG", default_i_val=6)
335 variants=[
"LMAXRHO1"], &
336 description=
"GAPW : max L number for expansion of the atomic densities in spherical gaussians", &
337 usage=
"LMAXN1 integer", &
343 variants=[
"LMAXRHO0"], &
344 description=
"GAPW : max L number for the expansion compensation densities in spherical gaussians", &
345 usage=
"LMAXN0 integer", &
351 description=
"GAPW : integer added to the max L of the basis set, used to determine the "// &
352 "maximum value of L for the compensation charge density.", &
353 usage=
"LADDN0 integer", &
360 description=
"GAPW: algorithm to construct the atomic radial grids", &
361 usage=
"QUADRATURE GC_SIMPLE", &
362 enum_c_vals=
s2a(
"GC_SIMPLE",
"GC_TRANSFORMED",
"GC_LOG"), &
364 enum_desc=
s2a(
"Gauss-Chebyshev quadrature", &
365 "Transformed Gauss-Chebyshev quadrature", &
366 "Logarithmic transformed Gauss-Chebyshev quadrature"), &
372 description=
"What kind of PW_GRID should be employed", &
373 usage=
"PW_GRID NS-FULLSPACE", &
374 enum_c_vals=
s2a(
"SPHERICAL",
"NS-FULLSPACE",
"NS-HALFSPACE"), &
375 enum_desc=
s2a(
"- not tested",
" tested",
" - not tested"), &
381 CALL keyword_create(keyword, __location__, name=
"PW_GRID_LAYOUT", &
382 description=
"Force a particular real-space layout for the plane waves grids. "// &
383 "Numbers ≤ 0 mean that this dimension is free, incorrect layouts will be ignored. "// &
384 "The default (/-1,-1/) causes CP2K to select a good value, "// &
385 "i.e. plane distributed for large grids, more general distribution for small grids.", &
386 usage=
"PW_GRID_LAYOUT 4 16", &
387 repeats=.false., n_var=2, &
388 default_i_vals=[-1, -1])
392 CALL keyword_create(keyword, __location__, name=
"PW_GRID_BLOCKED", &
393 description=
"Can be used to set the distribution in g-space for the pw grids and their FFT.", &
394 usage=
"PW_GRID_BLOCKED FREE", &
395 enum_c_vals=
s2a(
"FREE",
"TRUE",
"FALSE"), &
396 enum_desc=
s2a(
"CP2K will select an appropriate value",
"blocked",
"not blocked"), &
403 keyword, __location__, name=
"EXTRAPOLATION", &
404 variants=
s2a(
"INTERPOLATION",
"WF_INTERPOLATION"), &
405 description=
"Extrapolation strategy for the wavefunction during e.g. MD. "// &
406 "Not all options are available for all simulation methods. "// &
407 "PS and ASPC are recommended, see also EXTRAPOLATION_ORDER.", &
409 usage=
"EXTRAPOLATION PS", &
410 enum_c_vals=
s2a(
"USE_GUESS",
"USE_PREV_P",
"USE_PREV_RHO_R",
"LINEAR_WF", &
411 "LINEAR_P",
"LINEAR_PS",
"USE_PREV_WF",
"PS",
"FROZEN",
"ASPC"), &
413 "Use the method specified with SCF_GUESS, i.e. no extrapolation", &
414 "Use the previous density matrix", &
415 "Use the previous density in real space", &
416 "Linear extrapolation of the wavefunction (not available for k-points)", &
417 "Linear extrapolation of the density matrix", &
418 "Linear extrapolation of the density matrix times the overlap matrix (not available for k-points)", &
419 "Use the previous wavefunction", &
420 "Higher order extrapolation of the density matrix times the overlap matrix", &
421 "Frozen ... (not available for k-points)", &
422 "Always stable predictor corrector, similar to PS, but going for MD stability instead of initial guess accuracy."), &
438 CALL keyword_create(keyword, __location__, name=
"EXTRAPOLATION_ORDER", &
439 description=
"Order for the PS or ASPC extrapolation (typically 2-4). "// &
440 "Higher order might bring more accuracy, but comes, "// &
441 "for large systems, also at some cost. "// &
442 "In some cases, a high order extrapolation is not stable,"// &
443 " and the order needs to be reduced.", &
444 usage=
"EXTRAPOLATION_ORDER {integer}", default_i_val=3)
449 description=
"Specifies the electronic structure method that should be employed", &
450 usage=
"METHOD GAPW", &
451 enum_c_vals=
s2a(
"GAPW",
"GAPW_XC",
"GPW",
"LRIGPW",
"RIGPW", &
452 "MNDO",
"MNDOD",
"AM1",
"PM3",
"PM6",
"PM6-FM",
"PDG",
"RM1",
"PNNL",
"DFTB",
"xTB",
"OFGPW"), &
453 enum_desc=
s2a(
"Gaussian and augmented plane waves method", &
454 "Gaussian and augmented plane waves method only for XC", &
455 "Gaussian and plane waves method", &
456 "Local resolution of identity method", &
457 "Resolution of identity method for HXC terms", &
458 "MNDO semiempirical",
"MNDO-d semiempirical",
"AM1 semiempirical", &
459 "PM3 semiempirical",
"PM6 semiempirical",
"PM6-FM semiempirical",
"PDG semiempirical", &
460 "RM1 semiempirical", &
461 "PNNL semiempirical", &
462 "DFTB Density Functional based Tight-Binding", &
463 "GFN-xTB Extended Tight-Binding", &
464 "OFGPW Orbital-free GPW method"), &
477 description=
"Specifies the method used to calculate the local pseudopotential contribution.", &
478 usage=
"CORE_PPL ANALYTIC", &
479 enum_c_vals=
s2a(
"ANALYTIC",
"GRID"), &
480 enum_desc=
s2a(
"Analytic integration of integrals", &
481 "Numerical integration on real space grid. Lumped together with core charge"), &
487 CALL keyword_create(keyword, __location__, name=
"EMBED_RESTART_FILE_NAME", &
488 description=
"Root of the file name where to read the embedding "// &
489 "potential guess.", &
490 usage=
"EMBED_RESTART_FILE_NAME <FILENAME>", &
495 CALL keyword_create(keyword, __location__, name=
"EMBED_CUBE_FILE_NAME", &
496 description=
"Root of the file name where to read the embedding "// &
497 "potential (guess) as a cube.", &
498 usage=
"EMBED_CUBE_FILE_NAME <FILENAME>", &
503 CALL keyword_create(keyword, __location__, name=
"EMBED_SPIN_CUBE_FILE_NAME", &
504 description=
"Root of the file name where to read the spin part "// &
505 "of the embedding potential (guess) as a cube.", &
506 usage=
"EMBED_SPIN_CUBE_FILE_NAME <FILENAME>", &
527 CALL create_mulliken_section(subsection)
539 CALL create_s2_restraint_section(subsection)
566 SUBROUTINE create_mulliken_section(section)
572 cpassert(.NOT.
ASSOCIATED(section))
573 CALL section_create(section, __location__, name=
"MULLIKEN_RESTRAINT", &
574 description=
"Use mulliken charges in a restraint (check code for details)", &
575 n_keywords=7, n_subsections=0, repeats=.false.)
578 description=
"force constant of the restraint", &
579 usage=
"STRENGTH {real} ", default_r_val=0.1_dp)
584 description=
"target value of the restraint", &
585 usage=
"TARGET {real} ", default_r_val=1._dp)
590 description=
"Specifies the list of atoms that is summed in the restraint", &
591 usage=
"ATOMS {integer} {integer} .. {integer}", &
596 END SUBROUTINE create_mulliken_section
605 CHARACTER(len=*),
INTENT(in) :: section_name
610 NULLIFY (keyword, print_key)
611 cpassert(.NOT.
ASSOCIATED(section))
612 CALL section_create(section, __location__, name=trim(adjustl(section_name)), &
613 description=
"Use DDAPC charges in a restraint (check code for details)", &
614 n_keywords=7, n_subsections=0, repeats=.true.)
616 CALL keyword_create(keyword, __location__, name=
"TYPE_OF_DENSITY", &
617 description=
"Specifies the type of density used for the fitting", &
618 usage=
"TYPE_OF_DENSITY (FULL|SPIN)", &
619 enum_c_vals=
s2a(
"FULL",
"SPIN"), &
621 enum_desc=
s2a(
"Full density",
"Spin density"), &
627 description=
"force constant of the restraint", &
628 usage=
"STRENGTH {real} ", default_r_val=0.1_dp)
633 description=
"target value of the restraint", &
634 usage=
"TARGET {real} ", default_r_val=1._dp)
639 description=
"Specifies the list of atoms that is summed in the restraint", &
640 usage=
"ATOMS {integer} {integer} .. {integer}", &
646 description=
"Defines the the coefficient of the atom in the atom list (default is one) ", &
647 usage=
"COEFF 1.0 -1.0", &
648 type_of_var=
real_t, n_var=-1)
652 CALL keyword_create(keyword, __location__, name=
"FUNCTIONAL_FORM", &
653 description=
"Specifies the functional form of the term added", &
654 usage=
"FUNCTIONAL_FORM RESTRAINT", &
655 enum_c_vals=
s2a(
"RESTRAINT",
"CONSTRAINT"), &
657 enum_desc=
s2a(
"Harmonic potential: s*(q-t)**2",
"Constraint form: s*(q-t)"), &
663 description=
"Controls the printing basic info about the method", &
674 SUBROUTINE create_s2_restraint_section(section)
680 cpassert(.NOT.
ASSOCIATED(section))
683 description=
"Use S2 in a re/constraint (OT only)", &
684 n_keywords=7, n_subsections=0, repeats=.false.)
687 description=
"force constant of the restraint", &
688 usage=
"STRENGTH {real} ", default_r_val=0.1_dp)
693 description=
"target value of the restraint", &
694 usage=
"TARGET {real} ", default_r_val=1._dp)
698 CALL keyword_create(keyword, __location__, name=
"FUNCTIONAL_FORM", &
699 description=
"Specifies the functional form of the term added", &
700 usage=
"FUNCTIONAL_FORM RESTRAINT", &
701 enum_c_vals=
s2a(
"RESTRAINT",
"CONSTRAINT"), &
703 enum_desc=
s2a(
"Harmonic potential: s*(q-t)**2",
"Constraint form: s*(q-t)"), &
708 END SUBROUTINE create_s2_restraint_section
721 cpassert(.NOT.
ASSOCIATED(section))
723 description=
"This section specifies optional parameters for LRIGPW.", &
724 n_keywords=3, n_subsections=0, repeats=.false., citations=[
golze2017b])
728 CALL keyword_create(keyword, __location__, name=
"LRI_OVERLAP_MATRIX", &
729 description=
"Specifies whether to calculate the inverse or the "// &
730 "pseudoinverse of the overlap matrix of the auxiliary "// &
731 "basis set. Calculating the pseudoinverse is necessary "// &
732 "for very large auxiliary basis sets, but more expensive. "// &
733 "Using the pseudoinverse, consistent forces are not "// &
735 usage=
"LRI_OVERLAP_MATRIX INVERSE", &
736 enum_c_vals=
s2a(
"INVERSE",
"PSEUDO_INVERSE_SVD",
"PSEUDO_INVERSE_DIAG", &
738 enum_desc=
s2a(
"Calculate inverse of the overlap matrix.", &
739 "Calculate the pseuodinverse of the overlap matrix "// &
740 "using singular value decomposition.", &
741 "Calculate the pseudoinverse of the overlap matrix "// &
742 "by prior diagonalization.", &
743 "Choose automatically for each pair whether to "// &
744 "calculate the inverse or pseudoinverse based on the "// &
745 "condition number of the overlap matrix for each pair. "// &
746 "Calculating the pseudoinverse is much more expensive."), &
753 CALL keyword_create(keyword, __location__, name=
"MAX_CONDITION_NUM", &
754 description=
"If AUTOSELECT is chosen for LRI_OVERLAP_MATRIX, this "// &
755 "keyword specifies that the pseudoinverse is calculated "// &
756 "only if the LOG of the condition number of the lri "// &
757 "overlap matrix is larger than the given value.", &
758 usage=
"MAX_CONDITION_NUM 20.0", default_r_val=20.0_dp)
763 description=
"Threshold for ABA and ABB integrals in LRI. "// &
764 "This is used for screening in the KS and "// &
765 "force calculations (tensor contractions).", &
766 usage=
"EPS_O3_INT 1.e-10", default_r_val=1.0e-14_dp)
770 CALL keyword_create(keyword, __location__, name=
"DEBUG_LRI_INTEGRALS", &
771 description=
"Debug the integrals needed for LRIGPW.", &
772 usage=
"DEBUG_LRI_INTEGRALS TRUE", &
773 default_l_val=.false., lone_keyword_l_val=.true.)
777 CALL keyword_create(keyword, __location__, name=
"EXACT_1C_TERMS", &
778 description=
"Don't use LRI for one center densities.", &
779 usage=
"EXACT_1C_TERMS TRUE", &
780 default_l_val=.false., lone_keyword_l_val=.true.)
785 description=
"Use LRI/RI for local pseudopotential.", &
786 usage=
"PPL_RI TRUE", &
787 default_l_val=.false., lone_keyword_l_val=.true.)
792 description=
"Print statistical information on the RI calculation.", &
793 usage=
"RI_STATISTIC TRUE", &
794 default_l_val=.false., lone_keyword_l_val=.true.)
798 CALL keyword_create(keyword, __location__, name=
"DISTANT_PAIR_APPROXIMATION", &
799 description=
"Calculate distant pairs using an independent atom approximation.", &
800 usage=
"DISTANT_PAIR_APPROXIMATION TRUE", &
801 default_l_val=.false., lone_keyword_l_val=.true.)
805 CALL keyword_create(keyword, __location__, name=
"DISTANT_PAIR_METHOD", &
806 description=
"Method used to separate pair density for distant pairs. "// &
807 "Options: EW (equal weights); AW (atomic weights); SW (set weights); "// &
808 "LW (shell function weights)", &
809 usage=
"DISTANT_PAIR_METHOD {method}", &
814 CALL keyword_create(keyword, __location__, name=
"DISTANT_PAIR_RADII", &
815 description=
"Inner and outer radii used in distant "// &
816 "pair separation. Smooth interpolation between inner and outer "// &
818 usage=
"DISTANT_PAIR_RADII r_inner {real} r_outer {real} ", &
819 n_var=2, default_r_vals=[8._dp, 12._dp], unit_str=
'bohr', &
824 CALL keyword_create(keyword, __location__, name=
"SHG_LRI_INTEGRALS", &
825 description=
"Uses the SHG (solid harmonic Gaussian) integral "// &
826 "scheme instead of Obara-Saika", &
827 usage=
"SHG_LRI_INTEGRALS TRUE", &
828 default_l_val=.true., lone_keyword_l_val=.true., &
834 description=
"Approximation to be used for the inverse of the "// &
835 "RI overlap matrix. INVF, INVS: exact inverse, apply directly "// &
836 "for solver (F:full matrix, S:sparsematrix). AINV approximate inverse, use with PCG. "// &
837 "NONE: no approximation used with CG solver.", &
838 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.