(git:34ef472)
input_cp2k_poisson.F
Go to the documentation of this file.
1 !--------------------------------------------------------------------------------------------------!
2 ! CP2K: A general program to perform molecular dynamics simulations !
3 ! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4 ! !
5 ! SPDX-License-Identifier: GPL-2.0-or-later !
6 !--------------------------------------------------------------------------------------------------!
7 
8 ! **************************************************************************************************
9 !> \brief function that build the poisson section of the input
10 !> \par History
11 !> 03.2006 fusing of poisson_dft and poisson_mm
12 !> \author fawzi
13 ! **************************************************************************************************
15  USE bibliography, ONLY: &
18  USE cell_types, ONLY: use_perd_none,&
19  use_perd_x,&
20  use_perd_xy,&
21  use_perd_xyz,&
22  use_perd_xz,&
23  use_perd_y,&
24  use_perd_yz,&
30  USE cp_units, ONLY: cp_unit_to_cp2k
31  USE dct, ONLY: neumannx,&
32  neumannxy,&
33  neumannxyz,&
34  neumannxz,&
35  neumanny,&
36  neumannyz,&
37  neumannz
38  USE dielectric_types, ONLY: &
42  inscribed,&
43  x_axis,&
44  xy_plane,&
45  xz_plane,&
46  y_axis,&
47  yz_plane,&
48  z_axis
49  USE input_constants, ONLY: do_fist_pol_cg,&
55  keyword_type
60  section_type
61  USE input_val_types, ONLY: enum_t,&
62  integer_t,&
63  real_t
64  USE kinds, ONLY: dp
69  USE ps_implicit_types, ONLY: mixed_bc,&
71  neumann_bc,&
73  USE pw_poisson_types, ONLY: &
77  USE pw_spline_utils, ONLY: no_precond,&
83  USE string_utilities, ONLY: s2a
84 #include "./base/base_uses.f90"
85 
86  IMPLICIT NONE
87  PRIVATE
88 
89  LOGICAL, PRIVATE, PARAMETER :: debug_this_module = .true.
90  CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'input_cp2k_poisson'
91 
92  PUBLIC :: create_poisson_section, &
95 !***
96 CONTAINS
97 
98 ! **************************************************************************************************
99 !> \brief Creates the Poisson section
100 !> \param section the section to create
101 !> \author teo
102 ! **************************************************************************************************
103  SUBROUTINE create_poisson_section(section)
104  TYPE(section_type), POINTER :: section
105 
106  TYPE(keyword_type), POINTER :: keyword
107  TYPE(section_type), POINTER :: subsection
108 
109  cpassert(.NOT. ASSOCIATED(section))
110  CALL section_create(section, __location__, name="poisson", &
111  description="Sets up the poisson resolutor.", &
112  n_keywords=1, n_subsections=0, repeats=.false.)
113 
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."), &
134  citations=(/blochl1995, martyna1999, genovese2006, genovese2007/), &
135  default_i_val=pw_poisson_periodic)
136  CALL section_add_keyword(section, keyword)
137  CALL keyword_release(keyword)
138 
139  CALL keyword_create(keyword, __location__, name="PERIODIC", &
140  description="Specify the directions in which PBC apply. Important notice,"// &
141  " this only applies to the electrostatics."// &
142  " See the CELL section to specify the periodicity used for e.g. the pair lists."// &
143  " Typically the settings should be the same.", &
144  usage="PERIODIC (x|y|z|xy|xz|yz|xyz|none)", &
145  enum_c_vals=s2a("x", "y", "z", "xy", "xz", "yz", "xyz", "none"), &
146  enum_i_vals=(/use_perd_x, use_perd_y, use_perd_z, &
149  default_i_val=use_perd_xyz)
150  CALL section_add_keyword(section, keyword)
151  CALL keyword_release(keyword)
152 
153  CALL create_mt_section(subsection)
154  CALL section_add_subsection(section, subsection)
155  CALL section_release(subsection)
156 
157  CALL create_wavelet_section(subsection)
158  CALL section_add_subsection(section, subsection)
159  CALL section_release(subsection)
160 
161  CALL create_multipole_section(subsection)
162  CALL section_add_subsection(section, subsection)
163  CALL section_release(subsection)
164 
165  CALL create_ewald_section(subsection)
166  CALL section_add_subsection(section, subsection)
167  CALL section_release(subsection)
168 
169  CALL create_implicit_ps_section(subsection)
170  CALL section_add_subsection(section, subsection)
171  CALL section_release(subsection)
172  END SUBROUTINE create_poisson_section
173 
174 ! **************************************************************************************************
175 !> \brief Section to set-up parameters for decoupling using the Bloechl scheme
176 !> \param section the section to create
177 !> \author teo
178 ! **************************************************************************************************
179  SUBROUTINE create_multipole_section(section)
180  TYPE(section_type), POINTER :: section
181 
182  TYPE(keyword_type), POINTER :: keyword
183  TYPE(section_type), POINTER :: subsection
184 
185  cpassert(.NOT. ASSOCIATED(section))
186 
187  CALL section_create(section, __location__, name="MULTIPOLE", &
188  description="This section is used to set up the decoupling of QM periodic images with "// &
189  "the use of density derived atomic point charges.", &
190  n_keywords=1, n_subsections=0, repeats=.false.)
191 
192  NULLIFY (keyword, subsection)
193  CALL keyword_create(keyword, __location__, name="RCUT", &
194  description="Real space cutoff for the Ewald sum.", &
195  usage="RCUT {real}", n_var=1, type_of_var=real_t, &
196  unit_str="angstrom")
197  CALL section_add_keyword(section, keyword)
198  CALL keyword_release(keyword)
199 
200  CALL keyword_create(keyword, __location__, name="EWALD_PRECISION", &
201  description="Precision achieved in the Ewald sum.", &
202  usage="EWALD_PRECISION {real}", n_var=1, type_of_var=real_t, &
203  unit_str="hartree", default_r_val=1.0e-6_dp)
204  CALL section_add_keyword(section, keyword)
205  CALL keyword_release(keyword)
206 
207  CALL keyword_create(keyword, __location__, name="ANALYTICAL_GTERM", &
208  description="Evaluates the Gterm in the Ewald Scheme analytically instead of using Splines.", &
209  usage="ANALYTICAL_GTERM <LOGICAL>", &
210  default_l_val=.false., lone_keyword_l_val=.true.)
211  CALL section_add_keyword(section, keyword)
212  CALL keyword_release(keyword)
213 
214  CALL keyword_create(keyword, __location__, name="NGRIDS", &
215  description="Specifies the number of grid points used for the Interpolation of the G-space term", &
216  usage="NGRIDS <integer> <integer> <integer> ", n_var=3, default_i_vals=(/50, 50, 50/))
217  CALL section_add_keyword(section, keyword)
218  CALL keyword_release(keyword)
219 
220  CALL create_gspace_interp_section(subsection)
221  CALL section_add_subsection(section, subsection)
222  CALL section_release(subsection)
223 
224  CALL cp_print_key_section_create(subsection, __location__, "check_spline", &
225  description="Controls the checking of the G-space term Spline Interpolation.", &
226  print_level=medium_print_level, filename="GSpace-SplInterp")
227  CALL section_add_subsection(section, subsection)
228  CALL section_release(subsection)
229 
230  CALL cp_print_key_section_create(subsection, __location__, "program_run_info", &
231  description="Controls the printing of basic information during the run", &
232  print_level=low_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
233  CALL section_add_subsection(section, subsection)
234  CALL section_release(subsection)
235 
236  END SUBROUTINE create_multipole_section
237 
238 ! **************************************************************************************************
239 !> \brief Creates the Martyna-Tuckerman section
240 !> \param section the section to create
241 !> \author teo
242 ! **************************************************************************************************
243  SUBROUTINE create_mt_section(section)
244  TYPE(section_type), POINTER :: section
245 
246  TYPE(keyword_type), POINTER :: keyword
247 
248  cpassert(.NOT. ASSOCIATED(section))
249  CALL section_create(section, __location__, name="mt", &
250  description="Sets up parameters of Martyna-Tuckerman poisson solver. "// &
251  "Note that exact results are only guaranteed if the unit cell is "// &
252  "twice as large as charge density (and serious artefacts can result "// &
253  "if the cell is much smaller).", &
254  n_keywords=1, n_subsections=0, repeats=.false., &
255  citations=(/martyna1999/))
256 
257  NULLIFY (keyword)
258 
259  CALL keyword_create(keyword, __location__, name="ALPHA", &
260  description="Convergence parameter ALPHA*RMIN. Default value 7.0", &
261  usage="ALPHA real", &
262  n_var=1, default_r_val=7.0_dp)
263  CALL section_add_keyword(section, keyword)
264  CALL keyword_release(keyword)
265 
266  CALL keyword_create(keyword, __location__, name="REL_CUTOFF", &
267  description="Specify the multiplicative factor for the CUTOFF keyword in MULTI_GRID"// &
268  " section. The result gives the cutoff at which the 1/r non-periodic FFT3D is evaluated."// &
269  " Default is 2.0", &
270  usage="REL_CUTOFF real", &
271  n_var=1, default_r_val=2.0_dp)
272  CALL section_add_keyword(section, keyword)
273  CALL keyword_release(keyword)
274 
275  END SUBROUTINE create_mt_section
276 
277 ! **************************************************************************************************
278 !> \brief ...
279 !> \param section will contain the ewald section
280 !> \author fawzi
281 ! **************************************************************************************************
282  SUBROUTINE create_ewald_section(section)
283  TYPE(section_type), POINTER :: section
284 
285  TYPE(keyword_type), POINTER :: keyword
286  TYPE(section_type), POINTER :: print_key, subsection
287 
288  cpassert(.NOT. ASSOCIATED(section))
289  CALL section_create(section, __location__, name="ewald", &
290  description="Ewald parameters controlling electrostatic only for CLASSICAL MM.", &
291  n_keywords=7, n_subsections=0, repeats=.false., &
293 
294  NULLIFY (keyword, print_key, subsection)
295  CALL keyword_create( &
296  keyword, __location__, name="EWALD_TYPE", &
297  description="The type of ewald you want to perform.", &
298  citations=(/ewald1921, essmann1995, darden1993/), &
299  usage="EWALD_TYPE (NONE|EWALD|PME|SPME)", &
300  default_i_val=do_ewald_ewald, &
301  enum_c_vals=(/"none ", &
302  "ewald ", &
303  "pme ", &
304  "spme "/), &
305  enum_i_vals=(/do_ewald_none, &
306  do_ewald_ewald, &
307  do_ewald_pme, &
308  do_ewald_spme/), &
309  enum_desc=s2a("NONE standard real-space coulomb potential is computed together with the non-bonded contributions", &
310  "EWALD is the standard non-fft based ewald", &
311  "PME is the particle mesh using fft interpolation", &
312  "SPME is the smooth particle mesh using beta-Euler splines (recommended)"))
313  CALL section_add_keyword(section, keyword)
314  CALL keyword_release(keyword)
315 
316  CALL keyword_create(keyword, __location__, name="EWALD_ACCURACY", &
317  description="Expected accuracy in the Ewald sum. This number affects only the calculation of "// &
318  "the cutoff for the real-space term of the ewald summation (EWALD|PME|SPME) as well as the "// &
319  "construction of the neighbor lists (if the cutoff for non-bonded terms is smaller than the "// &
320  "value employed to compute the EWALD real-space term). This keyword has no "// &
321  "effect on the reciprocal space term (which can be tuned independently).", &
322  usage="EWALD_ACCURACY {real}", n_var=1, type_of_var=real_t, &
323  unit_str="hartree", default_r_val=1.0e-6_dp)
324  CALL section_add_keyword(section, keyword)
325  CALL keyword_release(keyword)
326 
327  CALL keyword_create(keyword, __location__, name="RCUT", &
328  description="Explicitly provide the real-space cutoff of the ewald summation (EWALD|PME|SPME). "// &
329  "This value is ignored in Tight-binding applications (rcut from basis overlap is used). "// &
330  "If present, overwrites the estimate of EWALD_ACCURACY and may affect the "// &
331  "construction of the neighbor lists for non-bonded terms (in FIST), if the value "// &
332  "specified is larger than the cutoff for non-bonded interactions.", &
333  usage="RCUT 5.0", n_var=1, type_of_var=real_t, unit_str="angstrom")
334  CALL section_add_keyword(section, keyword)
335  CALL keyword_release(keyword)
336 
337  CALL keyword_create(keyword, __location__, name="alpha", &
338  description="alpha parameter associated with Ewald (EWALD|PME|SPME). "// &
339  "Recommended for small systems is alpha = 3.5 / r_cut. "// &
340  "For Tight-binding application a recommended value is alpha = 1.0. "// &
341  "Tuning alpha, r_cut and gmax is needed to obtain O(N**1.5) scaling for ewald.", &
342  usage="alpha .30", &
343  default_r_val=cp_unit_to_cp2k(value=0.35_dp, unit_str="angstrom^-1"), &
344  unit_str='angstrom^-1')
345  CALL section_add_keyword(section, keyword)
346  CALL keyword_release(keyword)
347 
348  CALL keyword_create(keyword, __location__, name="gmax", &
349  description="number of grid points (SPME and EWALD). If a single number is specified, "// &
350  "the same number of points is used for all three directions on the grid. "// &
351  "If three numbers are given, each direction can have a different number of points. "// &
352  "The number of points needs to be FFTable (which depends on the library used) and odd for EWALD. "// &
353  "The optimal number depends e.g. on alpha and the size of the cell. 1 point per Angstrom is common.", &
354  usage="gmax 25 25 25", n_var=-1, type_of_var=integer_t)
355  CALL section_add_keyword(section, keyword)
356  CALL keyword_release(keyword)
357 
358  CALL keyword_create(keyword, __location__, name="ns_max", &
359  description="number of grid points on small mesh (PME only), should be odd.", &
360  usage="ns_max 11", default_i_val=11)
361  CALL section_add_keyword(section, keyword)
362  CALL keyword_release(keyword)
363 
364  CALL keyword_create(keyword, __location__, name="o_spline", &
365  description="order of the beta-Euler spline (SPME only)", &
366  usage="o_spline 6", default_i_val=6)
367  CALL section_add_keyword(section, keyword)
368  CALL keyword_release(keyword)
369 
370  CALL keyword_create(keyword, __location__, name="epsilon", &
371  description="tolerance of gaussians for fft interpolation (PME only)", &
372  usage="epsilon 1e-6", default_r_val=1.e-6_dp)
373  CALL section_add_keyword(section, keyword)
374  CALL keyword_release(keyword)
375 
376  NULLIFY (subsection)
377  CALL create_rsgrid_section(subsection)
378  CALL section_add_subsection(section, subsection)
379  CALL section_release(subsection)
380 
381  NULLIFY (subsection)
382  CALL section_create(subsection, __location__, name="MULTIPOLES", &
383  description="Enables the use of multipoles in the treatment of the electrostatics.", &
384  n_keywords=0, n_subsections=1, repeats=.false., &
385  citations=(/aguado2003, laino2008/))
386 
387  CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
388  description="Controls the activation of the Multipoles", &
389  usage="&MULTIPOLES T", default_l_val=.false., lone_keyword_l_val=.true.)
390  CALL section_add_keyword(subsection, keyword)
391  CALL keyword_release(keyword)
392 
393  CALL keyword_create(keyword, __location__, name="MAX_MULTIPOLE_EXPANSION", &
394  description="Specify the maximum level of multipoles expansion used "// &
395  "for the electrostatics.", &
396  usage="MAX_MULTIPOLE_EXPANSION DIPOLE", &
397  enum_c_vals=s2a("NONE", "CHARGE", "DIPOLE", "QUADRUPOLE"), &
398  enum_desc=s2a("No multipolar terms! Check the codes providing a zero contribution.", &
399  "Use up to the Charge term", &
400  "Use up to the Dipole term", &
401  "Use up to the Quadrupole term"), &
403  do_multipole_quadrupole/), type_of_var=enum_t)
404  CALL section_add_keyword(subsection, keyword)
405  CALL keyword_release(keyword)
406 
407  CALL keyword_create(keyword, __location__, name="POL_SCF", &
408  description="Specify the method to obtain self consistent induced "// &
409  "multipole moments.", &
410  usage="POL_SCF CONJUGATE_GRADIENT", &
411  enum_c_vals=s2a("NONE", "SELF_CONSISTENT", "CONJUGATE_GRADIENT"), &
412  enum_desc=s2a("No inducible multipoles.", &
413  "Conventional self-consistent iteration.", &
414  "Linear conjugate-gradient optimization of the sum "// &
415  "of the electrostatic and induction energy. This "// &
416  "method does not support non-linear polarization "// &
417  "but is sometimes faster."), &
418  enum_i_vals=(/do_fist_pol_none, do_fist_pol_sc, do_fist_pol_cg/), &
419  type_of_var=enum_t, default_i_val=do_fist_pol_none)
420  CALL section_add_keyword(subsection, keyword)
421  CALL keyword_release(keyword)
422 
423  CALL keyword_create(keyword, __location__, name="MAX_IPOL_ITER", &
424  description="Specify the maximum number of iterations for induced "// &
425  "dipoles", &
426  usage="MAX_IPOL_ITER {int}", type_of_var=integer_t, &
427  n_var=1, default_i_val=0)
428  CALL section_add_keyword(subsection, keyword)
429  CALL keyword_release(keyword)
430 
431  CALL keyword_create(keyword, __location__, name="EPS_POL", &
432  description="Specify the rmsd threshold for the derivatives "// &
433  "of the energy towards the Cartesian dipoles components", &
434  usage="EPS_POL {real}", type_of_var=real_t, &
435  n_var=1, default_r_val=0.5e-07_dp)
436  CALL section_add_keyword(subsection, keyword)
437  CALL keyword_release(keyword)
438 
439  CALL section_add_subsection(section, subsection)
440  CALL section_release(subsection)
441 
442  NULLIFY (subsection)
443  CALL section_create(subsection, __location__, name="PRINT", &
444  description="Controls printing of Ewald properties", &
445  n_keywords=0, n_subsections=1, repeats=.false.)
446  NULLIFY (print_key)
447  CALL cp_print_key_section_create(print_key, __location__, "PROGRAM_RUN_INFO", &
448  description="controls the printing of ewald setup", &
449  print_level=low_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
450  CALL section_add_subsection(subsection, print_key)
451  CALL section_release(print_key)
452  CALL section_add_subsection(section, subsection)
453  CALL section_release(subsection)
454 
455  END SUBROUTINE create_ewald_section
456 
457 ! **************************************************************************************************
458 !> \brief creates the interpolation section for the periodic QM/MM
459 !> \param section ...
460 !> \author tlaino
461 ! **************************************************************************************************
462  SUBROUTINE create_gspace_interp_section(section)
463  TYPE(section_type), POINTER :: section
464 
465  TYPE(keyword_type), POINTER :: keyword
466  TYPE(section_type), POINTER :: print_key
467 
468  cpassert(.NOT. ASSOCIATED(section))
469  CALL section_create(section, __location__, name="interpolator", &
470  description="controls the interpolation for the G-space term", &
471  n_keywords=5, n_subsections=0, repeats=.false.)
472 
473  NULLIFY (keyword, print_key)
474 
475  CALL keyword_create(keyword, __location__, name="aint_precond", &
476  description="the approximate inverse to use to get the starting point"// &
477  " for the linear solver of the spline3 methods", &
478  usage="kind spline3", &
479  default_i_val=precond_spl3_aint, &
480  enum_c_vals=s2a("copy", "spl3_nopbc_aint1", "spl3_nopbc_precond1", &
481  "spl3_nopbc_aint2", "spl3_nopbc_precond2", "spl3_nopbc_precond3"), &
482  enum_i_vals=(/no_precond, precond_spl3_aint, precond_spl3_1, &
484  CALL section_add_keyword(section, keyword)
485  CALL keyword_release(keyword)
486 
487  CALL keyword_create(keyword, __location__, name="precond", &
488  description="The preconditioner used"// &
489  " for the linear solver of the spline3 methods", &
490  usage="kind spline3", &
491  default_i_val=precond_spl3_3, &
492  enum_c_vals=s2a("copy", "spl3_nopbc_aint1", "spl3_nopbc_precond1", &
493  "spl3_nopbc_aint2", "spl3_nopbc_precond2", "spl3_nopbc_precond3"), &
494  enum_i_vals=(/no_precond, precond_spl3_aint, precond_spl3_1, &
496  CALL section_add_keyword(section, keyword)
497  CALL keyword_release(keyword)
498 
499  CALL keyword_create(keyword, __location__, name="eps_x", &
500  description="accuracy on the solution for spline3 the interpolators", &
501  usage="eps_x 1.e-15", default_r_val=1.e-10_dp)
502  CALL section_add_keyword(section, keyword)
503  CALL keyword_release(keyword)
504 
505  CALL keyword_create(keyword, __location__, name="eps_r", &
506  description="accuracy on the residual for spline3 the interpolators", &
507  usage="eps_r 1.e-15", default_r_val=1.e-10_dp)
508  CALL section_add_keyword(section, keyword)
509  CALL keyword_release(keyword)
510 
511  CALL keyword_create(keyword, __location__, name="max_iter", &
512  variants=(/'maxiter'/), &
513  description="the maximum number of iterations", &
514  usage="max_iter 200", default_i_val=100)
515  CALL section_add_keyword(section, keyword)
516  CALL keyword_release(keyword)
517 
518  NULLIFY (print_key)
519  CALL cp_print_key_section_create(print_key, __location__, "conv_info", &
520  description="if convergence information about the linear solver"// &
521  " of the spline methods should be printed", &
522  print_level=medium_print_level, each_iter_names=s2a("SPLINE_FIND_COEFFS"), &
523  each_iter_values=(/10/), filename="__STD_OUT__", &
524  add_last=add_last_numeric)
525  CALL section_add_subsection(section, print_key)
526  CALL section_release(print_key)
527 
528  END SUBROUTINE create_gspace_interp_section
529 
530 ! **************************************************************************************************
531 !> \brief Creates the wavelet section
532 !> \param section the section to create
533 !> \author fschiff
534 !> \note
535 !> this approach is based on the development of T. Deutsch and S. Goedecker
536 ! **************************************************************************************************
537  SUBROUTINE create_wavelet_section(section)
538  TYPE(section_type), POINTER :: section
539 
540  TYPE(keyword_type), POINTER :: keyword
541 
542  cpassert(.NOT. ASSOCIATED(section))
543  CALL section_create( &
544  section, __location__, name="wavelet", &
545  description="Sets up parameters of wavelet based poisson solver.", &
546  n_keywords=1, n_subsections=0, repeats=.false., &
547  citations=(/genovese2006, genovese2007/))
548 
549  NULLIFY (keyword)
550 
551  CALL keyword_create( &
552  keyword, __location__, name="SCF_TYPE", &
553  description="Type of scaling function used in the wavelet approach, the total energy depends on this choice, "// &
554  "and the convergence with respect to cutoff depends on the selected scaling functions. "// &
555  "Possible values are 8,14,16,20,24,30,40,50,60,100", &
556  usage="SCF_TYPE integer", &
557  n_var=1, default_i_val=40)
558  CALL section_add_keyword(section, keyword)
559  CALL keyword_release(keyword)
560 
561  END SUBROUTINE create_wavelet_section
562 
563 ! **************************************************************************************************
564 !> \brief Creates the section for the implicit (generalized) poisson solver
565 !> \param section the section to be created
566 !> \author Mohammad Hossein Bani-Hashemian
567 ! **************************************************************************************************
568  SUBROUTINE create_implicit_ps_section(section)
569  TYPE(section_type), POINTER :: section
570 
571  TYPE(keyword_type), POINTER :: keyword
572  TYPE(section_type), POINTER :: subsection
573 
574  cpassert(.NOT. ASSOCIATED(section))
575  CALL section_create(section, __location__, name="IMPLICIT", &
576  description="Parameters for the implicit (generalized) Poisson solver.", &
577  citations=(/banihashemian2016/), &
578  n_keywords=6, n_subsections=2, repeats=.false.)
579 
580  NULLIFY (subsection, keyword)
581 
582  CALL create_dielectric_section(subsection)
583  CALL section_add_subsection(section, subsection)
584  CALL section_release(subsection)
585 
586  CALL create_dbc_section(subsection)
587  CALL section_add_subsection(section, subsection)
588  CALL section_release(subsection)
589 
590  CALL keyword_create( &
591  keyword, __location__, name="BOUNDARY_CONDITIONS", &
592  enum_c_vals=s2a('PERIODIC', 'MIXED', 'MIXED_PERIODIC', 'NEUMANN'), &
593  enum_desc=s2a('periodic boundary conditions', 'Dirichlet + homogeneous Neumann boundary conditions', &
594  'Dirichlet + periodic boundary conditions', 'homogeneous Neumann BC (zero-average solution)'), &
595  enum_i_vals=(/periodic_bc, mixed_bc, mixed_periodic_bc, neumann_bc/), &
596  description="Specifies the type of boundary conditions. Dirichlet=fixed value, Neumann=zero normal deriv. "// &
597  "Mixed and Neumann boundaries essentially requires FFTW3 so that all grid sizes are FFT-able.", &
598  usage="BOUNDARY_CONDITIONS <bc_type>", default_i_val=periodic_bc)
599  CALL section_add_keyword(section, keyword)
600  CALL keyword_release(keyword)
601 
602  CALL keyword_create(keyword, __location__, name="ZERO_INITIAL_GUESS", &
603  description="Whether or not to use zero potential as initial guess.", &
604  usage="ZERO_INITIAL_GUESS <logical>", default_l_val=.false., lone_keyword_l_val=.true.)
605  CALL section_add_keyword(section, keyword)
606  CALL keyword_release(keyword)
607 
608  CALL keyword_create(keyword, __location__, name="max_iter", &
609  description="Maximum number of iterations.", &
610  usage="max_iter <integer>", default_i_val=30)
611  CALL section_add_keyword(section, keyword)
612  CALL keyword_release(keyword)
613 
614  CALL keyword_create(keyword, __location__, name="tol", &
615  description="Stopping tolerance.", &
616  usage="tol <real>", default_r_val=1.0e-8_dp)
617  CALL section_add_keyword(section, keyword)
618  CALL keyword_release(keyword)
619 
620  CALL keyword_create(keyword, __location__, name="OR_PARAMETER", variants=s2a('omega'), &
621  description="Over-relaxation parameter (large epsilon requires smaller omega ~0.1).", &
622  usage="OR_PARAMETER <real>", default_r_val=1.0_dp)
623  CALL section_add_keyword(section, keyword)
624  CALL keyword_release(keyword)
625 
626  CALL keyword_create( &
627  keyword, __location__, name="NEUMANN_DIRECTIONS", &
628  enum_c_vals=s2a('XYZ', 'XY', 'XZ', 'YZ', 'X', 'Y', 'Z'), &
630  description="Directions in which homogeneous Neumann conditions are imposed. In the remaining directions "// &
631  "periodic conditions will be enforced. Having specified MIXED or NEUMANN as BOUNDARY_CONDITIONS, "// &
632  "the keyword is meant to be used to combine periodic and homogeneous Neumann conditions at the "// &
633  "boundaries of the simulation cell.", &
634  usage="NEUMANN_DIRECTIONS <direction>", default_i_val=neumannxyz)
635  CALL section_add_keyword(section, keyword)
636  CALL keyword_release(keyword)
637 
638  END SUBROUTINE create_implicit_ps_section
639 
640 ! **************************************************************************************************
641 !> \brief Creates the dielectric constant section.
642 !> The dielectric constant is defined as a function of electronic density.
643 !> [see O. Andreussi, I. Dabo, and N. Marzari, J. Chem. Phys., 136, 064102(2012)]
644 !> \param section the section to be created
645 !> \author Mohammad Hossein Bani-Hashemian
646 ! **************************************************************************************************
647  SUBROUTINE create_dielectric_section(section)
648  TYPE(section_type), POINTER :: section
649 
650  TYPE(keyword_type), POINTER :: keyword
651  TYPE(section_type), POINTER :: subsection
652 
653  cpassert(.NOT. ASSOCIATED(section))
654  CALL section_create(section, __location__, name="DIELECTRIC", &
655  description="Parameters for the dielectric constant function.", &
656  n_keywords=6, n_subsections=2, repeats=.false.)
657 
658  NULLIFY (keyword, subsection)
659 
660  CALL keyword_create(keyword, __location__, name="DIELECTRIC_CORE_CORRECTION", &
661  description="Avoid spurious values of the dielectric constant at the ionic core for pseudopotentials "// &
662  "where the electron density goes to zero at the core (e.g. GTH). "// &
663  "The correction is based on rho_core.", &
664  usage="DIELECTRIC_CORE_CORRECTION <logical>", default_l_val=.true., lone_keyword_l_val=.true.)
665  CALL section_add_keyword(section, keyword)
666  CALL keyword_release(keyword)
667 
668  CALL keyword_create( &
669  keyword, __location__, name="DIELECTRIC_FUNCTION_TYPE", &
670  enum_c_vals=s2a('density_dependent', 'spatially_dependent', 'spatially_rho_dependent'), &
672  enum_desc=s2a("Dielectric constant as a function of the electron density "// &
673  "as e.g. proposed within the SCCS model.", &
674  "Various regions with different dielectric constants.", &
675  "Various regions with different dielectric constants. The dielectric constant decays to 1.0, "// &
676  "wherever the electron density is present."), &
677  description="Preferred type for the dielectric constant function.", &
678  usage="DIELECTRIC_FUNCTION_TYPE <method>", default_i_val=rho_dependent)
679  CALL section_add_keyword(section, keyword)
680  CALL keyword_release(keyword)
681 
682  CALL keyword_create(keyword, __location__, name="dielectric_constant", variants=s2a('epsilon'), &
683  description="Dielectric constant in the bulk of the solvent.", &
684  usage="dielectric_constant <real>", default_r_val=80.0_dp)
685  CALL section_add_keyword(section, keyword)
686  CALL keyword_release(keyword)
687 
688  CALL keyword_create(keyword, __location__, name="rho_min", &
689  description="Lower density threshold.", &
690  usage="rho_min <real>", default_r_val=1.0e-4_dp)
691  CALL section_add_keyword(section, keyword)
692  CALL keyword_release(keyword)
693 
694  CALL keyword_create(keyword, __location__, name="rho_max", &
695  description="Upper density threshold.", &
696  usage="rho_max <real>", default_r_val=1.0e-3_dp)
697  CALL section_add_keyword(section, keyword)
698  CALL keyword_release(keyword)
699 
700  CALL keyword_create( &
701  keyword, __location__, name="DERIVATIVE_METHOD", &
702  enum_c_vals=s2a('fft', 'fft_use_deps', 'fft_use_drho', 'cd3', 'cd5', 'cd7'), &
705  enum_desc=s2a("FFT based deriv of epsilon, without correction (high cutoff needed).", &
706  "FFT based deriv of epsilon, with correction using gradient of epsilon (high cutoff needed).", &
707  "FFT based deriv of epsilon, with correction using gradient of rho (high cutoff needed).", &
708  "3-point central difference derivative.", &
709  "5-point central difference derivative.", &
710  "7-point central difference derivative (recommended)."), &
711  description="Preferred method for evaluating the gradient of ln(eps).", &
712  usage="DERIVATIVE_METHOD <method>", default_i_val=derivative_cd7)
713  CALL section_add_keyword(section, keyword)
714  CALL keyword_release(keyword)
715 
716  CALL create_dielec_aa_cuboidal_section(subsection)
717  CALL section_add_subsection(section, subsection)
718  CALL section_release(subsection)
719 
720  CALL create_dielec_xaa_annular_section(subsection)
721  CALL section_add_subsection(section, subsection)
722  CALL section_release(subsection)
723 
724  END SUBROUTINE create_dielectric_section
725 
726 ! **************************************************************************************************
727 !> \brief Creates the section for creating axis-aligned cuboidal dielectric region.
728 !> \param section the section to be created
729 !> \author Mohammad Hossein Bani-Hashemian
730 ! **************************************************************************************************
731  SUBROUTINE create_dielec_aa_cuboidal_section(section)
732  TYPE(section_type), POINTER :: section
733 
734  TYPE(keyword_type), POINTER :: keyword
735 
736  cpassert(.NOT. ASSOCIATED(section))
737  CALL section_create(section, __location__, name="DIELEC_AA_CUBOIDAL", &
738  description="Parameters for creating axis-aligned cuboidal dielectric region. "// &
739  "Note that once such a region is defined, the 'background' dielectric constant "// &
740  "would be the default (80.0), unless a different value is specified using the "// &
741  "keyword IMPLICIT%DIELECTRIC%DIELECTRIC_CONSTANT.", &
742  n_keywords=5, n_subsections=0, repeats=.true.)
743 
744  NULLIFY (keyword)
745 
746  CALL keyword_create(keyword, __location__, name="dielectric_constant", variants=s2a('epsilon'), &
747  description="value of the dielectric constant inside the region.", &
748  usage="dielectric_constant <real>", default_r_val=80.0_dp)
749  CALL section_add_keyword(section, keyword)
750  CALL keyword_release(keyword)
751 
752  CALL keyword_create(keyword, __location__, name="X_xtnt", &
753  description="The X extents of the cuboid.", &
754  usage="X_xtnt <xmin(real)> <xmax(real)>", unit_str="angstrom", &
755  n_var=2, type_of_var=real_t)
756  CALL section_add_keyword(section, keyword)
757  CALL keyword_release(keyword)
758 
759  CALL keyword_create(keyword, __location__, name="Y_xtnt", &
760  description="The Y extents of the cuboid.", &
761  usage="Y_xtnt <ymin(real)> <ymax(real)>", unit_str="angstrom", &
762  n_var=2, type_of_var=real_t)
763  CALL section_add_keyword(section, keyword)
764  CALL keyword_release(keyword)
765 
766  CALL keyword_create(keyword, __location__, name="Z_xtnt", &
767  description="The Z extents of the cuboid.", &
768  usage="Z_xtnt <zmin(real)> <zmax(real)>", unit_str="angstrom", &
769  n_var=2, type_of_var=real_t)
770  CALL section_add_keyword(section, keyword)
771  CALL keyword_release(keyword)
772 
773  CALL keyword_create(keyword, __location__, name="SMOOTHING_WIDTH", variants=s2a('zeta'), &
774  description="The width of the standard mollifier.", &
775  usage="SMOOTHING_WIDTH <real>", unit_str="angstrom", &
776  default_r_val=cp_unit_to_cp2k(value=0.2_dp, unit_str="angstrom"))
777  CALL section_add_keyword(section, keyword)
778  CALL keyword_release(keyword)
779 
780  END SUBROUTINE create_dielec_aa_cuboidal_section
781 
782 ! **************************************************************************************************
783 !> \brief Creates the section for creating x-axis-aligned annular dielectric region.
784 !> \param section the section to be created
785 !> \author Mohammad Hossein Bani-Hashemian
786 ! **************************************************************************************************
787  SUBROUTINE create_dielec_xaa_annular_section(section)
788  TYPE(section_type), POINTER :: section
789 
790  TYPE(keyword_type), POINTER :: keyword
791 
792  cpassert(.NOT. ASSOCIATED(section))
793  CALL section_create(section, __location__, name="DIELEC_XAA_ANNULAR", &
794  description="Parameters for creating x-axis-aligned annular dielectric region. "// &
795  "Note that once such a region is defined, the 'background' dielectric constant "// &
796  "would be the default (80.0), unless a different value is specified using the "// &
797  "keyword IMPLICIT%DIELECTRIC%DIELECTRIC_CONSTANT.", &
798  n_keywords=5, n_subsections=0, repeats=.true.)
799 
800  NULLIFY (keyword)
801 
802  CALL keyword_create(keyword, __location__, name="dielectric_constant", variants=s2a('epsilon'), &
803  description="value of the dielectric constant inside the region.", &
804  usage="dielectric_constant <real>", default_r_val=80.0_dp)
805  CALL section_add_keyword(section, keyword)
806  CALL keyword_release(keyword)
807 
808  CALL keyword_create(keyword, __location__, name="X_xtnt", &
809  description="The X extents of the annulus.", &
810  usage="X_xtnt <xmin(real)> <xmax(real)>", unit_str="angstrom", &
811  n_var=2, type_of_var=real_t)
812  CALL section_add_keyword(section, keyword)
813  CALL keyword_release(keyword)
814 
815  CALL keyword_create(keyword, __location__, name="base_center", &
816  description="The y and z coordinates of the annulus' base center.", &
817  usage="base_center <y(real)> <z(real)>", unit_str="angstrom", &
818  n_var=2, type_of_var=real_t)
819  CALL section_add_keyword(section, keyword)
820  CALL keyword_release(keyword)
821 
822  CALL keyword_create(keyword, __location__, name="base_radii", &
823  description="The base radius of the annulus.", &
824  usage="base_radius <r1(real)> <r2(real)>", unit_str="angstrom", &
825  n_var=2, type_of_var=real_t)
826  CALL section_add_keyword(section, keyword)
827  CALL keyword_release(keyword)
828 
829  CALL keyword_create(keyword, __location__, name="smoothing_width", variants=s2a('zeta'), &
830  description="The width of the standard mollifier.", &
831  usage="smoothing_width <real>", unit_str="angstrom", &
832  default_r_val=cp_unit_to_cp2k(value=0.2_dp, unit_str="angstrom"))
833  CALL section_add_keyword(section, keyword)
834  CALL keyword_release(keyword)
835 
836  END SUBROUTINE create_dielec_xaa_annular_section
837 
838 ! **************************************************************************************************
839 !> \brief Creates the section for Dirichlet boundary conditions
840 !> \param section the section to be created
841 !> \author Mohammad Hossein Bani-Hashemian
842 ! **************************************************************************************************
843  SUBROUTINE create_dbc_section(section)
844  TYPE(section_type), POINTER :: section
845 
846  TYPE(keyword_type), POINTER :: keyword
847  TYPE(section_type), POINTER :: subsection
848 
849  cpassert(.NOT. ASSOCIATED(section))
850  CALL section_create(section, __location__, name="DIRICHLET_BC", &
851  description="Parameters for creating Dirichlet type boundary conditions.", &
852  n_keywords=1, n_subsections=4, repeats=.false.)
853 
854  NULLIFY (keyword)
855 
856  CALL keyword_create(keyword, __location__, name="VERBOSE_OUTPUT", &
857  description="Print out the coordinates of the vertices defining Dirichlet regions and their "// &
858  "tessellations (in Angstrom), the values of the electrostatic potential at the regions (in a.u.), "// &
859  "and their corresponding evaluated Lagrange multipliers.", &
860  usage="VERBOSE_OUTPUT <logical>", default_l_val=.false., lone_keyword_l_val=.true.)
861  CALL section_add_keyword(section, keyword)
862  CALL keyword_release(keyword)
863 
864  NULLIFY (subsection)
865 
866  CALL create_aa_planar_section(subsection)
867  CALL section_add_subsection(section, subsection)
868  CALL section_release(subsection)
869 
870  CALL create_planar_section(subsection)
871  CALL section_add_subsection(section, subsection)
872  CALL section_release(subsection)
873 
874  CALL create_aa_cylindrical_section(subsection)
875  CALL section_add_subsection(section, subsection)
876  CALL section_release(subsection)
877 
878  CALL create_aa_cuboidal_section(subsection)
879  CALL section_add_subsection(section, subsection)
880  CALL section_release(subsection)
881 
882  END SUBROUTINE create_dbc_section
883 
884 ! **************************************************************************************************
885 !> \brief Creates the section for creating axis-aligned planar Dirichlet BC.
886 !> \param section the section to be created
887 !> \author Mohammad Hossein Bani-Hashemian
888 ! **************************************************************************************************
889  SUBROUTINE create_aa_planar_section(section)
890  TYPE(section_type), POINTER :: section
891 
892  TYPE(keyword_type), POINTER :: keyword
893 
894  cpassert(.NOT. ASSOCIATED(section))
895  CALL section_create(section, __location__, name="AA_PLANAR", &
896  description="Parameters for creating axis-aligned planar (rectangular) Dirichlet boundary regions.", &
897  n_keywords=10, n_subsections=0, repeats=.true.)
898 
899  NULLIFY (keyword)
900 
901  CALL keyword_create(keyword, __location__, name="v_D", &
902  description="The value of the fixed potential to be imposed at the axis-aligned Dirichlet boundary.", &
903  usage="v_D <real>", unit_str="volt", type_of_var=real_t)
904  CALL section_add_keyword(section, keyword)
905  CALL keyword_release(keyword)
906 
907  CALL keyword_create(keyword, __location__, name="OSCILLATING_FRACTION", &
908  description="A fraction of the field can be set to oscilate over time.", &
909  usage="OSCILLATING_FRACTION <real>", default_r_val=0.0_dp, type_of_var=real_t)
910  CALL section_add_keyword(section, keyword)
911  CALL keyword_release(keyword)
912 
913  CALL keyword_create(keyword, __location__, name="FREQUENCY", &
914  description="The frequency with which the oscillating fraction oscillates.", &
915  usage="FREQUENCY <real>", default_r_val=0.0_dp, unit_str="s^-1", type_of_var=real_t)
916  CALL section_add_keyword(section, keyword)
917  CALL keyword_release(keyword)
918 
919  CALL keyword_create(keyword, __location__, name="PHASE", &
920  description="The phase of the oscillattion. A phase of zero corresponds to a cosine function. ", &
921  usage="PHASE <real>", default_r_val=0.0_dp, type_of_var=real_t)
922  CALL section_add_keyword(section, keyword)
923  CALL keyword_release(keyword)
924 
925  CALL keyword_create(keyword, __location__, name="PARALLEL_PLANE", &
926  enum_c_vals=s2a('XY', 'YZ', 'XZ'), &
927  enum_i_vals=(/xy_plane, yz_plane, xz_plane/), &
928  description="The coordinate plane that the region is parallel to.", &
929  usage="PARALLEL_PLANE <plane>", &
930  type_of_var=enum_t)
931  CALL section_add_keyword(section, keyword)
932  CALL keyword_release(keyword)
933 
934  CALL keyword_create(keyword, __location__, name="INTERCEPT", &
935  description="The intercept of the rectangle's plane.", &
936  usage="INTERCEPT <real>", unit_str="angstrom", &
937  type_of_var=real_t)
938  CALL section_add_keyword(section, keyword)
939  CALL keyword_release(keyword)
940 
941  CALL keyword_create(keyword, __location__, name="X_xtnt", &
942  description="The X extents of the rectangle.", &
943  usage="X_xtnt <xmin(real)> <xmax(real)>", unit_str="angstrom", &
944  n_var=2, type_of_var=real_t)
945  CALL section_add_keyword(section, keyword)
946  CALL keyword_release(keyword)
947 
948  CALL keyword_create(keyword, __location__, name="Y_xtnt", &
949  description="The Y extents of the rectangle.", &
950  usage="Y_xtnt <ymin(real)> <ymax(real)>", unit_str="angstrom", &
951  n_var=2, type_of_var=real_t)
952  CALL section_add_keyword(section, keyword)
953  CALL keyword_release(keyword)
954 
955  CALL keyword_create(keyword, __location__, name="Z_xtnt", &
956  description="The Z extents of the rectangle.", &
957  usage="Z_xtnt <zmin(real)> <zmax(real)>", unit_str="angstrom", &
958  n_var=2, type_of_var=real_t)
959  CALL section_add_keyword(section, keyword)
960  CALL keyword_release(keyword)
961 
962  CALL keyword_create(keyword, __location__, name="N_PRTN", &
963  description="The number of partitions in the directions of the unit vectors generating the "// &
964  "corresponding PARALLEL_PLANE (e1, e2 or e3) for tiling the rectangluar region.", &
965  usage="N_PRTN <integer> <integer>", &
966  n_var=2, default_i_vals=(/1, 1/))
967  CALL section_add_keyword(section, keyword)
968  CALL keyword_release(keyword)
969 
970  CALL keyword_create(keyword, __location__, name="THICKNESS", &
971  description="The thickness of the planar Dirichlet region.", &
972  usage="THICKNESS <real>", unit_str="angstrom", &
973  default_r_val=cp_unit_to_cp2k(value=0.75_dp, unit_str="angstrom"))
974  CALL section_add_keyword(section, keyword)
975  CALL keyword_release(keyword)
976 
977  CALL keyword_create(keyword, __location__, name="SMOOTHING_WIDTH", variants=s2a('SIGMA'), &
978  description="The width of the transition region between the Dirichlet subdomain and its surrounding.", &
979  usage="SMOOTHING_WIDTH <real>", unit_str="angstrom", &
980  default_r_val=cp_unit_to_cp2k(value=0.2_dp, unit_str="angstrom"))
981  CALL section_add_keyword(section, keyword)
982  CALL keyword_release(keyword)
983 
984  CALL keyword_create(keyword, __location__, name="PERIODIC_REGION", &
985  description="Whether or not to take into consideration the effects of the periodicity of the "// &
986  "simulation cell (MIXED_PERIODIC bc) for regions defined sufficiently close to the boundaries.", &
987  usage="PERIODIC_REGION <logical>", default_l_val=.false., lone_keyword_l_val=.true.)
988  CALL section_add_keyword(section, keyword)
989  CALL keyword_release(keyword)
990 
991  END SUBROUTINE create_aa_planar_section
992 
993 ! **************************************************************************************************
994 !> \brief Creates the section for creating axis-aligned planar Dirichlet BC.
995 !> \param section the section to be created
996 !> \author Mohammad Hossein Bani-Hashemian
997 ! **************************************************************************************************
998  SUBROUTINE create_planar_section(section)
999  TYPE(section_type), POINTER :: section
1000 
1001  TYPE(keyword_type), POINTER :: keyword
1002 
1003  cpassert(.NOT. ASSOCIATED(section))
1004  CALL section_create( &
1005  section, __location__, name="PLANAR", &
1006  description="Parameters for creating planar (rectangular) Dirichlet boundary regions with given vertices.", &
1007  n_keywords=7, n_subsections=0, repeats=.true.)
1008 
1009  NULLIFY (keyword)
1010 
1011  CALL keyword_create(keyword, __location__, name="v_D", &
1012  description="The value of the fixed potential to be imposed at the planar Dirichlet boundary.", &
1013  usage="v_D <real>", unit_str="volt", type_of_var=real_t)
1014  CALL section_add_keyword(section, keyword)
1015  CALL keyword_release(keyword)
1016 
1017  CALL keyword_create(keyword, __location__, name="OSCILLATING_FRACTION", &
1018  description="A fraction of the field can be set to oscilate over time.", &
1019  usage="OSCILLATING_FRACTION <real>", default_r_val=0.0_dp, type_of_var=real_t)
1020  CALL section_add_keyword(section, keyword)
1021  CALL keyword_release(keyword)
1022 
1023  CALL keyword_create(keyword, __location__, name="FREQUENCY", &
1024  description="The frequency with which the oscillating fraction oscillates.", &
1025  usage="FREQUENCY <real>", default_r_val=0.0_dp, unit_str="s^-1", type_of_var=real_t)
1026  CALL section_add_keyword(section, keyword)
1027  CALL keyword_release(keyword)
1028 
1029  CALL keyword_create(keyword, __location__, name="PHASE", &
1030  description="The phase of the oscillattion. A phase of zero corresponds to a cosine function. ", &
1031  usage="PHASE <real>", default_r_val=0.0_dp, type_of_var=real_t)
1032  CALL section_add_keyword(section, keyword)
1033  CALL keyword_release(keyword)
1034 
1035  CALL keyword_create(keyword, __location__, name="A", &
1036  description="Coordinates of the vertex A.", &
1037  usage="A <x(real)> <y(real)> <z(real)>", unit_str="angstrom", &
1038  n_var=3, type_of_var=real_t)
1039  CALL section_add_keyword(section, keyword)
1040  CALL keyword_release(keyword)
1041 
1042  CALL keyword_create(keyword, __location__, name="B", &
1043  description="Coordinates of the vertex B.", &
1044  usage="B <x(real)> <y(real)> <z(real)>", unit_str="angstrom", &
1045  n_var=3, type_of_var=real_t)
1046  CALL section_add_keyword(section, keyword)
1047  CALL keyword_release(keyword)
1048 
1049  CALL keyword_create(keyword, __location__, name="C", &
1050  description="Coordinates of the vertex C.", &
1051  usage="C <x(real)> <y(real)> <z(real)>", unit_str="angstrom", &
1052  n_var=3, type_of_var=real_t)
1053  CALL section_add_keyword(section, keyword)
1054  CALL keyword_release(keyword)
1055 
1056  CALL keyword_create( &
1057  keyword, __location__, name="N_PRTN", &
1058  description="The number of partitions along the edges for tiling the rectangular region. If the edges "// &
1059  "have different lengths, from the two given values, the larger one will be assigned to the longer edge.", &
1060  usage="N_PRTN <integer> <integer>", &
1061  n_var=2, default_i_vals=(/1, 1/))
1062  CALL section_add_keyword(section, keyword)
1063  CALL keyword_release(keyword)
1064 
1065  CALL keyword_create(keyword, __location__, name="THICKNESS", &
1066  description="The thickness of the planar Dirichlet region.", &
1067  usage="THICKNESS <real>", unit_str="angstrom", &
1068  default_r_val=cp_unit_to_cp2k(value=0.75_dp, unit_str="angstrom"))
1069  CALL section_add_keyword(section, keyword)
1070  CALL keyword_release(keyword)
1071 
1072  CALL keyword_create(keyword, __location__, name="SMOOTHING_WIDTH", variants=s2a('SIGMA'), &
1073  description="The width of the transition region between the Dirichlet subdomain and its surrounding.", &
1074  usage="SMOOTHING_WIDTH <real>", unit_str="angstrom", &
1075  default_r_val=cp_unit_to_cp2k(value=0.2_dp, unit_str="angstrom"))
1076  CALL section_add_keyword(section, keyword)
1077  CALL keyword_release(keyword)
1078 
1079  END SUBROUTINE create_planar_section
1080 
1081 ! **************************************************************************************************
1082 !> \brief Creates the section for creating x-axis-aligned cylindrical Dirichlet BC.
1083 !> \param section the section to be created
1084 !> \author Mohammad Hossein Bani-Hashemian
1085 ! **************************************************************************************************
1086  SUBROUTINE create_aa_cylindrical_section(section)
1087  TYPE(section_type), POINTER :: section
1088 
1089  TYPE(keyword_type), POINTER :: keyword
1090 
1091  cpassert(.NOT. ASSOCIATED(section))
1092  CALL section_create(section, __location__, name="AA_CYLINDRICAL", &
1093  description="Parameters for creating axis-aligned cylindrical Dirichlet boundary regions.", &
1094  n_keywords=11, n_subsections=0, repeats=.true.)
1095 
1096  NULLIFY (keyword)
1097 
1098  CALL keyword_create(keyword, __location__, name="v_D", &
1099  description="The value of the fixed potential to be imposed at the cylindrical Dirichlet boundary.", &
1100  usage="v_D <real>", unit_str="volt", type_of_var=real_t)
1101  CALL section_add_keyword(section, keyword)
1102  CALL keyword_release(keyword)
1103 
1104  CALL keyword_create(keyword, __location__, name="OSCILLATING_FRACTION", &
1105  description="A fraction of the field can be set to oscilate over time.", &
1106  usage="OSCILLATING_FRACTION <real>", default_r_val=0.0_dp, type_of_var=real_t)
1107  CALL section_add_keyword(section, keyword)
1108  CALL keyword_release(keyword)
1109 
1110  CALL keyword_create(keyword, __location__, name="FREQUENCY", &
1111  description="The frequency with which the oscillating fraction oscillates.", &
1112  usage="FREQUENCY <real>", default_r_val=0.0_dp, unit_str="s^-1", type_of_var=real_t)
1113  CALL section_add_keyword(section, keyword)
1114  CALL keyword_release(keyword)
1115 
1116  CALL keyword_create(keyword, __location__, name="PHASE", &
1117  description="The phase of the oscillattion. A phase of zero corresponds to a cosine function. ", &
1118  usage="PHASE <real>", default_r_val=0.0_dp, type_of_var=real_t)
1119  CALL section_add_keyword(section, keyword)
1120  CALL keyword_release(keyword)
1121 
1122  CALL keyword_create(keyword, __location__, name="PARALLEL_AXIS", &
1123  enum_c_vals=s2a('X', 'Y', 'Z'), &
1124  enum_i_vals=(/x_axis, y_axis, z_axis/), &
1125  description="The coordinate axis that the cylindrical region extends along.", &
1126  usage="PARALLEL_AXIS <axis>", &
1127  type_of_var=enum_t)
1128  CALL section_add_keyword(section, keyword)
1129  CALL keyword_release(keyword)
1130 
1131  CALL keyword_create(keyword, __location__, name="xtnt", &
1132  description="The extents of the cylinder along its central axis.", &
1133  usage="xtnt <xmin(real)> <xmax(real)>", unit_str="angstrom", &
1134  n_var=2, type_of_var=real_t)
1135  CALL section_add_keyword(section, keyword)
1136  CALL keyword_release(keyword)
1137 
1138  CALL keyword_create(keyword, __location__, name="BASE_CENTER", &
1139  description="The y and z coordinates (x and z or x and y coordinates, "// &
1140  "depending on the choice of the parallel axis) of the cylinder's base center.", &
1141  usage="BASE_CENTER <y(real)> <z(real)>", unit_str="angstrom", &
1142  n_var=2, type_of_var=real_t)
1143  CALL section_add_keyword(section, keyword)
1144  CALL keyword_release(keyword)
1145 
1146  CALL keyword_create(keyword, __location__, name="BASE_RADIUS", &
1147  description="The base radius of the cylinder.", &
1148  usage="BASE_RADIUS <real>", unit_str="angstrom", &
1149  default_r_val=cp_unit_to_cp2k(value=1.0_dp, unit_str="angstrom"))
1150  CALL section_add_keyword(section, keyword)
1151  CALL keyword_release(keyword)
1152 
1153  CALL keyword_create(keyword, __location__, name="N_SIDES", &
1154  description="The number of sides (faces) of the n-gonal prism approximating the cylinder.", &
1155  usage="N_SIDES <integer>", default_i_val=5)
1156  CALL section_add_keyword(section, keyword)
1157  CALL keyword_release(keyword)
1158 
1159  CALL keyword_create(keyword, __location__, name="APX_TYPE", &
1160  enum_c_vals=s2a('CIRCUMSCRIBED', 'INSCRIBED'), &
1161  enum_i_vals=(/circumscribed, inscribed/), &
1162  description="Specifies the type of the n-gonal prism approximating the cylinder.", &
1163  usage="APX_TYPE <apx_type>", default_i_val=circumscribed)
1164  CALL section_add_keyword(section, keyword)
1165  CALL keyword_release(keyword)
1166 
1167  CALL keyword_create( &
1168  keyword, __location__, name="N_PRTN", &
1169  description="The number of partitions along the face edges of the prism for tiling. If the edges "// &
1170  "have different lengths, from the two given values, the larger one will be assigned to the longer edge.", &
1171  usage="N_PRTN <integer> <integer>", &
1172  n_var=2, default_i_vals=(/1, 1/))
1173  CALL section_add_keyword(section, keyword)
1174  CALL keyword_release(keyword)
1175 
1176  CALL keyword_create(keyword, __location__, name="THICKNESS", &
1177  description="The thickness of the cylinder.", &
1178  usage="THICKNESS <real>", unit_str="angstrom", &
1179  default_r_val=cp_unit_to_cp2k(value=0.75_dp, unit_str="angstrom"))
1180  CALL section_add_keyword(section, keyword)
1181  CALL keyword_release(keyword)
1182 
1183  CALL keyword_create(keyword, __location__, name="SMOOTHING_WIDTH", variants=s2a('SIGMA'), &
1184  description="The width of the transition region between the Dirichlet subdomain and its surrounding.", &
1185  usage="SMOOTHING_WIDTH <real>", unit_str="angstrom", &
1186  default_r_val=cp_unit_to_cp2k(value=0.2_dp, unit_str="angstrom"))
1187  CALL section_add_keyword(section, keyword)
1188  CALL keyword_release(keyword)
1189 
1190  CALL keyword_create( &
1191  keyword, __location__, name="delta_alpha", &
1192  description="A central angle specifying the gap between the faces of the n-gonal prism. To avoide overlap "// &
1193  "between the cuboids (of the given thickness) built on top of the faces, a larger value is required if the"// &
1194  " number of faces (N_SIDES) is quite few and/or the base radius is fairly small.", &
1195  usage="delta_alpha <real>", default_r_val=0.05_dp, unit_str="rad")
1196  CALL section_add_keyword(section, keyword)
1197  CALL keyword_release(keyword)
1198 
1199  END SUBROUTINE create_aa_cylindrical_section
1200 
1201 ! **************************************************************************************************
1202 !> \brief Creates the section for creating axis-aligned cuboidal Dirichlet region.
1203 !> \param section the section to be created
1204 !> \author Mohammad Hossein Bani-Hashemian
1205 ! **************************************************************************************************
1206  SUBROUTINE create_aa_cuboidal_section(section)
1207  TYPE(section_type), POINTER :: section
1208 
1209  TYPE(keyword_type), POINTER :: keyword
1210 
1211  cpassert(.NOT. ASSOCIATED(section))
1212  CALL section_create(section, __location__, name="AA_CUBOIDAL", &
1213  description="Parameters for creating axis-aligned cuboidal Dirichlet regions.", &
1214  n_keywords=7, n_subsections=0, repeats=.true.)
1215 
1216  NULLIFY (keyword)
1217 
1218  CALL keyword_create(keyword, __location__, name="v_D", &
1219  description="The value of the fixed potential to be imposed at the region.", &
1220  usage="v_D <real>", unit_str="volt", type_of_var=real_t)
1221  CALL section_add_keyword(section, keyword)
1222  CALL keyword_release(keyword)
1223 
1224  CALL keyword_create(keyword, __location__, name="OSCILLATING_FRACTION", &
1225  description="A fraction of the field can be set to oscilate over time.", &
1226  usage="OSCILLATING_FRACTION <real>", default_r_val=0.0_dp, type_of_var=real_t)
1227  CALL section_add_keyword(section, keyword)
1228  CALL keyword_release(keyword)
1229 
1230  CALL keyword_create(keyword, __location__, name="FREQUENCY", &
1231  description="The frequency with which the oscillating fraction oscillates.", &
1232  usage="FREQUENCY <real>", default_r_val=0.0_dp, unit_str="s^-1", type_of_var=real_t)
1233  CALL section_add_keyword(section, keyword)
1234  CALL keyword_release(keyword)
1235 
1236  CALL keyword_create(keyword, __location__, name="PHASE", &
1237  description="The phase of the oscillattion. A phase of zero corresponds to a cosine function. ", &
1238  usage="PHASE <real>", default_r_val=0.0_dp, type_of_var=real_t)
1239  CALL section_add_keyword(section, keyword)
1240  CALL keyword_release(keyword)
1241 
1242  CALL keyword_create(keyword, __location__, name="X_xtnt", &
1243  description="The X extents of the cuboid.", &
1244  usage="X_xtnt <xmin(real)> <xmax(real)>", unit_str="angstrom", &
1245  n_var=2, type_of_var=real_t)
1246  CALL section_add_keyword(section, keyword)
1247  CALL keyword_release(keyword)
1248 
1249  CALL keyword_create(keyword, __location__, name="Y_xtnt", &
1250  description="The Y extents of the cuboid.", &
1251  usage="Y_xtnt <ymin(real)> <ymax(real)>", unit_str="angstrom", &
1252  n_var=2, type_of_var=real_t)
1253  CALL section_add_keyword(section, keyword)
1254  CALL keyword_release(keyword)
1255 
1256  CALL keyword_create(keyword, __location__, name="Z_xtnt", &
1257  description="The Z extents of the cuboid.", &
1258  usage="Z_xtnt <zmin(real)> <zmax(real)>", unit_str="angstrom", &
1259  n_var=2, type_of_var=real_t)
1260  CALL section_add_keyword(section, keyword)
1261  CALL keyword_release(keyword)
1262 
1263  CALL keyword_create(keyword, __location__, name="N_PRTN", &
1264  description="The number of partitions in the x, y and z directions for partitioning the cuboid.", &
1265  usage="N_PRTN <integer> <integer> <integer>", &
1266  n_var=3, default_i_vals=(/1, 1, 1/))
1267  CALL section_add_keyword(section, keyword)
1268  CALL keyword_release(keyword)
1269 
1270  CALL keyword_create(keyword, __location__, name="SMOOTHING_WIDTH", variants=s2a('SIGMA'), &
1271  description="The width of the transition region between the Dirichlet subdomain and its surrounding.", &
1272  usage="SMOOTHING_WIDTH <real>", unit_str="angstrom", &
1273  default_r_val=cp_unit_to_cp2k(value=0.2_dp, unit_str="angstrom"))
1274  CALL section_add_keyword(section, keyword)
1275  CALL keyword_release(keyword)
1276 
1277  CALL keyword_create(keyword, __location__, name="PERIODIC_REGION", &
1278  description="Whether or not to take into consideration the effects of the periodicity of the "// &
1279  "simulation cell (MIXED_PERIODIC bc) for regions defined sufficiently close to the boundaries.", &
1280  usage="PERIODIC_REGION <logical>", default_l_val=.false., lone_keyword_l_val=.true.)
1281  CALL section_add_keyword(section, keyword)
1282  CALL keyword_release(keyword)
1283 
1284  END SUBROUTINE create_aa_cuboidal_section
1285 
1286 END MODULE input_cp2k_poisson
collects all references to literature in CP2K as new algorithms / method are included from literature...
Definition: bibliography.F:28
integer, save, public blochl1995
Definition: bibliography.F:43
integer, save, public toukmaji1996
Definition: bibliography.F:43
integer, save, public aguado2003
Definition: bibliography.F:43
integer, save, public genovese2006
Definition: bibliography.F:43
integer, save, public laino2008
Definition: bibliography.F:43
integer, save, public genovese2007
Definition: bibliography.F:43
integer, save, public essmann1995
Definition: bibliography.F:43
integer, save, public darden1993
Definition: bibliography.F:43
integer, save, public ewald1921
Definition: bibliography.F:43
integer, save, public martyna1999
Definition: bibliography.F:43
integer, save, public banihashemian2016
Definition: bibliography.F:43
Handles all functions related to the CELL.
Definition: cell_types.F:15
integer, parameter, public use_perd_xyz
Definition: cell_types.F:42
integer, parameter, public use_perd_y
Definition: cell_types.F:42
integer, parameter, public use_perd_xz
Definition: cell_types.F:42
integer, parameter, public use_perd_x
Definition: cell_types.F:42
integer, parameter, public use_perd_z
Definition: cell_types.F:42
integer, parameter, public use_perd_yz
Definition: cell_types.F:42
integer, parameter, public use_perd_none
Definition: cell_types.F:42
integer, parameter, public use_perd_xy
Definition: cell_types.F:42
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
unit conversion facility
Definition: cp_units.F:30
real(kind=dp) function, public cp_unit_to_cp2k(value, unit_str, defaults, power)
converts to the internal cp2k units to the given unit
Definition: cp_units.F:1150
the type I Discrete Cosine Transform (DCT-I)
Definition: dct.F:16
integer, parameter, public neumannx
Definition: dct.F:64
integer, parameter, public neumannxy
Definition: dct.F:64
integer, parameter, public neumannxz
Definition: dct.F:64
integer, parameter, public neumannxyz
Definition: dct.F:64
integer, parameter, public neumannz
Definition: dct.F:64
integer, parameter, public neumannyz
Definition: dct.F:64
integer, parameter, public neumanny
Definition: dct.F:64
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
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public do_fist_pol_cg
integer, parameter, public do_fist_pol_none
integer, parameter, public do_fist_pol_sc
function that build the poisson section of the input
subroutine, public create_poisson_section(section)
Creates the Poisson section.
subroutine, public create_gspace_interp_section(section)
creates the interpolation section for the periodic QM/MM
subroutine, public create_ewald_section(section)
...
subroutine, public create_rsgrid_section(section)
...
represents keywords in an input
subroutine, public keyword_release(keyword)
releases the given keyword (see doc/ReferenceCounting.html)
subroutine, public keyword_create(keyword, location, name, description, usage, type_of_var, n_var, repeats, variants, default_val, default_l_val, default_r_val, default_lc_val, default_c_val, default_i_val, default_l_vals, default_r_vals, default_c_vals, default_i_vals, lone_keyword_val, lone_keyword_l_val, lone_keyword_r_val, lone_keyword_c_val, lone_keyword_i_val, lone_keyword_l_vals, lone_keyword_r_vals, lone_keyword_c_vals, lone_keyword_i_vals, enum_c_vals, enum_i_vals, enum, enum_strict, enum_desc, unit_str, citations, deprecation_notice, removed)
creates a keyword object
objects that represent the structure of input sections and the data contained in an input section
subroutine, public section_create(section, location, name, description, n_keywords, n_subsections, repeats, citations)
creates a list of keywords
subroutine, public section_add_keyword(section, keyword)
adds a keyword to the given section
subroutine, public section_add_subsection(section, subsection)
adds a subsection to the given section
recursive subroutine, public section_release(section)
releases the given keyword list (see doc/ReferenceCounting.html)
a wrapper for basic fortran types.
integer, parameter, public real_t
integer, parameter, public integer_t
integer, parameter, public enum_t
Defines the basic variable types.
Definition: kinds.F:23
integer, parameter, public dp
Definition: kinds.F:34
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.