Line data Source code
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 : ! **************************************************************************************************
14 : MODULE input_cp2k_poisson
15 : USE bibliography, ONLY: &
16 : Aguado2003, BaniHashemian2016, Blochl1995, Darden1993, Essmann1995, Ewald1921, &
17 : Genovese2006, Genovese2007, Laino2008, Martyna1999, Toukmaji1996
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,&
25 : use_perd_z
26 : USE cp_output_handling, ONLY: add_last_numeric,&
27 : cp_print_key_section_create,&
28 : low_print_level,&
29 : medium_print_level
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: &
39 : derivative_cd3, derivative_cd5, derivative_cd7, derivative_fft, derivative_fft_use_deps, &
40 : derivative_fft_use_drho, rho_dependent, spatially_dependent, spatially_rho_dependent
41 : USE dirichlet_bc_types, ONLY: CIRCUMSCRIBED,&
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,&
50 : do_fist_pol_none,&
51 : do_fist_pol_sc
52 : USE input_cp2k_rsgrid, ONLY: create_rsgrid_section
53 : USE input_keyword_types, ONLY: keyword_create,&
54 : keyword_release,&
55 : keyword_type
56 : USE input_section_types, ONLY: section_add_keyword,&
57 : section_add_subsection,&
58 : section_create,&
59 : section_release,&
60 : section_type
61 : USE input_val_types, ONLY: enum_t,&
62 : integer_t,&
63 : real_t
64 : USE kinds, ONLY: dp
65 : USE multipole_types, ONLY: do_multipole_charge,&
66 : do_multipole_dipole,&
67 : do_multipole_none,&
68 : do_multipole_quadrupole
69 : USE ps_implicit_types, ONLY: MIXED_BC,&
70 : MIXED_PERIODIC_BC,&
71 : NEUMANN_BC,&
72 : PERIODIC_BC
73 : USE pw_poisson_types, ONLY: &
74 : do_ewald_ewald, do_ewald_none, do_ewald_pme, do_ewald_spme, pw_poisson_analytic, &
75 : pw_poisson_implicit, pw_poisson_mt, pw_poisson_multipole, pw_poisson_periodic, &
76 : pw_poisson_wavelet
77 : USE pw_spline_utils, ONLY: no_precond,&
78 : precond_spl3_1,&
79 : precond_spl3_2,&
80 : precond_spl3_3,&
81 : precond_spl3_aint,&
82 : precond_spl3_aint2
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, &
93 : create_gspace_interp_section, &
94 : create_ewald_section
95 : !***
96 : CONTAINS
97 :
98 : ! **************************************************************************************************
99 : !> \brief Creates the Poisson section
100 : !> \param section the section to create
101 : !> \author teo
102 : ! **************************************************************************************************
103 25606 : 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 25606 : CPASSERT(.NOT. ASSOCIATED(section))
110 : CALL section_create(section, __LOCATION__, name="poisson", &
111 : description="Sets up the poisson resolutor.", &
112 25606 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
113 :
114 25606 : 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"), &
120 : enum_i_vals=(/pw_poisson_periodic, pw_poisson_analytic, pw_poisson_mt, pw_poisson_multipole, &
121 : pw_poisson_wavelet, pw_poisson_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 179242 : default_i_val=pw_poisson_periodic)
136 25606 : CALL section_add_keyword(section, keyword)
137 25606 : 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, &
147 : use_perd_xy, use_perd_xz, use_perd_yz, &
148 : use_perd_xyz, use_perd_none/), &
149 25606 : default_i_val=use_perd_xyz)
150 25606 : CALL section_add_keyword(section, keyword)
151 25606 : CALL keyword_release(keyword)
152 :
153 25606 : CALL create_mt_section(subsection)
154 25606 : CALL section_add_subsection(section, subsection)
155 25606 : CALL section_release(subsection)
156 :
157 25606 : CALL create_wavelet_section(subsection)
158 25606 : CALL section_add_subsection(section, subsection)
159 25606 : CALL section_release(subsection)
160 :
161 25606 : CALL create_multipole_section(subsection)
162 25606 : CALL section_add_subsection(section, subsection)
163 25606 : CALL section_release(subsection)
164 :
165 25606 : CALL create_ewald_section(subsection)
166 25606 : CALL section_add_subsection(section, subsection)
167 25606 : CALL section_release(subsection)
168 :
169 25606 : CALL create_implicit_ps_section(subsection)
170 25606 : CALL section_add_subsection(section, subsection)
171 25606 : CALL section_release(subsection)
172 25606 : 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 25606 : 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 25606 : 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 25606 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
191 :
192 25606 : 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 25606 : unit_str="angstrom")
197 25606 : CALL section_add_keyword(section, keyword)
198 25606 : 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 25606 : unit_str="hartree", default_r_val=1.0E-6_dp)
204 25606 : CALL section_add_keyword(section, keyword)
205 25606 : 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 25606 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
211 25606 : CALL section_add_keyword(section, keyword)
212 25606 : 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 25606 : usage="NGRIDS <integer> <integer> <integer> ", n_var=3, default_i_vals=(/50, 50, 50/))
217 25606 : CALL section_add_keyword(section, keyword)
218 25606 : CALL keyword_release(keyword)
219 :
220 25606 : CALL create_gspace_interp_section(subsection)
221 25606 : CALL section_add_subsection(section, subsection)
222 25606 : 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 25606 : print_level=medium_print_level, filename="GSpace-SplInterp")
227 25606 : CALL section_add_subsection(section, subsection)
228 25606 : 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 25606 : print_level=low_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
233 25606 : CALL section_add_subsection(section, subsection)
234 25606 : CALL section_release(subsection)
235 :
236 25606 : 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 25606 : SUBROUTINE create_mt_section(section)
244 : TYPE(section_type), POINTER :: section
245 :
246 : TYPE(keyword_type), POINTER :: keyword
247 :
248 25606 : 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 51212 : citations=(/Martyna1999/))
256 :
257 25606 : 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 25606 : n_var=1, default_r_val=7.0_dp)
263 25606 : CALL section_add_keyword(section, keyword)
264 25606 : 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 25606 : n_var=1, default_r_val=2.0_dp)
272 25606 : CALL section_add_keyword(section, keyword)
273 25606 : CALL keyword_release(keyword)
274 :
275 25606 : END SUBROUTINE create_mt_section
276 :
277 : ! **************************************************************************************************
278 : !> \brief ...
279 : !> \param section will contain the ewald section
280 : !> \author fawzi
281 : ! **************************************************************************************************
282 27510 : 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 27510 : 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., &
292 165060 : citations=(/Ewald1921, Darden1993, Essmann1995, Toukmaji1996, Laino2008/))
293 :
294 27510 : 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 220080 : "SPME is the smooth particle mesh using beta-Euler splines (recommended)"))
313 27510 : CALL section_add_keyword(section, keyword)
314 27510 : 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 27510 : unit_str="hartree", default_r_val=1.0E-6_dp)
324 27510 : CALL section_add_keyword(section, keyword)
325 27510 : 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 27510 : usage="RCUT 5.0", n_var=1, type_of_var=real_t, unit_str="angstrom")
334 27510 : CALL section_add_keyword(section, keyword)
335 27510 : 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 27510 : unit_str='angstrom^-1')
345 27510 : CALL section_add_keyword(section, keyword)
346 27510 : 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 27510 : usage="gmax 25 25 25", n_var=-1, type_of_var=integer_t)
355 27510 : CALL section_add_keyword(section, keyword)
356 27510 : 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 27510 : usage="ns_max 11", default_i_val=11)
361 27510 : CALL section_add_keyword(section, keyword)
362 27510 : 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 27510 : usage="o_spline 6", default_i_val=6)
367 27510 : CALL section_add_keyword(section, keyword)
368 27510 : CALL keyword_release(keyword)
369 :
370 : CALL keyword_create(keyword, __LOCATION__, name="epsilon", &
371 : description="tolerance of gaussians for fft interpolation (PME only)", &
372 27510 : usage="epsilon 1e-6", default_r_val=1.e-6_dp)
373 27510 : CALL section_add_keyword(section, keyword)
374 27510 : CALL keyword_release(keyword)
375 :
376 27510 : NULLIFY (subsection)
377 27510 : CALL create_rsgrid_section(subsection)
378 27510 : CALL section_add_subsection(section, subsection)
379 27510 : CALL section_release(subsection)
380 :
381 27510 : 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 82530 : citations=(/Aguado2003, Laino2008/))
386 :
387 : CALL keyword_create(keyword, __LOCATION__, name="_SECTION_PARAMETERS_", &
388 : description="Controls the activation of the Multipoles", &
389 27510 : usage="&MULTIPOLES T", default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
390 27510 : CALL section_add_keyword(subsection, keyword)
391 27510 : 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"), &
402 : enum_i_vals=(/do_multipole_none, do_multipole_charge, do_multipole_dipole, &
403 27510 : do_multipole_quadrupole/), type_of_var=enum_t)
404 27510 : CALL section_add_keyword(subsection, keyword)
405 27510 : 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 27510 : type_of_var=enum_t, default_i_val=do_fist_pol_none)
420 27510 : CALL section_add_keyword(subsection, keyword)
421 27510 : 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 27510 : n_var=1, default_i_val=0)
428 27510 : CALL section_add_keyword(subsection, keyword)
429 27510 : 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 27510 : n_var=1, default_r_val=0.5e-07_dp)
436 27510 : CALL section_add_keyword(subsection, keyword)
437 27510 : CALL keyword_release(keyword)
438 :
439 27510 : CALL section_add_subsection(section, subsection)
440 27510 : CALL section_release(subsection)
441 :
442 27510 : NULLIFY (subsection)
443 : CALL section_create(subsection, __LOCATION__, name="PRINT", &
444 : description="Controls printing of Ewald properties", &
445 27510 : n_keywords=0, n_subsections=1, repeats=.FALSE.)
446 27510 : 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 27510 : print_level=low_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
450 27510 : CALL section_add_subsection(subsection, print_key)
451 27510 : CALL section_release(print_key)
452 27510 : CALL section_add_subsection(section, subsection)
453 27510 : CALL section_release(subsection)
454 :
455 27510 : 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 42666 : 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 42666 : CPASSERT(.NOT. ASSOCIATED(section))
469 : CALL section_create(section, __LOCATION__, name="interpolator", &
470 : description="controls the interpolation for the G-space term", &
471 42666 : n_keywords=5, n_subsections=0, repeats=.FALSE.)
472 :
473 42666 : 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, &
483 42666 : precond_spl3_aint2, precond_spl3_2, precond_spl3_3/))
484 42666 : CALL section_add_keyword(section, keyword)
485 42666 : 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, &
495 42666 : precond_spl3_aint2, precond_spl3_2, precond_spl3_3/))
496 42666 : CALL section_add_keyword(section, keyword)
497 42666 : 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 42666 : usage="eps_x 1.e-15", default_r_val=1.e-10_dp)
502 42666 : CALL section_add_keyword(section, keyword)
503 42666 : 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 42666 : usage="eps_r 1.e-15", default_r_val=1.e-10_dp)
508 42666 : CALL section_add_keyword(section, keyword)
509 42666 : 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 85332 : usage="max_iter 200", default_i_val=100)
515 42666 : CALL section_add_keyword(section, keyword)
516 42666 : CALL keyword_release(keyword)
517 :
518 42666 : 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 42666 : add_last=add_last_numeric)
525 42666 : CALL section_add_subsection(section, print_key)
526 42666 : CALL section_release(print_key)
527 :
528 42666 : 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 25606 : SUBROUTINE create_wavelet_section(section)
538 : TYPE(section_type), POINTER :: section
539 :
540 : TYPE(keyword_type), POINTER :: keyword
541 :
542 25606 : 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 76818 : citations=(/Genovese2006, Genovese2007/))
548 :
549 25606 : 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 25606 : n_var=1, default_i_val=40)
558 25606 : CALL section_add_keyword(section, keyword)
559 25606 : CALL keyword_release(keyword)
560 :
561 25606 : 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 25606 : 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 25606 : 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 51212 : n_keywords=6, n_subsections=2, repeats=.FALSE.)
579 :
580 25606 : NULLIFY (subsection, keyword)
581 :
582 25606 : CALL create_dielectric_section(subsection)
583 25606 : CALL section_add_subsection(section, subsection)
584 25606 : CALL section_release(subsection)
585 :
586 25606 : CALL create_dbc_section(subsection)
587 25606 : CALL section_add_subsection(section, subsection)
588 25606 : 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 25606 : usage="BOUNDARY_CONDITIONS <bc_type>", default_i_val=PERIODIC_BC)
599 25606 : CALL section_add_keyword(section, keyword)
600 25606 : 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 25606 : usage="ZERO_INITIAL_GUESS <logical>", default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
605 25606 : CALL section_add_keyword(section, keyword)
606 25606 : CALL keyword_release(keyword)
607 :
608 : CALL keyword_create(keyword, __LOCATION__, name="max_iter", &
609 : description="Maximum number of iterations.", &
610 25606 : usage="max_iter <integer>", default_i_val=30)
611 25606 : CALL section_add_keyword(section, keyword)
612 25606 : CALL keyword_release(keyword)
613 :
614 : CALL keyword_create(keyword, __LOCATION__, name="tol", &
615 : description="Stopping tolerance.", &
616 25606 : usage="tol <real>", default_r_val=1.0E-8_dp)
617 25606 : CALL section_add_keyword(section, keyword)
618 25606 : 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 25606 : usage="OR_PARAMETER <real>", default_r_val=1.0_dp)
623 25606 : CALL section_add_keyword(section, keyword)
624 25606 : 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'), &
629 : enum_i_vals=(/neumannXYZ, neumannXY, neumannXZ, neumannYZ, neumannX, neumannY, neumannZ/), &
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 25606 : usage="NEUMANN_DIRECTIONS <direction>", default_i_val=neumannXYZ)
635 25606 : CALL section_add_keyword(section, keyword)
636 25606 : CALL keyword_release(keyword)
637 :
638 25606 : 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 25606 : 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 25606 : CPASSERT(.NOT. ASSOCIATED(section))
654 : CALL section_create(section, __LOCATION__, name="DIELECTRIC", &
655 : description="Parameters for the dielectric constant function.", &
656 25606 : n_keywords=6, n_subsections=2, repeats=.FALSE.)
657 :
658 25606 : 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 25606 : usage="DIELECTRIC_CORE_CORRECTION <logical>", default_l_val=.TRUE., lone_keyword_l_val=.TRUE.)
665 25606 : CALL section_add_keyword(section, keyword)
666 25606 : 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'), &
671 : enum_i_vals=(/rho_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 25606 : usage="DIELECTRIC_FUNCTION_TYPE <method>", default_i_val=rho_dependent)
679 25606 : CALL section_add_keyword(section, keyword)
680 25606 : 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 25606 : usage="dielectric_constant <real>", default_r_val=80.0_dp)
685 25606 : CALL section_add_keyword(section, keyword)
686 25606 : CALL keyword_release(keyword)
687 :
688 : CALL keyword_create(keyword, __LOCATION__, name="rho_min", &
689 : description="Lower density threshold.", &
690 25606 : usage="rho_min <real>", default_r_val=1.0E-4_dp)
691 25606 : CALL section_add_keyword(section, keyword)
692 25606 : CALL keyword_release(keyword)
693 :
694 : CALL keyword_create(keyword, __LOCATION__, name="rho_max", &
695 : description="Upper density threshold.", &
696 25606 : usage="rho_max <real>", default_r_val=1.0E-3_dp)
697 25606 : CALL section_add_keyword(section, keyword)
698 25606 : 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'), &
703 : enum_i_vals=(/derivative_fft, derivative_fft_use_deps, derivative_fft_use_drho, &
704 : derivative_cd3, derivative_cd5, derivative_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 25606 : usage="DERIVATIVE_METHOD <method>", default_i_val=derivative_cd7)
713 25606 : CALL section_add_keyword(section, keyword)
714 25606 : CALL keyword_release(keyword)
715 :
716 25606 : CALL create_dielec_aa_cuboidal_section(subsection)
717 25606 : CALL section_add_subsection(section, subsection)
718 25606 : CALL section_release(subsection)
719 :
720 25606 : CALL create_dielec_xaa_annular_section(subsection)
721 25606 : CALL section_add_subsection(section, subsection)
722 25606 : CALL section_release(subsection)
723 :
724 25606 : 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 25606 : SUBROUTINE create_dielec_aa_cuboidal_section(section)
732 : TYPE(section_type), POINTER :: section
733 :
734 : TYPE(keyword_type), POINTER :: keyword
735 :
736 25606 : 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 25606 : n_keywords=5, n_subsections=0, repeats=.TRUE.)
743 :
744 25606 : 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 25606 : usage="dielectric_constant <real>", default_r_val=80.0_dp)
749 25606 : CALL section_add_keyword(section, keyword)
750 25606 : 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 25606 : n_var=2, type_of_var=real_t)
756 25606 : CALL section_add_keyword(section, keyword)
757 25606 : 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 25606 : n_var=2, type_of_var=real_t)
763 25606 : CALL section_add_keyword(section, keyword)
764 25606 : 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 25606 : n_var=2, type_of_var=real_t)
770 25606 : CALL section_add_keyword(section, keyword)
771 25606 : 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 25606 : default_r_val=cp_unit_to_cp2k(value=0.2_dp, unit_str="angstrom"))
777 25606 : CALL section_add_keyword(section, keyword)
778 25606 : CALL keyword_release(keyword)
779 :
780 25606 : 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 25606 : SUBROUTINE create_dielec_xaa_annular_section(section)
788 : TYPE(section_type), POINTER :: section
789 :
790 : TYPE(keyword_type), POINTER :: keyword
791 :
792 25606 : 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 25606 : n_keywords=5, n_subsections=0, repeats=.TRUE.)
799 :
800 25606 : 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 25606 : usage="dielectric_constant <real>", default_r_val=80.0_dp)
805 25606 : CALL section_add_keyword(section, keyword)
806 25606 : 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 25606 : n_var=2, type_of_var=real_t)
812 25606 : CALL section_add_keyword(section, keyword)
813 25606 : 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 25606 : n_var=2, type_of_var=real_t)
819 25606 : CALL section_add_keyword(section, keyword)
820 25606 : 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 25606 : n_var=2, type_of_var=real_t)
826 25606 : CALL section_add_keyword(section, keyword)
827 25606 : 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 25606 : default_r_val=cp_unit_to_cp2k(value=0.2_dp, unit_str="angstrom"))
833 25606 : CALL section_add_keyword(section, keyword)
834 25606 : CALL keyword_release(keyword)
835 :
836 25606 : 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 25606 : 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 25606 : CPASSERT(.NOT. ASSOCIATED(section))
850 : CALL section_create(section, __LOCATION__, name="DIRICHLET_BC", &
851 : description="Parameters for creating Dirichlet type boundary conditions.", &
852 25606 : n_keywords=1, n_subsections=4, repeats=.FALSE.)
853 :
854 25606 : 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 25606 : usage="VERBOSE_OUTPUT <logical>", default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
861 25606 : CALL section_add_keyword(section, keyword)
862 25606 : CALL keyword_release(keyword)
863 :
864 25606 : NULLIFY (subsection)
865 :
866 25606 : CALL create_aa_planar_section(subsection)
867 25606 : CALL section_add_subsection(section, subsection)
868 25606 : CALL section_release(subsection)
869 :
870 25606 : CALL create_planar_section(subsection)
871 25606 : CALL section_add_subsection(section, subsection)
872 25606 : CALL section_release(subsection)
873 :
874 25606 : CALL create_aa_cylindrical_section(subsection)
875 25606 : CALL section_add_subsection(section, subsection)
876 25606 : CALL section_release(subsection)
877 :
878 25606 : CALL create_aa_cuboidal_section(subsection)
879 25606 : CALL section_add_subsection(section, subsection)
880 25606 : CALL section_release(subsection)
881 :
882 25606 : 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 25606 : SUBROUTINE create_aa_planar_section(section)
890 : TYPE(section_type), POINTER :: section
891 :
892 : TYPE(keyword_type), POINTER :: keyword
893 :
894 25606 : 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 25606 : n_keywords=10, n_subsections=0, repeats=.TRUE.)
898 :
899 25606 : 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 25606 : usage="v_D <real>", unit_str="volt", type_of_var=real_t)
904 25606 : CALL section_add_keyword(section, keyword)
905 25606 : 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 25606 : usage="OSCILLATING_FRACTION <real>", default_r_val=0.0_dp, type_of_var=real_t)
910 25606 : CALL section_add_keyword(section, keyword)
911 25606 : CALL keyword_release(keyword)
912 :
913 : CALL keyword_create(keyword, __LOCATION__, name="FREQUENCY", &
914 : description="The frequency with which the oscillating fraction oscillates.", &
915 25606 : usage="FREQUENCY <real>", default_r_val=0.0_dp, unit_str="s^-1", type_of_var=real_t)
916 25606 : CALL section_add_keyword(section, keyword)
917 25606 : 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 25606 : usage="PHASE <real>", default_r_val=0.0_dp, type_of_var=real_t)
922 25606 : CALL section_add_keyword(section, keyword)
923 25606 : 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 25606 : type_of_var=enum_t)
931 25606 : CALL section_add_keyword(section, keyword)
932 25606 : 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 25606 : type_of_var=real_t)
938 25606 : CALL section_add_keyword(section, keyword)
939 25606 : 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 25606 : n_var=2, type_of_var=real_t)
945 25606 : CALL section_add_keyword(section, keyword)
946 25606 : 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 25606 : n_var=2, type_of_var=real_t)
952 25606 : CALL section_add_keyword(section, keyword)
953 25606 : 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 25606 : n_var=2, type_of_var=real_t)
959 25606 : CALL section_add_keyword(section, keyword)
960 25606 : 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 25606 : n_var=2, default_i_vals=(/1, 1/))
967 25606 : CALL section_add_keyword(section, keyword)
968 25606 : 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 25606 : default_r_val=cp_unit_to_cp2k(value=0.75_dp, unit_str="angstrom"))
974 25606 : CALL section_add_keyword(section, keyword)
975 25606 : 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 25606 : default_r_val=cp_unit_to_cp2k(value=0.2_dp, unit_str="angstrom"))
981 25606 : CALL section_add_keyword(section, keyword)
982 25606 : 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 25606 : usage="PERIODIC_REGION <logical>", default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
988 25606 : CALL section_add_keyword(section, keyword)
989 25606 : CALL keyword_release(keyword)
990 :
991 25606 : 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 25606 : SUBROUTINE create_planar_section(section)
999 : TYPE(section_type), POINTER :: section
1000 :
1001 : TYPE(keyword_type), POINTER :: keyword
1002 :
1003 25606 : 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 25606 : n_keywords=7, n_subsections=0, repeats=.TRUE.)
1008 :
1009 25606 : 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 25606 : usage="v_D <real>", unit_str="volt", type_of_var=real_t)
1014 25606 : CALL section_add_keyword(section, keyword)
1015 25606 : 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 25606 : usage="OSCILLATING_FRACTION <real>", default_r_val=0.0_dp, type_of_var=real_t)
1020 25606 : CALL section_add_keyword(section, keyword)
1021 25606 : CALL keyword_release(keyword)
1022 :
1023 : CALL keyword_create(keyword, __LOCATION__, name="FREQUENCY", &
1024 : description="The frequency with which the oscillating fraction oscillates.", &
1025 25606 : usage="FREQUENCY <real>", default_r_val=0.0_dp, unit_str="s^-1", type_of_var=real_t)
1026 25606 : CALL section_add_keyword(section, keyword)
1027 25606 : 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 25606 : usage="PHASE <real>", default_r_val=0.0_dp, type_of_var=real_t)
1032 25606 : CALL section_add_keyword(section, keyword)
1033 25606 : 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 25606 : n_var=3, type_of_var=real_t)
1039 25606 : CALL section_add_keyword(section, keyword)
1040 25606 : 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 25606 : n_var=3, type_of_var=real_t)
1046 25606 : CALL section_add_keyword(section, keyword)
1047 25606 : 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 25606 : n_var=3, type_of_var=real_t)
1053 25606 : CALL section_add_keyword(section, keyword)
1054 25606 : 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 25606 : n_var=2, default_i_vals=(/1, 1/))
1062 25606 : CALL section_add_keyword(section, keyword)
1063 25606 : 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 25606 : default_r_val=cp_unit_to_cp2k(value=0.75_dp, unit_str="angstrom"))
1069 25606 : CALL section_add_keyword(section, keyword)
1070 25606 : 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 25606 : default_r_val=cp_unit_to_cp2k(value=0.2_dp, unit_str="angstrom"))
1076 25606 : CALL section_add_keyword(section, keyword)
1077 25606 : CALL keyword_release(keyword)
1078 :
1079 25606 : 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 25606 : SUBROUTINE create_aa_cylindrical_section(section)
1087 : TYPE(section_type), POINTER :: section
1088 :
1089 : TYPE(keyword_type), POINTER :: keyword
1090 :
1091 25606 : 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 25606 : n_keywords=11, n_subsections=0, repeats=.TRUE.)
1095 :
1096 25606 : 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 25606 : usage="v_D <real>", unit_str="volt", type_of_var=real_t)
1101 25606 : CALL section_add_keyword(section, keyword)
1102 25606 : 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 25606 : usage="OSCILLATING_FRACTION <real>", default_r_val=0.0_dp, type_of_var=real_t)
1107 25606 : CALL section_add_keyword(section, keyword)
1108 25606 : CALL keyword_release(keyword)
1109 :
1110 : CALL keyword_create(keyword, __LOCATION__, name="FREQUENCY", &
1111 : description="The frequency with which the oscillating fraction oscillates.", &
1112 25606 : usage="FREQUENCY <real>", default_r_val=0.0_dp, unit_str="s^-1", type_of_var=real_t)
1113 25606 : CALL section_add_keyword(section, keyword)
1114 25606 : 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 25606 : usage="PHASE <real>", default_r_val=0.0_dp, type_of_var=real_t)
1119 25606 : CALL section_add_keyword(section, keyword)
1120 25606 : 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 25606 : type_of_var=enum_t)
1128 25606 : CALL section_add_keyword(section, keyword)
1129 25606 : 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 25606 : n_var=2, type_of_var=real_t)
1135 25606 : CALL section_add_keyword(section, keyword)
1136 25606 : 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 25606 : n_var=2, type_of_var=real_t)
1143 25606 : CALL section_add_keyword(section, keyword)
1144 25606 : 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 25606 : default_r_val=cp_unit_to_cp2k(value=1.0_dp, unit_str="angstrom"))
1150 25606 : CALL section_add_keyword(section, keyword)
1151 25606 : 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 25606 : usage="N_SIDES <integer>", default_i_val=5)
1156 25606 : CALL section_add_keyword(section, keyword)
1157 25606 : 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 25606 : usage="APX_TYPE <apx_type>", default_i_val=CIRCUMSCRIBED)
1164 25606 : CALL section_add_keyword(section, keyword)
1165 25606 : 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 25606 : n_var=2, default_i_vals=(/1, 1/))
1173 25606 : CALL section_add_keyword(section, keyword)
1174 25606 : 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 25606 : default_r_val=cp_unit_to_cp2k(value=0.75_dp, unit_str="angstrom"))
1180 25606 : CALL section_add_keyword(section, keyword)
1181 25606 : 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 25606 : default_r_val=cp_unit_to_cp2k(value=0.2_dp, unit_str="angstrom"))
1187 25606 : CALL section_add_keyword(section, keyword)
1188 25606 : 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 25606 : usage="delta_alpha <real>", default_r_val=0.05_dp, unit_str="rad")
1196 25606 : CALL section_add_keyword(section, keyword)
1197 25606 : CALL keyword_release(keyword)
1198 :
1199 25606 : 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 25606 : SUBROUTINE create_aa_cuboidal_section(section)
1207 : TYPE(section_type), POINTER :: section
1208 :
1209 : TYPE(keyword_type), POINTER :: keyword
1210 :
1211 25606 : CPASSERT(.NOT. ASSOCIATED(section))
1212 : CALL section_create(section, __LOCATION__, name="AA_CUBOIDAL", &
1213 : description="Parameters for creating axis-aligned cuboidal Dirichlet regions.", &
1214 25606 : n_keywords=7, n_subsections=0, repeats=.TRUE.)
1215 :
1216 25606 : 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 25606 : usage="v_D <real>", unit_str="volt", type_of_var=real_t)
1221 25606 : CALL section_add_keyword(section, keyword)
1222 25606 : 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 25606 : usage="OSCILLATING_FRACTION <real>", default_r_val=0.0_dp, type_of_var=real_t)
1227 25606 : CALL section_add_keyword(section, keyword)
1228 25606 : CALL keyword_release(keyword)
1229 :
1230 : CALL keyword_create(keyword, __LOCATION__, name="FREQUENCY", &
1231 : description="The frequency with which the oscillating fraction oscillates.", &
1232 25606 : usage="FREQUENCY <real>", default_r_val=0.0_dp, unit_str="s^-1", type_of_var=real_t)
1233 25606 : CALL section_add_keyword(section, keyword)
1234 25606 : 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 25606 : usage="PHASE <real>", default_r_val=0.0_dp, type_of_var=real_t)
1239 25606 : CALL section_add_keyword(section, keyword)
1240 25606 : 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 25606 : n_var=2, type_of_var=real_t)
1246 25606 : CALL section_add_keyword(section, keyword)
1247 25606 : 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 25606 : n_var=2, type_of_var=real_t)
1253 25606 : CALL section_add_keyword(section, keyword)
1254 25606 : 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 25606 : n_var=2, type_of_var=real_t)
1260 25606 : CALL section_add_keyword(section, keyword)
1261 25606 : 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 25606 : n_var=3, default_i_vals=(/1, 1, 1/))
1267 25606 : CALL section_add_keyword(section, keyword)
1268 25606 : 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 25606 : default_r_val=cp_unit_to_cp2k(value=0.2_dp, unit_str="angstrom"))
1274 25606 : CALL section_add_keyword(section, keyword)
1275 25606 : 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 25606 : usage="PERIODIC_REGION <logical>", default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1281 25606 : CALL section_add_keyword(section, keyword)
1282 25606 : CALL keyword_release(keyword)
1283 :
1284 25606 : END SUBROUTINE create_aa_cuboidal_section
1285 :
1286 : END MODULE input_cp2k_poisson
|