84#include "./base/base_uses.f90"
89 LOGICAL,
PRIVATE,
PARAMETER :: debug_this_module = .true.
90 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'input_cp2k_poisson'
109 cpassert(.NOT.
ASSOCIATED(section))
111 description=
"Controls the Poisson solver and electrostatic boundary conditions used by DFT.", &
112 n_keywords=1, n_subsections=0, repeats=.false.)
114 NULLIFY (keyword, subsection)
115 CALL keyword_create(keyword, __location__, name=
"POISSON_SOLVER", &
116 variants=[
"POISSON",
"PSOLVER"], &
117 description=
"Specify which kind of solver to use to solve the Poisson equation.", &
118 usage=
"POISSON_SOLVER char", &
119 enum_c_vals=
s2a(
"PERIODIC",
"ANALYTIC",
"MT",
"MULTIPOLE",
"WAVELET",
"IMPLICIT"), &
122 enum_desc=
s2a(
"PERIODIC is only available for fully (3D) periodic systems.", &
123 "ANALYTIC is available for 0D, 1D and 2D periodic solutions using analytical green "// &
124 "functions in the g space (slow convergence).", &
125 "MT (Martyna Tuckermann) decoupling that interacts only with the nearest "// &
126 "neighbor. Beware results are completely wrong if the cell is smaller than twice the "// &
127 "cluster size (with electronic density). Available for 0D and 2D systems.", &
128 "MULTIPOLE uses a scheme that fits the total charge with one gaussian per atom. "// &
129 "Available only for cluster (0D) systems.", &
130 "WAVELET allows for 0D, 2D (but only PERIODIC XZ) and 3D systems. It does not "// &
131 "require very large unit cells, only that the density goes to zero on the faces of "// &
132 "the cell. The use of PREFERRED_FFT_LIBRARY FFTSG is required.", &
133 "IMPLICIT allows for 0D, 1D, 2D and 3D systems."), &
140 description=
"Specifies the directions in which periodic boundary conditions apply to electrostatics. "// &
141 "See the CELL section for the periodicity used by geometry and pair lists; "// &
142 "the settings are usually the same.", &
143 usage=
"PERIODIC (x|y|z|xy|xz|yz|xyz|none)", &
144 enum_c_vals=
s2a(
"x",
"y",
"z",
"xy",
"xz",
"yz",
"xyz",
"none"), &
152 CALL create_mt_section(subsection)
156 CALL create_wavelet_section(subsection)
160 CALL create_multipole_section(subsection)
168 CALL create_implicit_ps_section(subsection)
178 SUBROUTINE create_multipole_section(section)
184 cpassert(.NOT.
ASSOCIATED(section))
187 description=
"This section is used to set up the decoupling of QM periodic images with "// &
188 "the use of density derived atomic point charges.", &
189 n_keywords=1, n_subsections=0, repeats=.false.)
191 NULLIFY (keyword, subsection)
193 description=
"Real space cutoff for the Ewald sum.", &
194 usage=
"RCUT {real}", n_var=1, type_of_var=
real_t, &
199 CALL keyword_create(keyword, __location__, name=
"EWALD_PRECISION", &
200 description=
"Precision achieved in the Ewald sum.", &
201 usage=
"EWALD_PRECISION {real}", n_var=1, type_of_var=
real_t, &
202 unit_str=
"hartree", default_r_val=1.0e-6_dp)
206 CALL keyword_create(keyword, __location__, name=
"ANALYTICAL_GTERM", &
207 description=
"Evaluates the Gterm in the Ewald Scheme analytically instead of using Splines.", &
208 usage=
"ANALYTICAL_GTERM <LOGICAL>", &
209 default_l_val=.false., lone_keyword_l_val=.true.)
214 description=
"Specifies the number of grid points used for the Interpolation of the G-space term", &
215 usage=
"NGRIDS <integer> <integer> <integer> ", n_var=3, default_i_vals=[50, 50, 50])
224 description=
"Controls the checking of the G-space term Spline Interpolation.", &
230 description=
"Controls the printing of basic information during the run", &
235 END SUBROUTINE create_multipole_section
242 SUBROUTINE create_mt_section(section)
247 cpassert(.NOT.
ASSOCIATED(section))
249 description=
"Sets up parameters of Martyna-Tuckerman poisson solver. "// &
250 "Note that exact results are only guaranteed if the unit cell is "// &
251 "twice as large as charge density (and serious artefacts can result "// &
252 "if the cell is much smaller).", &
253 n_keywords=1, n_subsections=0, repeats=.false., &
259 description=
"Convergence parameter ALPHA*RMIN. Default value 7.0", &
260 usage=
"ALPHA real", &
261 n_var=1, default_r_val=7.0_dp)
266 description=
"Specify the multiplicative factor for the CUTOFF keyword in MULTI_GRID"// &
267 " section. The result gives the cutoff at which the 1/r non-periodic FFT3D is evaluated."// &
269 usage=
"REL_CUTOFF real", &
270 n_var=1, default_r_val=2.0_dp)
274 END SUBROUTINE create_mt_section
287 cpassert(.NOT.
ASSOCIATED(section))
289 description=
"Ewald parameters controlling electrostatic only for CLASSICAL MM.", &
290 n_keywords=7, n_subsections=0, repeats=.false., &
293 NULLIFY (keyword, print_key, subsection)
295 keyword, __location__, name=
"EWALD_TYPE", &
296 description=
"The type of ewald you want to perform.", &
298 usage=
"EWALD_TYPE (NONE|EWALD|PME|SPME)", &
300 enum_c_vals=[
"none ", &
308 enum_desc=
s2a(
"NONE standard real-space coulomb potential is computed together with the non-bonded contributions", &
309 "EWALD is the standard non-fft based ewald", &
310 "PME is the particle mesh using fft interpolation", &
311 "SPME is the smooth particle mesh using beta-Euler splines (recommended)"))
315 CALL keyword_create(keyword, __location__, name=
"EWALD_ACCURACY", &
316 description=
"Expected accuracy in the Ewald sum. This number affects only the calculation of "// &
317 "the cutoff for the real-space term of the ewald summation (EWALD|PME|SPME) as well as the "// &
318 "construction of the neighbor lists (if the cutoff for non-bonded terms is smaller than the "// &
319 "value employed to compute the EWALD real-space term). This keyword has no "// &
320 "effect on the reciprocal space term (which can be tuned independently).", &
321 usage=
"EWALD_ACCURACY {real}", n_var=1, type_of_var=
real_t, &
322 unit_str=
"hartree", default_r_val=1.0e-6_dp)
327 description=
"Explicitly provide the real-space cutoff of the ewald summation (EWALD|PME|SPME). "// &
328 "This value is ignored in Tight-binding applications (rcut from basis overlap is used). "// &
329 "If present, overwrites the estimate of EWALD_ACCURACY and may affect the "// &
330 "construction of the neighbor lists for non-bonded terms (in FIST), if the value "// &
331 "specified is larger than the cutoff for non-bonded interactions.", &
332 usage=
"RCUT 5.0", n_var=1, type_of_var=
real_t, unit_str=
"angstrom")
337 description=
"alpha parameter associated with Ewald (EWALD|PME|SPME). "// &
338 "Recommended for small systems is alpha = 3.5 / r_cut. "// &
339 "For Tight-binding application a recommended value is alpha = 1.0. "// &
340 "Tuning alpha, r_cut and gmax is needed to obtain O(N**1.5) scaling for ewald.", &
342 default_r_val=
cp_unit_to_cp2k(
value=0.35_dp, unit_str=
"angstrom^-1"), &
343 unit_str=
'angstrom^-1')
348 description=
"number of grid points (SPME and EWALD). If a single number is specified, "// &
349 "the same number of points is used for all three directions on the grid. "// &
350 "If three numbers are given, each direction can have a different number of points. "// &
351 "The number of points needs to be FFTable (which depends on the library used) and odd for EWALD. "// &
352 "The optimal number depends e.g. on alpha and the size of the cell. 1 point per Angstrom is common.", &
353 usage=
"gmax 25 25 25", n_var=-1, type_of_var=
integer_t)
358 description=
"number of grid points on small mesh (PME only), should be odd.", &
359 usage=
"ns_max 11", default_i_val=11)
364 description=
"order of the beta-Euler spline (SPME only)", &
365 usage=
"o_spline 6", default_i_val=6)
370 description=
"tolerance of gaussians for fft interpolation (PME only)", &
371 usage=
"epsilon 1e-6", default_r_val=1.e-6_dp)
381 CALL section_create(subsection, __location__, name=
"MULTIPOLES", &
382 description=
"Enables the use of multipoles in the treatment of the electrostatics.", &
383 n_keywords=0, n_subsections=1, repeats=.false., &
386 CALL keyword_create(keyword, __location__, name=
"_SECTION_PARAMETERS_", &
387 description=
"Controls the activation of the Multipoles", &
388 usage=
"&MULTIPOLES T", default_l_val=.false., lone_keyword_l_val=.true.)
392 CALL keyword_create(keyword, __location__, name=
"MAX_MULTIPOLE_EXPANSION", &
393 description=
"Specify the maximum level of multipoles expansion used "// &
394 "for the electrostatics.", &
395 usage=
"MAX_MULTIPOLE_EXPANSION DIPOLE", &
396 enum_c_vals=
s2a(
"NONE",
"CHARGE",
"DIPOLE",
"QUADRUPOLE"), &
397 enum_desc=
s2a(
"No multipolar terms! Check the codes providing a zero contribution.", &
398 "Use up to the Charge term", &
399 "Use up to the Dipole term", &
400 "Use up to the Quadrupole term"), &
407 description=
"Specify the method to obtain self consistent induced "// &
408 "multipole moments.", &
409 usage=
"POL_SCF CONJUGATE_GRADIENT", &
410 enum_c_vals=
s2a(
"NONE",
"SELF_CONSISTENT",
"CONJUGATE_GRADIENT"), &
411 enum_desc=
s2a(
"No inducible multipoles.", &
412 "Conventional self-consistent iteration.", &
413 "Linear conjugate-gradient optimization of the sum "// &
414 "of the electrostatic and induction energy. This "// &
415 "method does not support non-linear polarization "// &
416 "but is sometimes faster."), &
422 CALL keyword_create(keyword, __location__, name=
"MAX_IPOL_ITER", &
423 description=
"Specify the maximum number of iterations for induced "// &
425 usage=
"MAX_IPOL_ITER {int}", type_of_var=
integer_t, &
426 n_var=1, default_i_val=0)
431 description=
"Specify the rmsd threshold for the derivatives "// &
432 "of the energy towards the Cartesian dipoles components", &
433 usage=
"EPS_POL {real}", type_of_var=
real_t, &
434 n_var=1, default_r_val=0.5e-07_dp)
443 description=
"Controls printing of Ewald properties", &
444 n_keywords=0, n_subsections=1, repeats=.false.)
447 description=
"controls the printing of ewald setup", &
467 cpassert(.NOT.
ASSOCIATED(section))
469 description=
"controls the interpolation for the G-space term", &
470 n_keywords=5, n_subsections=0, repeats=.false.)
472 NULLIFY (keyword, print_key)
475 description=
"the approximate inverse to use to get the starting point"// &
476 " for the linear solver of the spline3 methods", &
477 usage=
"aint_precond copy", &
479 enum_c_vals=
s2a(
"copy",
"spl3_nopbc_aint1",
"spl3_nopbc_precond1", &
480 "spl3_nopbc_aint2",
"spl3_nopbc_precond2",
"spl3_nopbc_precond3"), &
487 description=
"The preconditioner used"// &
488 " for the linear solver of the spline3 methods", &
489 usage=
"precond copy", &
491 enum_c_vals=
s2a(
"copy",
"spl3_nopbc_aint1",
"spl3_nopbc_precond1", &
492 "spl3_nopbc_aint2",
"spl3_nopbc_precond2",
"spl3_nopbc_precond3"), &
499 description=
"accuracy on the solution for spline3 the interpolators", &
500 usage=
"eps_x 1.e-15", default_r_val=1.e-10_dp)
505 description=
"accuracy on the residual for spline3 the interpolators", &
506 usage=
"eps_r 1.e-15", default_r_val=1.e-10_dp)
511 variants=[
'maxiter'], &
512 description=
"the maximum number of iterations", &
513 usage=
"max_iter 200", default_i_val=100)
519 description=
"if convergence information about the linear solver"// &
520 " of the spline methods should be printed", &
522 each_iter_values=[10], filename=
"__STD_OUT__", &
536 SUBROUTINE create_wavelet_section(section)
541 cpassert(.NOT.
ASSOCIATED(section))
543 section, __location__, name=
"wavelet", &
544 description=
"Sets up parameters of wavelet based poisson solver.", &
545 n_keywords=1, n_subsections=0, repeats=.false., &
551 keyword, __location__, name=
"SCF_TYPE", &
552 description=
"Type of scaling function used in the wavelet approach, the total energy depends on this choice, "// &
553 "and the convergence with respect to cutoff depends on the selected scaling functions. "// &
554 "Possible values are 8,14,16,20,24,30,40,50,60,100", &
555 usage=
"SCF_TYPE integer", &
556 n_var=1, default_i_val=40)
560 END SUBROUTINE create_wavelet_section
567 SUBROUTINE create_implicit_ps_section(section)
573 cpassert(.NOT.
ASSOCIATED(section))
575 description=
"Parameters for the implicit (generalized) Poisson solver.", &
577 n_keywords=6, n_subsections=2, repeats=.false.)
579 NULLIFY (subsection, keyword)
581 CALL create_dielectric_section(subsection)
585 CALL create_dbc_section(subsection)
590 keyword, __location__, name=
"BOUNDARY_CONDITIONS", &
591 enum_c_vals=
s2a(
'PERIODIC',
'MIXED',
'MIXED_PERIODIC',
'NEUMANN'), &
592 enum_desc=
s2a(
'periodic boundary conditions',
'Dirichlet + homogeneous Neumann boundary conditions', &
593 'Dirichlet + periodic boundary conditions',
'homogeneous Neumann BC (zero-average solution)'), &
595 description=
"Specifies the type of boundary conditions. Dirichlet=fixed value, Neumann=zero normal deriv. "// &
596 "Mixed and Neumann boundaries essentially requires FFTW3 so that all grid sizes are FFT-able.", &
597 usage=
"BOUNDARY_CONDITIONS <bc_type>", default_i_val=
periodic_bc)
601 CALL keyword_create(keyword, __location__, name=
"ZERO_INITIAL_GUESS", &
602 description=
"Whether or not to use zero potential as initial guess.", &
603 usage=
"ZERO_INITIAL_GUESS <logical>", default_l_val=.false., lone_keyword_l_val=.true.)
608 description=
"Maximum number of iterations.", &
609 usage=
"max_iter <integer>", default_i_val=30)
614 description=
"Stopping tolerance.", &
615 usage=
"tol <real>", default_r_val=1.0e-8_dp)
619 CALL keyword_create(keyword, __location__, name=
"OR_PARAMETER", variants=
s2a(
'omega'), &
620 description=
"Over-relaxation parameter (large epsilon requires smaller omega ~0.1).", &
621 usage=
"OR_PARAMETER <real>", default_r_val=1.0_dp)
626 keyword, __location__, name=
"NEUMANN_DIRECTIONS", &
627 enum_c_vals=
s2a(
'XYZ',
'XY',
'XZ',
'YZ',
'X',
'Y',
'Z'), &
629 description=
"Directions in which homogeneous Neumann conditions are imposed. In the remaining directions "// &
630 "periodic conditions will be enforced. Having specified MIXED or NEUMANN as BOUNDARY_CONDITIONS, "// &
631 "the keyword is meant to be used to combine periodic and homogeneous Neumann conditions at the "// &
632 "boundaries of the simulation cell.", &
633 usage=
"NEUMANN_DIRECTIONS <direction>", default_i_val=
neumannxyz)
637 END SUBROUTINE create_implicit_ps_section
646 SUBROUTINE create_dielectric_section(section)
652 cpassert(.NOT.
ASSOCIATED(section))
654 description=
"Parameters for the dielectric constant function.", &
655 n_keywords=6, n_subsections=2, repeats=.false.)
657 NULLIFY (keyword, subsection)
659 CALL keyword_create(keyword, __location__, name=
"DIELECTRIC_CORE_CORRECTION", &
660 description=
"Avoid spurious values of the dielectric constant at the ionic core for pseudopotentials "// &
661 "where the electron density goes to zero at the core (e.g. GTH). "// &
662 "The correction is based on rho_core.", &
663 usage=
"DIELECTRIC_CORE_CORRECTION <logical>", default_l_val=.true., lone_keyword_l_val=.true.)
668 keyword, __location__, name=
"DIELECTRIC_FUNCTION_TYPE", &
669 enum_c_vals=
s2a(
'density_dependent',
'spatially_dependent',
'spatially_rho_dependent'), &
671 enum_desc=
s2a(
"Dielectric constant as a function of the electron density "// &
672 "as e.g. proposed within the SCCS model.", &
673 "Various regions with different dielectric constants.", &
674 "Various regions with different dielectric constants. The dielectric constant decays to 1.0, "// &
675 "wherever the electron density is present."), &
676 description=
"Preferred type for the dielectric constant function.", &
677 usage=
"DIELECTRIC_FUNCTION_TYPE <method>", default_i_val=
rho_dependent)
681 CALL keyword_create(keyword, __location__, name=
"dielectric_constant", variants=
s2a(
'epsilon'), &
682 description=
"Dielectric constant in the bulk of the solvent.", &
683 usage=
"dielectric_constant <real>", default_r_val=80.0_dp)
688 description=
"Lower density threshold.", &
689 usage=
"rho_min <real>", default_r_val=1.0e-4_dp)
694 description=
"Upper density threshold.", &
695 usage=
"rho_max <real>", default_r_val=1.0e-3_dp)
700 keyword, __location__, name=
"DERIVATIVE_METHOD", &
701 enum_c_vals=
s2a(
'fft',
'fft_use_deps',
'fft_use_drho',
'cd3',
'cd5',
'cd7'), &
704 enum_desc=
s2a(
"FFT based deriv of epsilon, without correction (high cutoff needed).", &
705 "FFT based deriv of epsilon, with correction using gradient of epsilon (high cutoff needed).", &
706 "FFT based deriv of epsilon, with correction using gradient of rho (high cutoff needed).", &
707 "3-point central difference derivative.", &
708 "5-point central difference derivative.", &
709 "7-point central difference derivative (recommended)."), &
710 description=
"Preferred method for evaluating the gradient of ln(eps).", &
715 CALL create_dielec_aa_cuboidal_section(subsection)
719 CALL create_dielec_xaa_annular_section(subsection)
723 END SUBROUTINE create_dielectric_section
730 SUBROUTINE create_dielec_aa_cuboidal_section(section)
735 cpassert(.NOT.
ASSOCIATED(section))
736 CALL section_create(section, __location__, name=
"DIELEC_AA_CUBOIDAL", &
737 description=
"Parameters for creating axis-aligned cuboidal dielectric region. "// &
738 "Note that once such a region is defined, the 'background' dielectric constant "// &
739 "would be the default (80.0), unless a different value is specified using the "// &
740 "keyword IMPLICIT%DIELECTRIC%DIELECTRIC_CONSTANT.", &
741 n_keywords=5, n_subsections=0, repeats=.true.)
745 CALL keyword_create(keyword, __location__, name=
"dielectric_constant", variants=
s2a(
'epsilon'), &
746 description=
"value of the dielectric constant inside the region.", &
747 usage=
"dielectric_constant <real>", default_r_val=80.0_dp)
752 description=
"The X extents of the cuboid.", &
753 usage=
"X_xtnt <xmin(real)> <xmax(real)>", unit_str=
"angstrom", &
754 n_var=2, type_of_var=
real_t)
759 description=
"The Y extents of the cuboid.", &
760 usage=
"Y_xtnt <ymin(real)> <ymax(real)>", unit_str=
"angstrom", &
761 n_var=2, type_of_var=
real_t)
766 description=
"The Z extents of the cuboid.", &
767 usage=
"Z_xtnt <zmin(real)> <zmax(real)>", unit_str=
"angstrom", &
768 n_var=2, type_of_var=
real_t)
772 CALL keyword_create(keyword, __location__, name=
"SMOOTHING_WIDTH", variants=
s2a(
'zeta'), &
773 description=
"The width of the standard mollifier.", &
774 usage=
"SMOOTHING_WIDTH <real>", unit_str=
"angstrom", &
779 END SUBROUTINE create_dielec_aa_cuboidal_section
786 SUBROUTINE create_dielec_xaa_annular_section(section)
791 cpassert(.NOT.
ASSOCIATED(section))
792 CALL section_create(section, __location__, name=
"DIELEC_XAA_ANNULAR", &
793 description=
"Parameters for creating x-axis-aligned annular dielectric region. "// &
794 "Note that once such a region is defined, the 'background' dielectric constant "// &
795 "would be the default (80.0), unless a different value is specified using the "// &
796 "keyword IMPLICIT%DIELECTRIC%DIELECTRIC_CONSTANT.", &
797 n_keywords=5, n_subsections=0, repeats=.true.)
801 CALL keyword_create(keyword, __location__, name=
"dielectric_constant", variants=
s2a(
'epsilon'), &
802 description=
"value of the dielectric constant inside the region.", &
803 usage=
"dielectric_constant <real>", default_r_val=80.0_dp)
808 description=
"The X extents of the annulus.", &
809 usage=
"X_xtnt <xmin(real)> <xmax(real)>", unit_str=
"angstrom", &
810 n_var=2, type_of_var=
real_t)
815 description=
"The y and z coordinates of the annulus' base center.", &
816 usage=
"base_center <y(real)> <z(real)>", unit_str=
"angstrom", &
817 n_var=2, type_of_var=
real_t)
822 description=
"The base radius of the annulus.", &
823 usage=
"base_radii <r1(real)> <r2(real)>", unit_str=
"angstrom", &
824 n_var=2, type_of_var=
real_t)
828 CALL keyword_create(keyword, __location__, name=
"smoothing_width", variants=
s2a(
'zeta'), &
829 description=
"The width of the standard mollifier.", &
830 usage=
"smoothing_width <real>", unit_str=
"angstrom", &
835 END SUBROUTINE create_dielec_xaa_annular_section
842 SUBROUTINE create_dbc_section(section)
848 cpassert(.NOT.
ASSOCIATED(section))
850 description=
"Parameters for creating Dirichlet type boundary conditions.", &
851 n_keywords=1, n_subsections=4, repeats=.false.)
855 CALL keyword_create(keyword, __location__, name=
"VERBOSE_OUTPUT", &
856 description=
"Print out the coordinates of the vertices defining Dirichlet regions and their "// &
857 "tessellations (in Angstrom), the values of the electrostatic potential at the regions (in a.u.), "// &
858 "and their corresponding evaluated Lagrange multipliers.", &
859 usage=
"VERBOSE_OUTPUT <logical>", default_l_val=.false., lone_keyword_l_val=.true.)
865 CALL create_aa_planar_section(subsection)
869 CALL create_planar_section(subsection)
873 CALL create_aa_cylindrical_section(subsection)
877 CALL create_aa_cuboidal_section(subsection)
881 END SUBROUTINE create_dbc_section
888 SUBROUTINE create_aa_planar_section(section)
893 cpassert(.NOT.
ASSOCIATED(section))
895 description=
"Parameters for creating axis-aligned planar (rectangular) Dirichlet boundary regions.", &
896 n_keywords=10, n_subsections=0, repeats=.true.)
901 description=
"The value of the fixed potential to be imposed at the axis-aligned Dirichlet boundary.", &
902 usage=
"v_D <real>", unit_str=
"volt", type_of_var=
real_t)
906 CALL keyword_create(keyword, __location__, name=
"OSCILLATING_FRACTION", &
907 description=
"A fraction of the field can be set to oscilate over time.", &
908 usage=
"OSCILLATING_FRACTION <real>", default_r_val=0.0_dp, type_of_var=
real_t)
913 description=
"The frequency with which the oscillating fraction oscillates.", &
914 usage=
"FREQUENCY <real>", default_r_val=0.0_dp, unit_str=
"s^-1", type_of_var=
real_t)
919 description=
"The phase of the oscillattion. A phase of zero corresponds to a cosine function. ", &
920 usage=
"PHASE <real>", default_r_val=0.0_dp, type_of_var=
real_t)
924 CALL keyword_create(keyword, __location__, name=
"PARALLEL_PLANE", &
925 enum_c_vals=
s2a(
'XY',
'YZ',
'XZ'), &
927 description=
"The coordinate plane that the region is parallel to.", &
928 usage=
"PARALLEL_PLANE <plane>", &
934 description=
"The intercept of the rectangle's plane.", &
935 usage=
"INTERCEPT <real>", unit_str=
"angstrom", &
941 description=
"The X extents of the rectangle.", &
942 usage=
"X_xtnt <xmin(real)> <xmax(real)>", unit_str=
"angstrom", &
943 n_var=2, type_of_var=
real_t)
948 description=
"The Y extents of the rectangle.", &
949 usage=
"Y_xtnt <ymin(real)> <ymax(real)>", unit_str=
"angstrom", &
950 n_var=2, type_of_var=
real_t)
955 description=
"The Z extents of the rectangle.", &
956 usage=
"Z_xtnt <zmin(real)> <zmax(real)>", unit_str=
"angstrom", &
957 n_var=2, type_of_var=
real_t)
962 description=
"The number of partitions in the directions of the unit vectors generating the "// &
963 "corresponding PARALLEL_PLANE (e1, e2 or e3) for tiling the rectangluar region.", &
964 usage=
"N_PRTN <integer> <integer>", &
965 n_var=2, default_i_vals=[1, 1])
970 description=
"The thickness of the planar Dirichlet region.", &
971 usage=
"THICKNESS <real>", unit_str=
"angstrom", &
976 CALL keyword_create(keyword, __location__, name=
"SMOOTHING_WIDTH", variants=
s2a(
'SIGMA'), &
977 description=
"The width of the transition region between the Dirichlet subdomain and its surrounding.", &
978 usage=
"SMOOTHING_WIDTH <real>", unit_str=
"angstrom", &
983 CALL keyword_create(keyword, __location__, name=
"PERIODIC_REGION", &
984 description=
"Whether or not to take into consideration the effects of the periodicity of the "// &
985 "simulation cell (MIXED_PERIODIC bc) for regions defined sufficiently close to the boundaries.", &
986 usage=
"PERIODIC_REGION <logical>", default_l_val=.false., lone_keyword_l_val=.true.)
990 END SUBROUTINE create_aa_planar_section
997 SUBROUTINE create_planar_section(section)
1002 cpassert(.NOT.
ASSOCIATED(section))
1004 section, __location__, name=
"PLANAR", &
1005 description=
"Parameters for creating planar (rectangular) Dirichlet boundary regions with given vertices.", &
1006 n_keywords=7, n_subsections=0, repeats=.true.)
1011 description=
"The value of the fixed potential to be imposed at the planar Dirichlet boundary.", &
1012 usage=
"v_D <real>", unit_str=
"volt", type_of_var=
real_t)
1016 CALL keyword_create(keyword, __location__, name=
"OSCILLATING_FRACTION", &
1017 description=
"A fraction of the field can be set to oscilate over time.", &
1018 usage=
"OSCILLATING_FRACTION <real>", default_r_val=0.0_dp, type_of_var=
real_t)
1023 description=
"The frequency with which the oscillating fraction oscillates.", &
1024 usage=
"FREQUENCY <real>", default_r_val=0.0_dp, unit_str=
"s^-1", type_of_var=
real_t)
1029 description=
"The phase of the oscillattion. A phase of zero corresponds to a cosine function. ", &
1030 usage=
"PHASE <real>", default_r_val=0.0_dp, type_of_var=
real_t)
1035 description=
"Coordinates of the vertex A.", &
1036 usage=
"A <x(real)> <y(real)> <z(real)>", unit_str=
"angstrom", &
1037 n_var=3, type_of_var=
real_t)
1042 description=
"Coordinates of the vertex B.", &
1043 usage=
"B <x(real)> <y(real)> <z(real)>", unit_str=
"angstrom", &
1044 n_var=3, type_of_var=
real_t)
1049 description=
"Coordinates of the vertex C.", &
1050 usage=
"C <x(real)> <y(real)> <z(real)>", unit_str=
"angstrom", &
1051 n_var=3, type_of_var=
real_t)
1056 keyword, __location__, name=
"N_PRTN", &
1057 description=
"The number of partitions along the edges for tiling the rectangular region. If the edges "// &
1058 "have different lengths, from the two given values, the larger one will be assigned to the longer edge.", &
1059 usage=
"N_PRTN <integer> <integer>", &
1060 n_var=2, default_i_vals=[1, 1])
1065 description=
"The thickness of the planar Dirichlet region.", &
1066 usage=
"THICKNESS <real>", unit_str=
"angstrom", &
1071 CALL keyword_create(keyword, __location__, name=
"SMOOTHING_WIDTH", variants=
s2a(
'SIGMA'), &
1072 description=
"The width of the transition region between the Dirichlet subdomain and its surrounding.", &
1073 usage=
"SMOOTHING_WIDTH <real>", unit_str=
"angstrom", &
1078 END SUBROUTINE create_planar_section
1085 SUBROUTINE create_aa_cylindrical_section(section)
1090 cpassert(.NOT.
ASSOCIATED(section))
1091 CALL section_create(section, __location__, name=
"AA_CYLINDRICAL", &
1092 description=
"Parameters for creating axis-aligned cylindrical Dirichlet boundary regions.", &
1093 n_keywords=11, n_subsections=0, repeats=.true.)
1098 description=
"The value of the fixed potential to be imposed at the cylindrical Dirichlet boundary.", &
1099 usage=
"v_D <real>", unit_str=
"volt", type_of_var=
real_t)
1103 CALL keyword_create(keyword, __location__, name=
"OSCILLATING_FRACTION", &
1104 description=
"A fraction of the field can be set to oscilate over time.", &
1105 usage=
"OSCILLATING_FRACTION <real>", default_r_val=0.0_dp, type_of_var=
real_t)
1110 description=
"The frequency with which the oscillating fraction oscillates.", &
1111 usage=
"FREQUENCY <real>", default_r_val=0.0_dp, unit_str=
"s^-1", type_of_var=
real_t)
1116 description=
"The phase of the oscillattion. A phase of zero corresponds to a cosine function. ", &
1117 usage=
"PHASE <real>", default_r_val=0.0_dp, type_of_var=
real_t)
1121 CALL keyword_create(keyword, __location__, name=
"PARALLEL_AXIS", &
1122 enum_c_vals=
s2a(
'X',
'Y',
'Z'), &
1124 description=
"The coordinate axis that the cylindrical region extends along.", &
1125 usage=
"PARALLEL_AXIS <axis>", &
1131 description=
"The extents of the cylinder along its central axis.", &
1132 usage=
"xtnt <xmin(real)> <xmax(real)>", unit_str=
"angstrom", &
1133 n_var=2, type_of_var=
real_t)
1138 description=
"The y and z coordinates (x and z or x and y coordinates, "// &
1139 "depending on the choice of the parallel axis) of the cylinder's base center.", &
1140 usage=
"BASE_CENTER <y(real)> <z(real)>", unit_str=
"angstrom", &
1141 n_var=2, type_of_var=
real_t)
1146 description=
"The base radius of the cylinder.", &
1147 usage=
"BASE_RADIUS <real>", unit_str=
"angstrom", &
1153 description=
"The number of sides (faces) of the n-gonal prism approximating the cylinder.", &
1154 usage=
"N_SIDES <integer>", default_i_val=5)
1159 enum_c_vals=
s2a(
'CIRCUMSCRIBED',
'INSCRIBED'), &
1161 description=
"Specifies the type of the n-gonal prism approximating the cylinder.", &
1167 keyword, __location__, name=
"N_PRTN", &
1168 description=
"The number of partitions along the face edges of the prism for tiling. If the edges "// &
1169 "have different lengths, from the two given values, the larger one will be assigned to the longer edge.", &
1170 usage=
"N_PRTN <integer> <integer>", &
1171 n_var=2, default_i_vals=[1, 1])
1176 description=
"The thickness of the cylinder.", &
1177 usage=
"THICKNESS <real>", unit_str=
"angstrom", &
1182 CALL keyword_create(keyword, __location__, name=
"SMOOTHING_WIDTH", variants=
s2a(
'SIGMA'), &
1183 description=
"The width of the transition region between the Dirichlet subdomain and its surrounding.", &
1184 usage=
"SMOOTHING_WIDTH <real>", unit_str=
"angstrom", &
1190 keyword, __location__, name=
"delta_alpha", &
1191 description=
"A central angle specifying the gap between the faces of the n-gonal prism. To avoide overlap "// &
1192 "between the cuboids (of the given thickness) built on top of the faces, a larger value is required if the"// &
1193 " number of faces (N_SIDES) is quite few and/or the base radius is fairly small.", &
1194 usage=
"delta_alpha <real>", default_r_val=0.05_dp, unit_str=
"rad")
1198 END SUBROUTINE create_aa_cylindrical_section
1205 SUBROUTINE create_aa_cuboidal_section(section)
1210 cpassert(.NOT.
ASSOCIATED(section))
1212 description=
"Parameters for creating axis-aligned cuboidal Dirichlet regions.", &
1213 n_keywords=7, n_subsections=0, repeats=.true.)
1218 description=
"The value of the fixed potential to be imposed at the region.", &
1219 usage=
"v_D <real>", unit_str=
"volt", type_of_var=
real_t)
1223 CALL keyword_create(keyword, __location__, name=
"OSCILLATING_FRACTION", &
1224 description=
"A fraction of the field can be set to oscilate over time.", &
1225 usage=
"OSCILLATING_FRACTION <real>", default_r_val=0.0_dp, type_of_var=
real_t)
1230 description=
"The frequency with which the oscillating fraction oscillates.", &
1231 usage=
"FREQUENCY <real>", default_r_val=0.0_dp, unit_str=
"s^-1", type_of_var=
real_t)
1236 description=
"The phase of the oscillattion. A phase of zero corresponds to a cosine function. ", &
1237 usage=
"PHASE <real>", default_r_val=0.0_dp, type_of_var=
real_t)
1242 description=
"The X extents of the cuboid.", &
1243 usage=
"X_xtnt <xmin(real)> <xmax(real)>", unit_str=
"angstrom", &
1244 n_var=2, type_of_var=
real_t)
1249 description=
"The Y extents of the cuboid.", &
1250 usage=
"Y_xtnt <ymin(real)> <ymax(real)>", unit_str=
"angstrom", &
1251 n_var=2, type_of_var=
real_t)
1256 description=
"The Z extents of the cuboid.", &
1257 usage=
"Z_xtnt <zmin(real)> <zmax(real)>", unit_str=
"angstrom", &
1258 n_var=2, type_of_var=
real_t)
1263 description=
"The number of partitions in the x, y and z directions for partitioning the cuboid.", &
1264 usage=
"N_PRTN <integer> <integer> <integer>", &
1265 n_var=3, default_i_vals=[1, 1, 1])
1269 CALL keyword_create(keyword, __location__, name=
"SMOOTHING_WIDTH", variants=
s2a(
'SIGMA'), &
1270 description=
"The width of the transition region between the Dirichlet subdomain and its surrounding.", &
1271 usage=
"SMOOTHING_WIDTH <real>", unit_str=
"angstrom", &
1276 CALL keyword_create(keyword, __location__, name=
"PERIODIC_REGION", &
1277 description=
"Whether or not to take into consideration the effects of the periodicity of the "// &
1278 "simulation cell (MIXED_PERIODIC bc) for regions defined sufficiently close to the boundaries.", &
1279 usage=
"PERIODIC_REGION <logical>", default_l_val=.false., lone_keyword_l_val=.true.)
1283 END SUBROUTINE create_aa_cuboidal_section
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public blochl1995
integer, save, public toukmaji1996
integer, save, public aguado2003
integer, save, public genovese2006
integer, save, public laino2008
integer, save, public genovese2007
integer, save, public essmann1995
integer, save, public darden1993
integer, save, public ewald1921
integer, save, public martyna1999
integer, save, public banihashemian2016
Handles all functions related to the CELL.
integer, parameter, public use_perd_xyz
integer, parameter, public use_perd_y
integer, parameter, public use_perd_xz
integer, parameter, public use_perd_x
integer, parameter, public use_perd_z
integer, parameter, public use_perd_yz
integer, parameter, public use_perd_none
integer, parameter, public use_perd_xy
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 medium_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
real(kind=dp) function, public cp_unit_to_cp2k(value, unit_str, defaults, power)
converts to the internal cp2k units to the given unit
the type I Discrete Cosine Transform (DCT-I)
integer, parameter, public neumannx
integer, parameter, public neumannxy
integer, parameter, public neumannxz
integer, parameter, public neumannxyz
integer, parameter, public neumannz
integer, parameter, public neumannyz
integer, parameter, public neumanny
dielectric constant data type
integer, parameter, public derivative_fft_use_drho
integer, parameter, public derivative_fft_use_deps
integer, parameter, public derivative_fft
integer, parameter, public derivative_cd5
integer, parameter, public spatially_rho_dependent
integer, parameter, public derivative_cd3
integer, parameter, public spatially_dependent
integer, parameter, public rho_dependent
integer, parameter, public derivative_cd7
Dirichlet boundary condition data types.
integer, parameter, public z_axis
integer, parameter, public inscribed
integer, parameter, public y_axis
integer, parameter, public circumscribed
integer, parameter, public x_axis
integer, parameter, public xz_plane
integer, parameter, public yz_plane
integer, parameter, public xy_plane
Defines the basic variable types.
integer, parameter, public dp
Multipole structure: for multipole (fixed and induced) in FF based MD.
integer, parameter, public do_multipole_quadrupole
integer, parameter, public do_multipole_dipole
integer, parameter, public do_multipole_charge
integer, parameter, public do_multipole_none
Types containing essential information for running implicit (iterative) Poisson solver.
integer, parameter, public neumann_bc
integer, parameter, public mixed_bc
integer, parameter, public mixed_periodic_bc
integer, parameter, public periodic_bc
functions related to the poisson solver on regular grids
integer, parameter, public pw_poisson_wavelet
integer, parameter, public pw_poisson_periodic
integer, parameter, public pw_poisson_mt
integer, parameter, public do_ewald_pme
integer, parameter, public do_ewald_ewald
integer, parameter, public pw_poisson_implicit
integer, parameter, public do_ewald_none
integer, parameter, public pw_poisson_analytic
integer, parameter, public do_ewald_spme
integer, parameter, public pw_poisson_multipole
different utils that are useful to manipulate splines on the regular grid of a pw
integer, parameter, public precond_spl3_3
integer, parameter, public precond_spl3_aint
integer, parameter, public no_precond
integer, parameter, public precond_spl3_2
integer, parameter, public precond_spl3_aint2
integer, parameter, public precond_spl3_1
Utilities for string manipulations.