(git:cccd2f3)
Loading...
Searching...
No Matches
input_cp2k_dft.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2026 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 dft section of the input
10!> \par History
11!> 10.2005 moved out of input_cp2k [fawzi]
12!> \author fawzi
13! **************************************************************************************************
17 USE bibliography, ONLY: &
27 USE cp_spline_utils, ONLY: pw_interp,&
30 USE cp_units, ONLY: cp_unit_to_cp2k
31 USE input_constants, ONLY: &
87 USE input_val_types, ONLY: char_t,&
88 integer_t,&
89 lchar_t,&
90 logical_t,&
91 real_t
92 USE kinds, ONLY: dp
93 USE physcon, ONLY: evolt,&
95 USE pw_spline_utils, ONLY: no_precond,&
101 USE string_utilities, ONLY: s2a
102#include "./base/base_uses.f90"
103
104 IMPLICIT NONE
105 PRIVATE
106
107 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'input_cp2k_dft'
108
109 PUBLIC :: create_dft_section
110 PUBLIC :: create_bsse_section
111 PUBLIC :: create_interp_section
112 PUBLIC :: create_mgrid_section
113
114CONTAINS
115
116! **************************************************************************************************
117!> \brief creates the dft section
118!> \param section the section to be created
119!> \author fawzi
120! **************************************************************************************************
121 SUBROUTINE create_dft_section(section)
122 TYPE(section_type), POINTER :: section
123
124 TYPE(keyword_type), POINTER :: keyword
125 TYPE(section_type), POINTER :: subsection
126
127 cpassert(.NOT. ASSOCIATED(section))
128 CALL section_create(section, __location__, name="DFT", &
129 description="Controls electronic-structure settings for Quickstep and related "// &
130 "Gaussian-basis DFT methods.", &
131 n_keywords=3, n_subsections=4, repeats=.false.)
132
133 NULLIFY (keyword)
134 CALL keyword_create(keyword, __location__, name="BASIS_SET_FILE_NAME", &
135 description="Name of a basis-set library file, optionally including a path. "// &
136 "This keyword can be repeated to search several basis-set files.", &
137 usage="BASIS_SET_FILE_NAME <FILENAME>", &
138 type_of_var=lchar_t, repeats=.true., &
139 default_lc_val="BASIS_SET", n_var=1)
140 CALL section_add_keyword(section, keyword)
141 CALL keyword_release(keyword)
142
143 CALL keyword_create(keyword, __location__, name="POTENTIAL_FILE_NAME", &
144 description="Name of the pseudopotential library file, optionally including a path. "// &
145 "The potential selected for each kind is set with KIND%POTENTIAL.", &
146 usage="POTENTIAL_FILE_NAME <FILENAME>", &
147 default_lc_val="POTENTIAL")
148 CALL section_add_keyword(section, keyword)
149 CALL keyword_release(keyword)
150
151 CALL keyword_create(keyword, __location__, name="WFN_RESTART_FILE_NAME", &
152 variants=["RESTART_FILE_NAME"], &
153 description="Name of the wavefunction restart file, may include a path."// &
154 " If no file is specified, the default is to open the file as generated by the wfn restart print key.", &
155 usage="WFN_RESTART_FILE_NAME <FILENAME>", &
156 type_of_var=lchar_t)
157 CALL section_add_keyword(section, keyword)
158 CALL keyword_release(keyword)
159
160 CALL keyword_create(keyword, __location__, &
161 name="UKS", &
162 variants=s2a("UNRESTRICTED_KOHN_SHAM", &
163 "LSD", &
164 "SPIN_POLARIZED"), &
165 description="Requests a spin-polarized calculation using alpha "// &
166 "and beta orbitals, i.e. no spin restriction is applied", &
167 usage="LSD", &
168 default_l_val=.false., &
169 lone_keyword_l_val=.true.)
170 CALL section_add_keyword(section, keyword)
171 CALL keyword_release(keyword)
172 CALL keyword_create(keyword, __location__, &
173 name="ROKS", &
174 variants=["RESTRICTED_OPEN_KOHN_SHAM"], &
175 description="Requests a restricted open Kohn-Sham calculation", &
176 usage="ROKS", &
177 default_l_val=.false., &
178 lone_keyword_l_val=.true.)
179 CALL section_add_keyword(section, keyword)
180 CALL keyword_release(keyword)
181 CALL keyword_create(keyword, __location__, &
182 name="MULTIPLICITY", &
183 variants=["MULTIP"], &
184 description="Two times the total spin plus one. "// &
185 "Specify 3 for a triplet, 4 for a quartet, "// &
186 "and so on. Default is 1 (singlet) for an "// &
187 "even number and 2 (doublet) for an odd number "// &
188 "of electrons.", &
189 usage="MULTIPLICITY 3", &
190 default_i_val=0) ! this default value is just a flag to get the above
191 CALL section_add_keyword(section, keyword)
192 CALL keyword_release(keyword)
193 CALL keyword_create(keyword, __location__, name="CHARGE", &
194 description="The total charge of the system", &
195 usage="CHARGE -1", &
196 default_i_val=0)
197 CALL section_add_keyword(section, keyword)
198 CALL keyword_release(keyword)
199
200 CALL keyword_create(keyword, __location__, &
201 name="PLUS_U_METHOD", &
202 description="Method employed for the calculation of the DFT+U contribution", &
203 repeats=.false., &
204 enum_c_vals=s2a("LOWDIN", "MULLIKEN", "MULLIKEN_CHARGES"), &
206 enum_desc=s2a("Method based on Lowdin population analysis "// &
207 "(computationally expensive, since the diagonalization of the "// &
208 "overlap matrix is required, but possibly more robust than Mulliken)", &
209 "Method based on Mulliken population analysis using the net AO and "// &
210 "overlap populations (computationally cheap method)", &
211 "Method based on Mulliken gross orbital populations (GOP)"), &
212 n_var=1, &
213 default_i_val=plus_u_mulliken, &
214 usage="PLUS_U_METHOD Lowdin")
215 CALL section_add_keyword(section, keyword)
216 CALL keyword_release(keyword)
217
218 CALL keyword_create(keyword, __location__, &
219 name="RELAX_MULTIPLICITY", &
220 variants=["RELAX_MULTIP"], &
221 description="Tolerance in Hartrees. Do not enforce the occupation "// &
222 "of alpha and beta MOs due to the initially "// &
223 "defined multiplicity, but rather follow the Aufbau principle. "// &
224 "A value greater than zero activates this option. "// &
225 "If alpha/beta MOs differ in energy less than this tolerance, "// &
226 "then alpha-MO occupation is preferred even if it is higher "// &
227 "in energy (within the tolerance). "// &
228 "Such spin-symmetry broken (spin-polarized) occupation is used "// &
229 "as SCF input, which (is assumed to) bias the SCF "// &
230 "towards a spin-polarized solution. "// &
231 "Thus, larger tolerance increases chances of ending up "// &
232 "with spin-polarization. "// &
233 "This option is only valid for unrestricted (i.e. spin polarised) "// &
234 "Kohn-Sham (UKS) calculations. It also needs non-zero "// &
235 "[ADDED_MOS](#CP2K_INPUT.FORCE_EVAL.DFT.SCF.ADDED_MOS) to actually affect the calculations, "// &
236 "which is why it is not expected to work with [OT](#CP2K_INPUT.FORCE_EVAL.DFT.SCF.OT) "// &
237 "and may raise errors when used with OT. "// &
238 "For more details see [this discussion](https://github.com/cp2k/cp2k/issues/4389).", &
239 usage="RELAX_MULTIPLICITY 0.00001", &
240 repeats=.false., &
241 default_r_val=0.0_dp)
242 CALL section_add_keyword(section, keyword)
243 CALL keyword_release(keyword)
244
245 CALL keyword_create(keyword, __location__, name="SUBCELLS", &
246 description="Read the grid size for subcell generation in the construction of "// &
247 "neighbor lists.", usage="SUBCELLS 1.5", &
248 n_var=1, default_r_val=2.0_dp)
249 CALL section_add_keyword(section, keyword)
250 CALL keyword_release(keyword)
251
252 CALL keyword_create(keyword, __location__, name="AUTO_BASIS", &
253 description="Specify size of automatically generated auxiliary (RI) basis sets: "// &
254 "Options={small,medium,large,huge}", &
255 usage="AUTO_BASIS {basis_type} {basis_size}", &
256 type_of_var=char_t, repeats=.true., n_var=-1, default_c_vals=["X", "X"])
257 CALL section_add_keyword(section, keyword)
258 CALL keyword_release(keyword)
259
260 CALL keyword_create(keyword, __location__, &
261 name="SURFACE_DIPOLE_CORRECTION", &
262 variants=s2a("SURFACE_DIPOLE", &
263 "SURF_DIP"), &
264 description="For slab calculations with asymmetric geometries, activate the correction of "// &
265 "the electrostatic potential with "// &
266 "by compensating for the surface dipole. Implemented only for slabs with normal "// &
267 "parallel to one Cartesian axis. The normal direction is given by the keyword SURF_DIP_DIR", &
268 usage="SURF_DIP", &
269 default_l_val=.false., &
270 lone_keyword_l_val=.true., &
271 citations=[bengtsson1999])
272 CALL section_add_keyword(section, keyword)
273 CALL keyword_release(keyword)
274
275 CALL keyword_create(keyword, __location__, &
276 name="SURF_DIP_DIR", &
277 description="Cartesian axis parallel to surface normal.", &
278 enum_c_vals=s2a("X", "Y", "Z"), &
279 enum_i_vals=[1, 2, 3], &
280 enum_desc=s2a("Along x", "Along y", "Along z"), &
281 n_var=1, &
282 default_i_val=3, &
283 usage="SURF_DIP_DIR Z")
284 CALL section_add_keyword(section, keyword)
285 CALL keyword_release(keyword)
286
287 CALL keyword_create(keyword, __location__, &
288 name="SURF_DIP_POS", &
289 description="This keyword assigns an user defined position in Angstroms "// &
290 "in the direction normal to the surface (given by SURF_DIP_DIR). "// &
291 "The default value is -1.0_dp which appplies the correction at a position "// &
292 "that has minimum electron density on the grid.", &
293 usage="SURF_DIP_POS -1.0_dp", &
294 default_r_val=-1.0_dp)
295 CALL section_add_keyword(section, keyword)
296 CALL keyword_release(keyword)
297
298 CALL keyword_create(keyword, __location__, &
299 name="SURF_DIP_SWITCH", &
300 description="WARNING: Experimental feature under development that will help the "// &
301 "user to switch parameters to facilitate SCF convergence. In its current form the "// &
302 "surface dipole correction is switched off if the calculation does not converge in "// &
303 "(0.5*MAX_SCF + 1) outer_scf steps. "// &
304 "The default value is .FALSE.", &
305 usage="SURF_DIP_SWITCH .TRUE.", &
306 default_l_val=.false., &
307 lone_keyword_l_val=.true.)
308 CALL section_add_keyword(section, keyword)
309 CALL keyword_release(keyword)
310
311 CALL keyword_create(keyword, __location__, &
312 name="CORE_CORR_DIP", &
313 description="If the total CORE_CORRECTION is non-zero and surface dipole "// &
314 "correction is switched on, presence of this keyword will adjust electron "// &
315 "density via MO occupation to reflect the total CORE_CORRECTION. "// &
316 "The default value is .FALSE.", &
317 usage="CORE_CORR_DIP .TRUE.", &
318 default_l_val=.false., &
319 lone_keyword_l_val=.true.)
320 CALL section_add_keyword(section, keyword)
321 CALL keyword_release(keyword)
322
323 CALL keyword_create(keyword, __location__, &
324 name="SORT_BASIS", &
325 description="Sorts basis functions according to a selected criterion. "// &
326 "Sorting by exponent can improve data locality for selected exact-exchange and RI workflows.", &
327 enum_c_vals=s2a("DEFAULT", "EXP"), &
328 enum_i_vals=[basis_sort_default, basis_sort_zet], &
329 enum_desc=s2a("don't sort", "sort w.r.t. exponent"), &
330 default_i_val=basis_sort_default, &
331 usage="SORT_BASIS EXP")
332 CALL section_add_keyword(section, keyword)
333 CALL keyword_release(keyword)
334
335 NULLIFY (subsection)
336 CALL create_scf_section(subsection)
337 CALL section_add_subsection(section, subsection)
338 CALL section_release(subsection)
339
340 CALL create_ls_scf_section(subsection)
341 CALL section_add_subsection(section, subsection)
342 CALL section_release(subsection)
343
344 CALL create_almo_scf_section(subsection)
345 CALL section_add_subsection(section, subsection)
346 CALL section_release(subsection)
347
348 CALL create_kg_section(subsection)
349 CALL section_add_subsection(section, subsection)
350 CALL section_release(subsection)
351
352 CALL create_harris_section(subsection)
353 CALL section_add_subsection(section, subsection)
354 CALL section_release(subsection)
355
356 CALL create_ec_section(subsection)
357 CALL section_add_subsection(section, subsection)
358 CALL section_release(subsection)
359
360 CALL create_exstate_section(subsection)
361 CALL section_add_subsection(section, subsection)
362 CALL section_release(subsection)
363
364 CALL create_admm_section(subsection)
365 CALL section_add_subsection(section, subsection)
366 CALL section_release(subsection)
367
368 CALL create_qs_section(subsection)
369 CALL section_add_subsection(section, subsection)
370 CALL section_release(subsection)
371
372 CALL create_mgrid_section(subsection, create_subsections=.true.)
373 CALL section_add_subsection(section, subsection)
374 CALL section_release(subsection)
375
376 CALL create_xc_section(subsection)
377 CALL section_add_subsection(section, subsection)
378 CALL section_release(subsection)
379
380 CALL create_relativistic_section(subsection)
381 CALL section_add_subsection(section, subsection)
382 CALL section_release(subsection)
383
384 CALL create_sic_section(subsection)
385 CALL section_add_subsection(section, subsection)
386 CALL section_release(subsection)
387
388 CALL create_low_spin_roks_section(subsection)
389 CALL section_add_subsection(section, subsection)
390 CALL section_release(subsection)
391
392 CALL create_efield_section(subsection)
393 CALL section_add_subsection(section, subsection)
394 CALL section_release(subsection)
395
396 CALL create_per_efield_section(subsection)
397 CALL section_add_subsection(section, subsection)
398 CALL section_release(subsection)
399
400 CALL create_ext_pot_section(subsection)
401 CALL section_add_subsection(section, subsection)
402 CALL section_release(subsection)
403
404 CALL create_transport_section(subsection)
405 CALL section_add_subsection(section, subsection)
406 CALL section_release(subsection)
407
408 ! ZMP sections to include the external density or v_xc potential
409 CALL create_ext_den_section(subsection)
410 CALL section_add_subsection(section, subsection)
411 CALL section_release(subsection)
412
413 CALL create_ext_vxc_section(subsection)
414 CALL section_add_subsection(section, subsection)
415 CALL section_release(subsection)
416
417 CALL create_poisson_section(subsection)
418 CALL section_add_subsection(section, subsection)
419 CALL section_release(subsection)
420
421 CALL create_kpoints_section(subsection)
422 CALL section_add_subsection(section, subsection)
423 CALL section_release(subsection)
424
425 CALL create_kpoint_set_section(subsection)
426 CALL section_add_subsection(section, subsection)
427 CALL section_release(subsection)
428
429 CALL create_implicit_solv_section(subsection)
430 CALL section_add_subsection(section, subsection)
431 CALL section_release(subsection)
432
433 CALL create_density_fitting_section(subsection)
434 CALL section_add_subsection(section, subsection)
435 CALL section_release(subsection)
436
437 CALL create_xas_section(subsection)
438 CALL section_add_subsection(section, subsection)
439 CALL section_release(subsection)
440
441 CALL create_xas_tdp_section(subsection)
442 CALL section_add_subsection(section, subsection)
443 CALL section_release(subsection)
444
445 CALL create_localize_section(subsection)
446 CALL section_add_subsection(section, subsection)
447 CALL section_release(subsection)
448
449 CALL create_rtp_section(subsection)
450 CALL section_add_subsection(section, subsection)
451 CALL section_release(subsection)
452
453 CALL create_print_dft_section(subsection)
454 CALL section_add_subsection(section, subsection)
455 CALL section_release(subsection)
456
457 CALL create_sccs_section(subsection)
458 CALL section_add_subsection(section, subsection)
459 CALL section_release(subsection)
460
461 CALL create_active_space_section(subsection)
462 CALL section_add_subsection(section, subsection)
463 CALL section_release(subsection)
464
465 CALL create_dft_smeagol_section(subsection)
466 CALL section_add_subsection(section, subsection)
467 CALL section_release(subsection)
468
469 CALL create_hairy_probes_section(subsection)
470 CALL section_add_subsection(section, subsection)
471 CALL section_release(subsection)
472
473 END SUBROUTINE create_dft_section
474
475! **************************************************************************************************
476!> \brief Hairy Probe DFT Model
477!> \param section ...
478!> \author Margherita Buraschi
479! **************************************************************************************************
480
481 SUBROUTINE create_hairy_probes_section(section)
482 TYPE(section_type), POINTER :: section
483
484 TYPE(keyword_type), POINTER :: keyword
485
486 NULLIFY (keyword)
487 cpassert(.NOT. ASSOCIATED(section))
488 CALL section_create(section, __location__, &
489 name="HAIRY_PROBES", &
490 description="Sets up a Hairy Probe calculation. ", &
491 n_keywords=0, n_subsections=0, repeats=.true.)
492
493 CALL keyword_create(keyword, __location__, &
494 name="_SECTION_PARAMETERS_", &
495 description="Controls the activation of hairy probe", &
496 usage="&HAIRY_PROBES ON", &
497 default_l_val=.false., &
498 lone_keyword_l_val=.true.)
499 CALL section_add_keyword(section, keyword)
500 CALL keyword_release(keyword)
501
502 CALL keyword_create(keyword, __location__, name="ATOM_IDS", &
503 description="Indexes of the atoms to which the probes are attached.", &
504 usage="ATOM_IDS <INTEGER> .. <INTEGER>", &
505 type_of_var=integer_t, n_var=-1)
506 CALL section_add_keyword(section, keyword)
507 CALL keyword_release(keyword)
508
509 CALL keyword_create(keyword, __location__, name="T", &
510 description="Electronic temperature [K]", &
511 usage="T <REAL>", &
512 default_r_val=cp_unit_to_cp2k(value=300.0_dp, unit_str="K"), &
513 unit_str="K")
514 CALL section_add_keyword(section, keyword)
515 CALL keyword_release(keyword)
516
517 CALL keyword_create(keyword, __location__, name="MU", &
518 description="Chemical potential of the electrons in the probes [eV] ", &
519 usage="MU <REAL>", &
520 default_r_val=cp_unit_to_cp2k(value=0.0_dp, unit_str="eV"), &
521 unit_str="eV")
522 CALL section_add_keyword(section, keyword)
523 CALL keyword_release(keyword)
524
525 CALL keyword_create(keyword, __location__, name="ALPHA", &
526 description="Parameter for solution probes ", &
527 usage="ALPHA <REAL>", &
528 default_r_val=1.0_dp)
529 CALL section_add_keyword(section, keyword)
530 CALL keyword_release(keyword)
531
532 CALL keyword_create(keyword, __location__, name="EPS_HP", &
533 description=" Tolerance for accuracy checks on occupation numbers "// &
534 "calculated using hair-probes. ", &
535 usage="EPS_HP <REAL>", &
536 default_r_val=1.0e-5_dp)
537 CALL section_add_keyword(section, keyword)
538 CALL keyword_release(keyword)
539 END SUBROUTINE create_hairy_probes_section
540!####################################################################################
541
542! **************************************************************************************************
543!> \brief Implicit Solvation Model
544!> \param section ...
545!> \author tlaino
546! **************************************************************************************************
547 SUBROUTINE create_implicit_solv_section(section)
548 TYPE(section_type), POINTER :: section
549
550 TYPE(keyword_type), POINTER :: keyword
551 TYPE(section_type), POINTER :: print_key, subsection
552
553 NULLIFY (keyword, subsection, print_key)
554 cpassert(.NOT. ASSOCIATED(section))
555 CALL section_create(section, __location__, name="SCRF", &
556 description="Adds an implicit solvation model to the DFT calculation."// &
557 " Know also as Self Consistent Reaction Field.", &
558 n_keywords=0, n_subsections=0, repeats=.false.)
559
560 CALL keyword_create(keyword, __location__, name="EPS_OUT", &
561 description="Value of the dielectric constant outside the sphere", &
562 usage="EPS_OUT <REAL>", &
563 default_r_val=1.0_dp)
564 CALL section_add_keyword(section, keyword)
565 CALL keyword_release(keyword)
566
567 CALL keyword_create(keyword, __location__, name="LMAX", &
568 description="Maximum value of L used in the multipole expansion", &
569 usage="LMAX <INTEGER>", &
570 default_i_val=3)
571 CALL section_add_keyword(section, keyword)
572 CALL keyword_release(keyword)
573
574 CALL create_sphere_section(subsection)
575 CALL section_add_subsection(section, subsection)
576 CALL section_release(subsection)
577
578 CALL cp_print_key_section_create(print_key, __location__, "program_run_info", &
579 description="Controls the printing basic info about the method", &
580 print_level=low_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
581 CALL section_add_subsection(section, print_key)
582 CALL section_release(print_key)
583
584 END SUBROUTINE create_implicit_solv_section
585
586! **************************************************************************************************
587!> \brief Create Sphere cavity
588!> \param section ...
589!> \author tlaino
590! **************************************************************************************************
591 SUBROUTINE create_sphere_section(section)
592 TYPE(section_type), POINTER :: section
593
594 TYPE(keyword_type), POINTER :: keyword
595 TYPE(section_type), POINTER :: subsection
596
597 NULLIFY (keyword, subsection)
598 cpassert(.NOT. ASSOCIATED(section))
599 CALL section_create(section, __location__, name="SPHERE", &
600 description="Treats the implicit solvent environment like a sphere", &
601 n_keywords=0, n_subsections=0, repeats=.false.)
602
603 CALL keyword_create(keyword, __location__, name="RADIUS", &
604 description="Value of the spherical cavity in the dielectric medium", &
605 usage="RADIUS <REAL>", &
606 unit_str="angstrom", &
607 type_of_var=real_t)
608 CALL section_add_keyword(section, keyword)
609 CALL keyword_release(keyword)
610
611 CALL create_center_section(subsection)
612 CALL section_add_subsection(section, subsection)
613 CALL section_release(subsection)
614
615 END SUBROUTINE create_sphere_section
616
617! **************************************************************************************************
618!> \brief ...
619!> \param section ...
620!> \author tlaino
621! **************************************************************************************************
622 SUBROUTINE create_center_section(section)
623 TYPE(section_type), POINTER :: section
624
625 TYPE(keyword_type), POINTER :: keyword
626
627 NULLIFY (keyword)
628 cpassert(.NOT. ASSOCIATED(section))
629 CALL section_create(section, __location__, name="CENTER", &
630 description="Defines the center of the sphere.", &
631 n_keywords=0, n_subsections=0, repeats=.false.)
632 CALL keyword_create(keyword, __location__, name="XYZ", &
633 description="Coordinates of the center of the sphere", &
634 usage="XYZ <REAL> <REAL> <REAL>", &
635 unit_str="angstrom", &
636 type_of_var=real_t, n_var=3)
637 CALL section_add_keyword(section, keyword)
638 CALL keyword_release(keyword)
639
640 CALL keyword_create(keyword, __location__, name="ATOM_LIST", &
641 description="Defines a list of atoms to define the center of the sphere", &
642 usage="ATOM_LIST <INTEGER> .. <INTEGER>", &
643 type_of_var=integer_t, n_var=-1)
644 CALL section_add_keyword(section, keyword)
645 CALL keyword_release(keyword)
646
647 CALL keyword_create(keyword, __location__, name="WEIGHT_TYPE", &
648 description="Defines the weight used to define the center of the sphere"// &
649 " (if ATOM_LIST is provided)", &
650 usage="WEIGHT_TYPE (UNIT|MASS)", &
651 enum_c_vals=["UNIT", "MASS"], &
652 enum_i_vals=[weight_type_unit, weight_type_mass], &
653 default_i_val=weight_type_unit)
654 CALL section_add_keyword(section, keyword)
655 CALL keyword_release(keyword)
656
657 CALL keyword_create(keyword, __location__, name="FIXED", &
658 description="Specify if the center of the sphere should be fixed or"// &
659 " allowed to move", &
660 usage="FIXED <LOGICAL>", &
661 default_l_val=.true.)
662 CALL section_add_keyword(section, keyword)
663 CALL keyword_release(keyword)
664
665 END SUBROUTINE create_center_section
666
667! **************************************************************************************************
668!> \brief ...
669!> \param section ...
670! **************************************************************************************************
671 SUBROUTINE create_admm_section(section)
672 TYPE(section_type), POINTER :: section
673
674 TYPE(keyword_type), POINTER :: keyword
675
676 NULLIFY (keyword)
677 cpassert(.NOT. ASSOCIATED(section))
678 CALL section_create(section, __location__, name="AUXILIARY_DENSITY_MATRIX_METHOD", &
679 description="Controls the auxiliary density matrix method (ADMM), which evaluates "// &
680 "Hartree-Fock exchange on a smaller auxiliary basis and adds an exchange correction.", &
681 n_keywords=1, n_subsections=1, repeats=.false., &
682 citations=[guidon2010])
683
684 CALL keyword_create( &
685 keyword, __location__, &
686 name="ADMM_TYPE", &
687 description="Named ADMM variant from the literature. This shortcut sets METHOD, "// &
688 "ADMM_PURIFICATION_METHOD, and EXCH_SCALING_MODEL consistently for the selected variant.", &
689 enum_c_vals=s2a("NONE", "ADMM1", "ADMM2", "ADMMS", "ADMMP", "ADMMQ"), &
690 enum_desc=s2a("No short name is used, use specific definitions (default)", &
691 "ADMM1 method from Guidon2010", &
692 "ADMM2 method from Guidon2010", &
693 "ADMMS method from Merlot2014", &
694 "ADMMP method from Merlot2014", &
695 "ADMMQ method from Merlot2014"), &
697 default_i_val=no_admm_type, &
698 citations=[guidon2010, merlot2014])
699 CALL section_add_keyword(section, keyword)
700 CALL keyword_release(keyword)
701
702 CALL keyword_create( &
703 keyword, __location__, &
704 name="ADMM_PURIFICATION_METHOD", &
705 description="Method that shall be used for wavefunction fitting. Use MO_DIAG for MD.", &
706 enum_c_vals=s2a("NONE", "CAUCHY", "CAUCHY_SUBSPACE", "MO_DIAG", "MO_NO_DIAG", "MCWEENY", "NONE_DM"), &
710 enum_desc=s2a("Do not apply any purification", &
711 "Perform purification via general Cauchy representation", &
712 "Perform purification via Cauchy representation in occupied subspace", &
713 "Calculate MO derivatives via Cauchy representation by diagonalization", &
714 "Calculate MO derivatives via Cauchy representation by inversion", &
715 "Perform original McWeeny purification via matrix multiplications", &
716 "Do not apply any purification, works directly with density matrix"), &
717 default_i_val=do_admm_purify_mo_diag)
718 CALL section_add_keyword(section, keyword)
719 CALL keyword_release(keyword)
720
721 CALL keyword_create( &
722 keyword, __location__, &
723 name="METHOD", &
724 description="Method that shall be used for wavefunction fitting. Use BASIS_PROJECTION for MD.", &
725 enum_c_vals=s2a("BASIS_PROJECTION", "BLOCKED_PROJECTION_PURIFY_FULL", "BLOCKED_PROJECTION", &
726 "CHARGE_CONSTRAINED_PROJECTION"), &
729 enum_desc=s2a("Construct auxiliary density matrix from auxiliary basis.", &
730 "Construct auxiliary density from a blocked Fock matrix,"// &
731 " but use the original matrix for purification.", &
732 "Construct auxiliary density from a blocked Fock matrix.", &
733 "Construct auxiliary density from auxiliary basis enforcing charge constrain."), &
734 default_i_val=do_admm_basis_projection)
735 CALL section_add_keyword(section, keyword)
736 CALL keyword_release(keyword)
737
738 CALL keyword_create( &
739 keyword, __location__, &
740 name="EXCH_SCALING_MODEL", &
741 description="Scaling of the exchange correction calculated by the auxiliary density matrix.", &
742 enum_c_vals=s2a("NONE", "MERLOT"), &
744 enum_desc=s2a("No scaling is enabled, refers to methods ADMM1, ADMM2 or ADMMQ.", &
745 "Exchange scaling according to Merlot (2014)"), &
746 default_i_val=do_admm_exch_scaling_none)
747 CALL section_add_keyword(section, keyword)
748 CALL keyword_release(keyword)
749
750 CALL keyword_create( &
751 keyword, __location__, &
752 name="EXCH_CORRECTION_FUNC", &
753 description="Exchange functional used for the ADMM correction. It should be chosen consistently "// &
754 "with the exchange functional in the main XC setup. LibXC implementations require linking with LibXC.", &
755 enum_c_vals=s2a("DEFAULT", "PBEX", "NONE", "OPTX", "BECKE88X", &
756 "PBEX_LIBXC", "BECKE88X_LIBXC", "OPTX_LIBXC", "DEFAULT_LIBXC", "LDA_X_LIBXC"), &
762 enum_desc=s2a("Use PBE-based corrections according to the chosen interaction operator.", &
763 "Use PBEX functional for exchange correction.", &
764 "No correction: X(D)-x(d)-> 0.", &
765 "Use OPTX functional for exchange correction.", &
766 "Use Becke88X functional for exchange correction.", &
767 "Use PBEX functional (LibXC implementation) for exchange correction.", &
768 "Use Becke88X functional (LibXC implementation) for exchange correction.", &
769 "Use OPTX functional (LibXC implementation) for exchange correction.", &
770 "Use PBE-based corrections (LibXC where possible) to the chosen interaction operator.", &
771 "Use Slater X functional (LibXC where possible) for exchange correction."), &
772 default_i_val=do_admm_aux_exch_func_default)
773 CALL section_add_keyword(section, keyword)
774 CALL keyword_release(keyword)
775
776 CALL keyword_create(keyword, __location__, name="optx_a1", &
777 description="OPTX a1 coefficient", &
778 default_r_val=1.05151_dp)
779 CALL section_add_keyword(section, keyword)
780 CALL keyword_release(keyword)
781 CALL keyword_create(keyword, __location__, name="optx_a2", &
782 description="OPTX a2 coefficient", &
783 default_r_val=1.43169_dp)
784 CALL section_add_keyword(section, keyword)
785 CALL keyword_release(keyword)
786 CALL keyword_create(keyword, __location__, name="optx_gamma", &
787 description="OPTX gamma coefficient", &
788 default_r_val=0.006_dp)
789 CALL section_add_keyword(section, keyword)
790 CALL keyword_release(keyword)
791
792 CALL keyword_create(keyword, __location__, name="BLOCK_LIST", &
793 description="Specifies a list of atoms.", &
794 usage="BLOCK_LIST {integer} {integer} .. {integer}", &
795 n_var=-1, type_of_var=integer_t, repeats=.true.)
796 CALL section_add_keyword(section, keyword)
797 CALL keyword_release(keyword)
798
799 CALL keyword_create(keyword, __location__, name="EPS_FILTER", &
800 description="Define accuracy of DBCSR operations", &
801 usage="EPS_FILTER", default_r_val=0.0_dp)
802 CALL section_add_keyword(section, keyword)
803 CALL keyword_release(keyword)
804
805 END SUBROUTINE create_admm_section
806
807! **************************************************************************************************
808!> \brief ...
809!> \param section ...
810! **************************************************************************************************
811 SUBROUTINE create_density_fitting_section(section)
812 TYPE(section_type), POINTER :: section
813
814 TYPE(keyword_type), POINTER :: keyword
815 TYPE(section_type), POINTER :: print_key
816
817 NULLIFY (keyword, print_key)
818 cpassert(.NOT. ASSOCIATED(section))
819 CALL section_create(section, __location__, name="DENSITY_FITTING", &
820 description="Setup parameters for density fitting (Bloechl charges or density derived "// &
821 "atomic point charges (DDAPC) charges)", &
822 n_keywords=7, n_subsections=0, repeats=.false., &
823 citations=[blochl1995])
824
825 CALL keyword_create(keyword, __location__, name="NUM_GAUSS", &
826 description="Specifies the numbers of gaussian used to fit the QM density for each atomic site.", &
827 usage="NUM_GAUSS {integer}", &
828 n_var=1, type_of_var=integer_t, default_i_val=3)
829 CALL section_add_keyword(section, keyword)
830 CALL keyword_release(keyword)
831
832 CALL keyword_create(keyword, __location__, name="PFACTOR", &
833 description="Specifies the progression factor for the gaussian exponent for each atomic site.", &
834 usage="PFACTOR {real}", &
835 n_var=1, type_of_var=real_t, default_r_val=1.5_dp)
836 CALL section_add_keyword(section, keyword)
837 CALL keyword_release(keyword)
838
839 CALL keyword_create(keyword, __location__, name="MIN_RADIUS", &
840 description="Specifies the smallest radius of the gaussian used in the fit. All other radius are"// &
841 " obtained with the progression factor.", &
842 usage="MIN_RADIUS {real}", &
843 unit_str="angstrom", n_var=1, type_of_var=real_t, default_r_val=0.5_dp)
844 CALL section_add_keyword(section, keyword)
845 CALL keyword_release(keyword)
846
847 CALL keyword_create(keyword, __location__, name="RADII", &
848 description="Specifies all the radius of the gaussian used in the fit for each atomic site. The use"// &
849 " of this keyword disables all other keywords of this section.", &
850 usage="RADII {real} {real} .. {real}", &
851 unit_str="angstrom", n_var=-1, type_of_var=real_t)
852 CALL section_add_keyword(section, keyword)
853 CALL keyword_release(keyword)
854
855 CALL keyword_create(keyword, __location__, name="GCUT", &
856 description="Cutoff for charge fit in G-space.", &
857 usage="GCUT {real}", &
858 n_var=1, type_of_var=real_t, default_r_val=sqrt(6.0_dp))
859 CALL section_add_keyword(section, keyword)
860 CALL keyword_release(keyword)
861
862 CALL cp_print_key_section_create(print_key, __location__, "program_run_info", &
863 description="Controls the printing of basic information during the run", &
864 print_level=high_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
865
866 CALL keyword_create(keyword, __location__, name="CONDITION_NUMBER", &
867 description="Prints information regarding the condition numbers of the A matrix (to be inverted)", &
868 usage="CONDITION_NUMBER <LOGICAL>", &
869 default_l_val=.false., lone_keyword_l_val=.true.)
870 CALL section_add_keyword(print_key, keyword)
871 CALL keyword_release(keyword)
872
873 CALL section_add_subsection(section, print_key)
874 CALL section_release(print_key)
875
876 END SUBROUTINE create_density_fitting_section
877
878! **************************************************************************************************
879!> \brief creates the input section for the relativistic part
880!> \param section the section to create
881!> \author jens
882! **************************************************************************************************
883 SUBROUTINE create_relativistic_section(section)
884 TYPE(section_type), POINTER :: section
885
886 TYPE(keyword_type), POINTER :: keyword
887
888 cpassert(.NOT. ASSOCIATED(section))
889 CALL section_create(section, __location__, name="relativistic", &
890 description="parameters needed and setup for relativistic calculations", &
891 n_keywords=5, n_subsections=0, repeats=.false.)
892
893 NULLIFY (keyword)
894
895 CALL keyword_create(keyword, __location__, name="method", &
896 description="type of relativistic correction used", &
897 usage="method (NONE|DKH|ZORA)", default_i_val=rel_none, &
898 enum_c_vals=s2a("NONE", "DKH", "ZORA"), &
899 enum_i_vals=[rel_none, rel_dkh, rel_zora], &
900 enum_desc=s2a("Use no relativistic correction", &
901 "Use Douglas-Kroll-Hess method", &
902 "Use ZORA method"))
903 CALL section_add_keyword(section, keyword)
904 CALL keyword_release(keyword)
905
906 CALL keyword_create(keyword, __location__, name="DKH_order", &
907 description="The order of the DKH transformation ", &
908 usage="DKH_order 2", default_i_val=2)
909 CALL section_add_keyword(section, keyword)
910 CALL keyword_release(keyword)
911
912 CALL keyword_create(keyword, __location__, name="ZORA_type", &
913 description="Type of ZORA method to be used", &
914 usage="ZORA_type scMP", default_i_val=rel_zora_full, &
915 enum_c_vals=s2a("FULL", "MP", "scMP"), &
916 enum_desc=s2a("Full ZORA method (not implemented)", &
917 "ZORA with atomic model potential", &
918 "Scaled ZORA with atomic model potential"), &
920 CALL section_add_keyword(section, keyword)
921 CALL keyword_release(keyword)
922
923 CALL keyword_create(keyword, __location__, name="transformation", &
924 description="Type of DKH transformation", &
925 usage="transformation (FULL|MOLECULE|ATOM)", default_i_val=rel_trans_atom, &
926 enum_c_vals=s2a("FULL", "MOLECULE", "ATOM"), &
928 enum_desc=s2a("Use full matrix transformation", &
929 "Use transformation blocked by molecule", &
930 "Use atomic blocks"))
931 CALL section_add_keyword(section, keyword)
932 CALL keyword_release(keyword)
933
934 CALL keyword_create(keyword, __location__, name="z_cutoff", &
935 description="The minimal atomic number considered for atom transformation", &
936 usage="z_cutoff 50", default_i_val=1)
937 CALL section_add_keyword(section, keyword)
938 CALL keyword_release(keyword)
939
940 CALL keyword_create(keyword, __location__, name="potential", &
941 description="External potential used in DKH transformation, full 1/r or erfc(r)/r", &
942 usage="POTENTIAL {FULL,ERFC}", default_i_val=rel_pot_erfc, &
943 enum_c_vals=s2a("FULL", "ERFC"), &
944 enum_i_vals=[rel_pot_full, rel_pot_erfc])
945 CALL section_add_keyword(section, keyword)
946 CALL keyword_release(keyword)
947
948 END SUBROUTINE create_relativistic_section
949
950! **************************************************************************************************
951!> \brief creates the KG section
952!> \param section ...
953!> \author Martin Haeufel [2012.07]
954! **************************************************************************************************
955 SUBROUTINE create_kg_section(section)
956 TYPE(section_type), POINTER :: section
957
958 TYPE(keyword_type), POINTER :: keyword
959 TYPE(section_type), POINTER :: print_key, subsection
960
961 cpassert(.NOT. ASSOCIATED(section))
962 CALL section_create(section, __location__, name="KG_METHOD", &
963 description="Specifies the parameters for a Kim-Gordon-like partitioning"// &
964 " into molecular subunits", &
965 n_keywords=0, n_subsections=1, repeats=.false., &
967
968 NULLIFY (keyword, subsection, print_key)
969
970 ! add a XC section
971 CALL create_xc_section(subsection)
972 CALL section_add_subsection(section, subsection)
973 CALL section_release(subsection)
974
975 ! add LRI section
976 CALL create_lrigpw_section(subsection)
977 CALL section_add_subsection(section, subsection)
978 CALL section_release(subsection)
979
980 CALL keyword_create(keyword, __location__, name="COLORING_METHOD", &
981 description="Which algorithm to use for coloring.", &
982 usage="COLORING_METHOD GREEDY", &
983 default_i_val=kg_color_dsatur, &
984 enum_c_vals=s2a("DSATUR", "GREEDY"), &
985 enum_desc=s2a("Maximum degree of saturation, relatively accurate", &
986 "Greedy, fast coloring, less accurate"), &
987 enum_i_vals=[kg_color_dsatur, kg_color_greedy])
988 CALL section_add_keyword(section, keyword)
989 CALL keyword_release(keyword)
990
991 CALL keyword_create(keyword, __location__, name="TNADD_METHOD", &
992 description="Algorithm to use for the calculation of the nonadditive kinetic energy.", &
993 usage="TNADD_METHOD ATOMIC", &
994 default_i_val=kg_tnadd_embed, &
995 enum_c_vals=s2a("EMBEDDING", "RI_EMBEDDING", "ATOMIC", "NONE"), &
996 enum_desc=s2a("Use full embedding potential (see Iannuzzi et al)", &
997 "Use full embedding potential with RI density fitting", &
998 "Use sum of atomic model potentials", &
999 "Do not use kinetic energy embedding"), &
1001 CALL section_add_keyword(section, keyword)
1002 CALL keyword_release(keyword)
1003
1004 CALL keyword_create(keyword, __location__, name="INTEGRATION_GRID", &
1005 description="Grid [small,medium,large,huge]to be used for the TNADD integration.", &
1006 usage="INTEGRATION_GRID MEDIUM", &
1007 default_c_val="MEDIUM")
1008 CALL section_add_keyword(section, keyword)
1009 CALL keyword_release(keyword)
1010
1011 CALL section_create(subsection, __location__, name="PRINT", &
1012 description="Print section", &
1013 n_keywords=0, n_subsections=1, repeats=.false.)
1014
1015 CALL cp_print_key_section_create(print_key, __location__, "NEIGHBOR_LISTS", &
1016 description="Controls the printing of the neighbor lists.", &
1017 print_level=low_print_level, filename="__STD_OUT__", unit_str="angstrom")
1018
1019 CALL keyword_create(keyword, __location__, &
1020 name="SAB_ORB_FULL", &
1021 description="Activates the printing of the full orbital "// &
1022 "orbital neighbor lists.", &
1023 default_l_val=.false., &
1024 lone_keyword_l_val=.true.)
1025 CALL section_add_keyword(print_key, keyword)
1026 CALL keyword_release(keyword)
1027
1028 CALL keyword_create(keyword, __location__, &
1029 name="SAB_ORB_MOLECULAR", &
1030 description="Activates the printing of the orbital "// &
1031 "orbital neighbor lists for molecular subsets.", &
1032 default_l_val=.false., &
1033 lone_keyword_l_val=.true.)
1034 CALL section_add_keyword(print_key, keyword)
1035 CALL keyword_release(keyword)
1036
1037 CALL keyword_create(keyword, __location__, &
1038 name="SAC_KIN", &
1039 description="Activates the printing of the orbital "// &
1040 "atomic potential neighbor list.", &
1041 default_l_val=.false., &
1042 lone_keyword_l_val=.true.)
1043 CALL section_add_keyword(print_key, keyword)
1044 CALL keyword_release(keyword)
1045
1046 CALL section_add_subsection(subsection, print_key)
1047 CALL section_release(print_key)
1048
1049 CALL section_add_subsection(section, subsection)
1050 CALL section_release(subsection)
1051
1052 END SUBROUTINE create_kg_section
1053
1054! **************************************************************************************************
1055!> \brief Create the BSSE section for counterpoise correction
1056!> \param section the section to create
1057!> \author teo
1058! **************************************************************************************************
1059 SUBROUTINE create_bsse_section(section)
1060 TYPE(section_type), POINTER :: section
1061
1062 TYPE(keyword_type), POINTER :: keyword
1063 TYPE(section_type), POINTER :: subsection
1064
1065 cpassert(.NOT. ASSOCIATED(section))
1066 CALL section_create(section, __location__, name="BSSE", &
1067 description="This section is used to set up the BSSE calculation. "// &
1068 "It also requires that for each atomic kind X a kind X_ghost is present, "// &
1069 "with the GHOST keyword specified, in addition to the other required fields.", &
1070 n_keywords=3, n_subsections=1, repeats=.false.)
1071
1072 NULLIFY (keyword, subsection)
1073 ! FRAGMENT SECTION
1074 CALL section_create(subsection, __location__, name="FRAGMENT", &
1075 description="Specify the atom number belonging to this fragment.", &
1076 n_keywords=2, n_subsections=0, repeats=.true.)
1077
1078 CALL keyword_create(keyword, __location__, name="LIST", &
1079 description="Specifies a list of atoms.", &
1080 usage="LIST {integer} {integer} .. {integer}", &
1081 repeats=.true., n_var=-1, type_of_var=integer_t)
1082 CALL section_add_keyword(subsection, keyword)
1083 CALL keyword_release(keyword)
1084
1085 CALL section_add_subsection(section, subsection)
1086 CALL section_release(subsection)
1087
1088 ! CONFIGURATION SECTION
1089 CALL section_create(subsection, __location__, name="CONFIGURATION", &
1090 description="Specify additional parameters for the combinatorial configurations. "// &
1091 "Use this section to manually specify charge and multiplicity of the fragments "// &
1092 "and their combinations.", &
1093 n_keywords=2, n_subsections=0, repeats=.true.)
1094
1095 CALL keyword_create(keyword, __location__, name="GLB_CONF", &
1096 description="Specifies the global configuration using 1 or 0 for each fragment. "// &
1097 "1 specifies the respective fragment as used, 0 as unused.", &
1098 usage="GLB_CONF {integer} {integer} .. {integer}", &
1099 n_var=-1, type_of_var=integer_t)
1100 CALL section_add_keyword(subsection, keyword)
1101 CALL keyword_release(keyword)
1102
1103 CALL keyword_create(keyword, __location__, name="SUB_CONF", &
1104 description="Specifies the subconfiguration using 1 or 0 belonging to the global configuration. "// &
1105 "1 specifies the respective fragment as real, 0 as ghost.", &
1106 usage="SUB_CONF {integer} {integer} .. {integer}", &
1107 n_var=-1, type_of_var=integer_t)
1108 CALL section_add_keyword(subsection, keyword)
1109 CALL keyword_release(keyword)
1110
1111 CALL keyword_create(keyword, __location__, &
1112 name="MULTIPLICITY", &
1113 variants=["MULTIP"], &
1114 description="Specify for each fragment the multiplicity. Two times the total spin plus one. "// &
1115 "Specify 3 for a triplet, 4 for a quartet,and so on. Default is 1 (singlet) for an "// &
1116 "even number and 2 (doublet) for an odd number of electrons.", &
1117 usage="MULTIPLICITY 3", &
1118 default_i_val=0) ! this default value is just a flag to get the above
1119 CALL section_add_keyword(subsection, keyword)
1120 CALL keyword_release(keyword)
1121
1122 CALL keyword_create(keyword, __location__, name="CHARGE", &
1123 description="The total charge for each fragment.", &
1124 usage="CHARGE -1", &
1125 default_i_val=0)
1126 CALL section_add_keyword(subsection, keyword)
1127 CALL keyword_release(keyword)
1128 CALL section_add_subsection(section, subsection)
1129 CALL section_release(subsection)
1130
1131 CALL section_create(subsection, __location__, name="FRAGMENT_ENERGIES", &
1132 description="This section contains the energies of the fragments already"// &
1133 " computed. It is useful as a summary and specifically for restarting BSSE runs.", &
1134 n_keywords=2, n_subsections=0, repeats=.true.)
1135 CALL keyword_create(keyword, __location__, name="_DEFAULT_KEYWORD_", &
1136 description="The energy computed for each fragment", repeats=.true., &
1137 usage="{REAL}", type_of_var=real_t)
1138 CALL section_add_keyword(subsection, keyword)
1139 CALL keyword_release(keyword)
1140 CALL section_add_subsection(section, subsection)
1141 CALL section_release(subsection)
1142
1143 CALL create_print_bsse_section(subsection)
1144 CALL section_add_subsection(section, subsection)
1145 CALL section_release(subsection)
1146
1147 END SUBROUTINE create_bsse_section
1148
1149! **************************************************************************************************
1150!> \brief Create the print bsse section
1151!> \param section the section to create
1152!> \author teo
1153! **************************************************************************************************
1154 SUBROUTINE create_print_bsse_section(section)
1155 TYPE(section_type), POINTER :: section
1156
1157 TYPE(section_type), POINTER :: print_key
1158
1159 cpassert(.NOT. ASSOCIATED(section))
1160 CALL section_create(section, __location__, name="print", &
1161 description="Section of possible print options in BSSE code.", &
1162 n_keywords=0, n_subsections=1, repeats=.false.)
1163
1164 NULLIFY (print_key)
1165 CALL cp_print_key_section_create(print_key, __location__, "PROGRAM_RUN_INFO", &
1166 description="Controls the printing of information regarding the run.", &
1167 print_level=low_print_level, filename="__STD_OUT__")
1168 CALL section_add_subsection(section, print_key)
1169 CALL section_release(print_key)
1170
1171 CALL cp_print_key_section_create(print_key, __location__, "RESTART", &
1172 description="Controls the dumping of the restart file during BSSE runs. "// &
1173 "By default the restart is updated after each configuration calculation is "// &
1174 "completed.", &
1175 print_level=silent_print_level, common_iter_levels=0, &
1176 add_last=add_last_numeric, filename="")
1177 CALL section_add_subsection(section, print_key)
1178 CALL section_release(print_key)
1179
1180 END SUBROUTINE create_print_bsse_section
1181
1182! **************************************************************************************************
1183!> \brief input section for optional parameters for RIGPW
1184!> \param section the section to create
1185!> \author JGH [06.2017]
1186! **************************************************************************************************
1187 SUBROUTINE create_rigpw_section(section)
1188 TYPE(section_type), POINTER :: section
1189
1190 cpassert(.NOT. ASSOCIATED(section))
1191 CALL section_create(section, __location__, name="RIGPW", &
1192 description="This section specifies optional parameters for RIGPW.", &
1193 n_keywords=1, n_subsections=0, repeats=.false.)
1194
1195! CALL keyword_create(keyword, __LOCATION__, name="RI_OVERLAP_MATRIX", &
1196! description="Specifies whether to calculate the inverse or the "// &
1197! "pseudoinverse of the overlap matrix of the auxiliary "// &
1198! "basis set. Calculating the pseudoinverse is necessary "// &
1199! "for very large auxiliary basis sets, but more expensive. "// &
1200! "Using the pseudoinverse, consistent forces are not "// &
1201! "guaranteed yet.", &
1202! usage="RI_OVERLAP_MATRIX INVERSE", &
1203! enum_c_vals=s2a("INVERSE", "PSEUDO_INVERSE_SVD", "PSEUDO_INVERSE_DIAG", &
1204! "AUTOSELECT"), &
1205! enum_desc=s2a("Calculate inverse of the overlap matrix.", &
1206! "Calculate the pseuodinverse of the overlap matrix "// &
1207! "using singular value decomposition.", &
1208! "Calculate the pseudoinverse of the overlap matrix "// &
1209! "by prior diagonalization.", &
1210! "Choose automatically for each pair whether to "// &
1211! "calculate the inverse or pseudoinverse based on the "// &
1212! "condition number of the overlap matrix for each pair. "// &
1213! "Calculating the pseudoinverse is much more expensive."), &
1214! enum_i_vals=(/do_lri_inv, do_lri_pseudoinv_svd, &
1215! do_lri_pseudoinv_diag, do_lri_inv_auto/), &
1216! default_i_val=do_lri_inv)
1217! CALL section_add_keyword(section, keyword)
1218! CALL keyword_release(keyword)
1219
1220 END SUBROUTINE create_rigpw_section
1221
1222! **************************************************************************************************
1223!> \brief creates the multigrid
1224!> \param section input section to create
1225!> \param create_subsections indicates whether or not subsections INTERPOLATOR and RS_GRID
1226!> should be created
1227!> \author fawzi
1228! **************************************************************************************************
1229 SUBROUTINE create_mgrid_section(section, create_subsections)
1230 TYPE(section_type), POINTER :: section
1231 LOGICAL, INTENT(in) :: create_subsections
1232
1233 TYPE(keyword_type), POINTER :: keyword
1234 TYPE(section_type), POINTER :: subsection
1235
1236 cpassert(.NOT. ASSOCIATED(section))
1237 CALL section_create(section, __location__, name="mgrid", &
1238 description="Controls the multigrid used by GPW/GAPW to represent densities, "// &
1239 "potentials, and Gaussian products on real-space grids.", &
1240 n_keywords=5, n_subsections=1, repeats=.false.)
1241 NULLIFY (keyword)
1242 CALL keyword_create(keyword, __location__, name="NGRIDS", &
1243 description="Number of multigrid levels. Smooth Gaussian products can be mapped to "// &
1244 "coarser levels, while sharper products require finer levels.", &
1245 usage="ngrids 1", default_i_val=4)
1246 CALL section_add_keyword(section, keyword)
1247 CALL keyword_release(keyword)
1248
1249 CALL keyword_create(keyword, __location__, name="cutoff", &
1250 description= &
1251 "Plane-wave cutoff of the finest real-space grid level. "// &
1252 "Increasing this value improves the grid representation, but it is "// &
1253 "not a substitute for converging the Gaussian basis set. "// &
1254 "Default value for SE or DFTB calculation is 1.0 [Ry].", &
1255 usage="cutoff 300", &
1256 default_r_val=cp_unit_to_cp2k(value=280.0_dp, unit_str="Ry"), &
1257 n_var=1, unit_str="Ry")
1258 CALL section_add_keyword(section, keyword)
1259 CALL keyword_release(keyword)
1260
1261 CALL keyword_create(keyword, __location__, name="progression_factor", &
1262 description="Factor used to derive the cutoff of coarser multigrid levels when "// &
1263 "they are not given explicitly.", &
1264 usage="progression_factor <integer>", default_r_val=3._dp)
1265 CALL section_add_keyword(section, keyword)
1266 CALL keyword_release(keyword)
1267
1268 CALL keyword_create(keyword, __location__, name="commensurate", &
1269 description="If the grids should be commensurate. If true overrides "// &
1270 "the progression factor and the cutoffs of the sub grids", &
1271 usage="commensurate", default_l_val=.false., &
1272 lone_keyword_l_val=.true.)
1273 CALL section_add_keyword(section, keyword)
1274 CALL keyword_release(keyword)
1275
1276 CALL keyword_create(keyword, __location__, name="realspace", &
1277 description="If both rho and rho_gspace are needed ", &
1278 usage="realspace", default_l_val=.false., &
1279 lone_keyword_l_val=.true.)
1280 CALL section_add_keyword(section, keyword)
1281 CALL keyword_release(keyword)
1282
1283 CALL keyword_create(keyword, __location__, name="REL_CUTOFF", &
1284 variants=["RELATIVE_CUTOFF"], &
1285 description="Controls to which multigrid level a Gaussian product is mapped. "// &
1286 "It is the reference cutoff for a Gaussian with exponent alpha=1. Larger values "// &
1287 "keep more Gaussian products on finer grids and can be important for accurate "// &
1288 "energies, forces, stress tensors, and variable-cell simulations.", &
1289 usage="RELATIVE_CUTOFF real", default_r_val=20.0_dp, &
1290 unit_str="Ry")
1291 CALL section_add_keyword(section, keyword)
1292 CALL keyword_release(keyword)
1293
1294 CALL keyword_create(keyword, __location__, name="MULTIGRID_SET", &
1295 description="Activate a manual setting of the multigrids", &
1296 usage="MULTIGRID_SET", default_l_val=.false.)
1297 CALL section_add_keyword(section, keyword)
1298 CALL keyword_release(keyword)
1299
1300 CALL keyword_create(keyword, __location__, &
1301 name="SKIP_LOAD_BALANCE_DISTRIBUTED", &
1302 description="Skips load balancing on distributed multigrids. "// &
1303 "Memory usage is O(p) so may be used "// &
1304 "for all but the very largest runs.", &
1305 usage="SKIP_LOAD_BALANCE_DISTRIBUTED", &
1306 default_l_val=.false., &
1307 lone_keyword_l_val=.true.)
1308! CALL keyword_create(keyword, __LOCATION__, name="SKIP_LOAD_BALANCE_DISTRIBUTED",&
1309! description="Skip load balancing on distributed multigrids, which might be memory intensive."//&
1310! "If not explicitly specified, runs using more than 1024 MPI tasks will default to .TRUE.",&
1311! usage="SKIP_LOAD_BALANCE_DISTRIBUTED", default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1312
1313 CALL section_add_keyword(section, keyword)
1314 CALL keyword_release(keyword)
1315
1316 CALL keyword_create(keyword, __location__, name="MULTIGRID_CUTOFF", &
1317 variants=["CUTOFF_LIST"], &
1318 description="List of cutoff values to set up multigrids manually", &
1319 usage="MULTIGRID_CUTOFF 200.0 100.0 ", &
1320 n_var=-1, &
1321 type_of_var=real_t, &
1322 unit_str="Ry")
1323 CALL section_add_keyword(section, keyword)
1324 CALL keyword_release(keyword)
1325
1326 IF (create_subsections) THEN
1327 NULLIFY (subsection)
1328 CALL create_rsgrid_section(subsection)
1329 CALL section_add_subsection(section, subsection)
1330 CALL section_release(subsection)
1331
1332 NULLIFY (subsection)
1333 CALL create_interp_section(subsection)
1334 CALL section_add_subsection(section, subsection)
1335 CALL section_release(subsection)
1336 END IF
1337 END SUBROUTINE create_mgrid_section
1338
1339! **************************************************************************************************
1340!> \brief creates the interpolation section
1341!> \param section ...
1342!> \author tlaino
1343! **************************************************************************************************
1344 SUBROUTINE create_interp_section(section)
1345 TYPE(section_type), POINTER :: section
1346
1347 TYPE(keyword_type), POINTER :: keyword
1348 TYPE(section_type), POINTER :: print_key
1349
1350 cpassert(.NOT. ASSOCIATED(section))
1351 CALL section_create(section, __location__, name="interpolator", &
1352 description="kind of interpolation used between the multigrids", &
1353 n_keywords=5, n_subsections=0, repeats=.false.)
1354
1355 NULLIFY (keyword, print_key)
1356
1357 CALL keyword_create(keyword, __location__, name="kind", &
1358 description="the interpolator to use", &
1359 usage="kind spline3", &
1360 default_i_val=pw_interp, &
1361 enum_c_vals=s2a("pw", "spline3_nopbc", "spline3"), &
1362 enum_i_vals=[pw_interp, &
1364 CALL section_add_keyword(section, keyword)
1365 CALL keyword_release(keyword)
1366
1367 CALL keyword_create(keyword, __location__, name="safe_computation", &
1368 description="if a non unrolled calculation is to be performed in parallel", &
1369 usage="safe_computation OFF", &
1370 default_l_val=.false., &
1371 lone_keyword_l_val=.true.)
1372 CALL section_add_keyword(section, keyword)
1373 CALL keyword_release(keyword)
1374
1375 CALL keyword_create(keyword, __location__, name="aint_precond", &
1376 description="the approximate inverse to use to get the starting point"// &
1377 " for the linear solver of the spline3 methods", &
1378 usage="aint_precond copy", &
1379 default_i_val=precond_spl3_aint, &
1380 enum_c_vals=s2a("copy", "spl3_nopbc_aint1", "spl3_nopbc_aint2", &
1381 "spl3_nopbc_precond1", "spl3_nopbc_precond2", "spl3_nopbc_precond3"), &
1384 CALL section_add_keyword(section, keyword)
1385 CALL keyword_release(keyword)
1386
1387 CALL keyword_create(keyword, __location__, name="precond", &
1388 description="The preconditioner used"// &
1389 " for the linear solver of the spline3 methods", &
1390 usage="PRECOND copy", &
1391 default_i_val=precond_spl3_3, &
1392 enum_c_vals=s2a("copy", "spl3_nopbc_aint1", "spl3_nopbc_aint2", &
1393 "spl3_nopbc_precond1", "spl3_nopbc_precond2", "spl3_nopbc_precond3"), &
1396 CALL section_add_keyword(section, keyword)
1397 CALL keyword_release(keyword)
1398
1399 CALL keyword_create(keyword, __location__, name="eps_x", &
1400 description="accuracy on the solution for spline3 the interpolators", &
1401 usage="eps_x 1.e-15", default_r_val=1.e-10_dp)
1402 CALL section_add_keyword(section, keyword)
1403 CALL keyword_release(keyword)
1404
1405 CALL keyword_create(keyword, __location__, name="eps_r", &
1406 description="accuracy on the residual for spline3 the interpolators", &
1407 usage="eps_r 1.e-15", default_r_val=1.e-10_dp)
1408 CALL section_add_keyword(section, keyword)
1409 CALL keyword_release(keyword)
1410
1411 CALL keyword_create(keyword, __location__, name="max_iter", &
1412 variants=['maxiter'], &
1413 description="the maximum number of iterations", &
1414 usage="max_iter 200", default_i_val=100)
1415 CALL section_add_keyword(section, keyword)
1416 CALL keyword_release(keyword)
1417
1418 NULLIFY (print_key)
1419 CALL cp_print_key_section_create(print_key, __location__, "conv_info", &
1420 description="if convergence information about the linear solver"// &
1421 " of the spline methods should be printed", &
1422 print_level=medium_print_level, each_iter_names=s2a("SPLINE_FIND_COEFFS"), &
1423 each_iter_values=[10], filename="__STD_OUT__", &
1424 add_last=add_last_numeric)
1425 CALL section_add_subsection(section, print_key)
1426 CALL section_release(print_key)
1427
1428 END SUBROUTINE create_interp_section
1429
1430! **************************************************************************************************
1431!> \brief creates the sic (self interaction correction) section
1432!> \param section ...
1433!> \author fawzi
1434! **************************************************************************************************
1435 SUBROUTINE create_sic_section(section)
1436 TYPE(section_type), POINTER :: section
1437
1438 TYPE(keyword_type), POINTER :: keyword
1439
1440 cpassert(.NOT. ASSOCIATED(section))
1441 CALL section_create(section, __location__, name="sic", &
1442 description="parameters for the self interaction correction", &
1443 n_keywords=6, n_subsections=0, repeats=.false., &
1445
1446 NULLIFY (keyword)
1447
1448 CALL keyword_create(keyword, __location__, name="SIC_SCALING_A", &
1449 description="Scaling of the coulomb term in sic [experimental]", &
1450 usage="SIC_SCALING_A 0.5", &
1451 citations=[vandevondele2005b], &
1452 default_r_val=1.0_dp)
1453 CALL section_add_keyword(section, keyword)
1454 CALL keyword_release(keyword)
1455
1456 CALL keyword_create(keyword, __location__, name="SIC_SCALING_B", &
1457 description="Scaling of the xc term in sic [experimental]", &
1458 usage="SIC_SCALING_B 0.5", &
1459 citations=[vandevondele2005b], &
1460 default_r_val=1.0_dp)
1461 CALL section_add_keyword(section, keyword)
1462 CALL keyword_release(keyword)
1463
1464 CALL keyword_create(keyword, __location__, name="SIC_METHOD", &
1465 description="Method used to remove the self interaction", &
1466 usage="SIC_METHOD MAURI_US", &
1467 default_i_val=sic_none, &
1468 enum_c_vals=s2a("NONE", "MAURI_US", "MAURI_SPZ", "AD", "EXPLICIT_ORBITALS"), &
1469 enum_i_vals=[sic_none, sic_mauri_us, sic_mauri_spz, sic_ad, sic_eo], &
1470 enum_desc=s2a("Do not apply a sic correction", &
1471 "Employ a (scaled) correction proposed by Mauri and co-workers"// &
1472 " on the spin density / doublet unpaired orbital", &
1473 "Employ a (scaled) Perdew-Zunger expression"// &
1474 " on the spin density / doublet unpaired orbital", &
1475 "The average density correction", &
1476 "(scaled) Perdew-Zunger correction explicitly on a set of orbitals."), &
1478 CALL section_add_keyword(section, keyword)
1479 CALL keyword_release(keyword)
1480
1481 CALL keyword_create(keyword, __location__, name="ORBITAL_SET", &
1482 description="Type of orbitals treated with the SIC", &
1483 usage="ORBITAL_SET ALL", &
1484 default_i_val=sic_list_unpaired, &
1485 enum_c_vals=s2a("UNPAIRED", "ALL"), &
1486 enum_desc=s2a("correction for the unpaired orbitals only, requires a restricted open shell calculation", &
1487 "correction for all orbitals, requires a LSD or ROKS calculation"), &
1488 enum_i_vals=[sic_list_unpaired, sic_list_all])
1489 CALL section_add_keyword(section, keyword)
1490 CALL keyword_release(keyword)
1491
1492 END SUBROUTINE create_sic_section
1493
1494! **************************************************************************************************
1495!> \brief creates the low spin roks section
1496!> \param section ...
1497!> \author Joost VandeVondele
1498! **************************************************************************************************
1499 SUBROUTINE create_low_spin_roks_section(section)
1500 TYPE(section_type), POINTER :: section
1501
1502 TYPE(keyword_type), POINTER :: keyword
1503
1504 cpassert(.NOT. ASSOCIATED(section))
1505 CALL section_create(section, __location__, name="LOW_SPIN_ROKS", &
1506 description="Specify the details of the low spin ROKS method. "// &
1507 "In particular, one can specify various terms added to the energy of the high spin roks configuration"// &
1508 " with a energy scaling factor, and a prescription of the spin state.", &
1509 n_keywords=6, n_subsections=0, repeats=.false.)
1510
1511 NULLIFY (keyword)
1512 CALL keyword_create(keyword, __location__, name="ENERGY_SCALING", &
1513 description="The scaling factors for each term added to the total energy. "// &
1514 "This list should contain one number for each term added to the total energy.", &
1515 usage="ENERGY_SCALING 1.0 -1.0 ", &
1516 n_var=-1, type_of_var=real_t, repeats=.false.)
1517 CALL section_add_keyword(section, keyword)
1518 CALL keyword_release(keyword)
1519 CALL keyword_create( &
1520 keyword, __location__, name="SPIN_CONFIGURATION", &
1521 description="For each singly occupied orbital, specify if this should be an alpha (=1) or a beta (=2) orbital. "// &
1522 "This keyword should be repeated, each repetition corresponding to an additional term.", &
1523 usage="SPIN_CONFIGURATION 1 2", &
1524 n_var=-1, type_of_var=integer_t, repeats=.true.)
1525 CALL section_add_keyword(section, keyword)
1526 CALL keyword_release(keyword)
1527
1528 END SUBROUTINE create_low_spin_roks_section
1529
1530! **************************************************************************************************
1531!> \brief ...
1532!> \param section ...
1533! **************************************************************************************************
1534 SUBROUTINE create_rtp_section(section)
1535 TYPE(section_type), POINTER :: section
1536
1537 TYPE(keyword_type), POINTER :: keyword
1538 TYPE(section_type), POINTER :: print_key, print_section, subsection
1539
1540 NULLIFY (keyword)
1541 cpassert(.NOT. ASSOCIATED(section))
1542 CALL section_create(section, __location__, name="REAL_TIME_PROPAGATION", &
1543 description="Parameters needed to set up the real time propagation"// &
1544 " for the electron dynamics. This currently works only in the NVE ensemble.", &
1545 n_keywords=4, n_subsections=4, repeats=.false., &
1546 citations=[kunert2003, andermatt2016])
1547
1548 CALL keyword_create(keyword, __location__, name="MAX_ITER", &
1549 description="Maximal number of iterations for the self consistent propagator loop.", &
1550 usage="MAX_ITER 10", &
1551 default_i_val=10)
1552 CALL section_add_keyword(section, keyword)
1553 CALL keyword_release(keyword)
1554
1555 CALL keyword_create(keyword, __location__, name="EPS_ITER", &
1556 description="Convergence criterion for the self consistent propagator loop.", &
1557 usage="EPS_ITER 1.0E-5", &
1558 default_r_val=1.0e-7_dp)
1559 CALL section_add_keyword(section, keyword)
1560 CALL keyword_release(keyword)
1561
1562 CALL keyword_create(keyword, __location__, name="ASPC_ORDER", &
1563 description="Speciefies how many steps will be used for extrapolation. "// &
1564 "One will be always used which is means X(t+dt)=X(t)", &
1565 usage="ASPC_ORDER 3", &
1566 default_i_val=3)
1567 CALL section_add_keyword(section, keyword)
1568 CALL keyword_release(keyword)
1569
1570 CALL keyword_create(keyword, __location__, name="MAT_EXP", &
1571 description="Which method should be used to calculate the exponential"// &
1572 " in the propagator. It is recommended to use BCH when employing density_propagation "// &
1573 "and ARNOLDI otherwise.", &
1574 usage="MAT_EXP TAYLOR", default_i_val=do_arnoldi, &
1575 enum_c_vals=s2a("TAYLOR", "PADE", "ARNOLDI", "BCH", "EXACT"), &
1576 enum_i_vals=[do_taylor, do_pade, do_arnoldi, do_bch, do_exact], &
1577 enum_desc=s2a("exponential is evaluated using scaling and squaring in combination"// &
1578 " with a taylor expansion of the exponential.", &
1579 "uses scaling and squaring together with the pade approximation", &
1580 "uses arnoldi subspace algorithm to compute exp(H)*MO directly, can't be used in "// &
1581 "combination with Crank Nicholson or density propagation", &
1582 "Uses a Baker-Campbell-Hausdorff expansion to propagate the density matrix,"// &
1583 " only works for density propagation", &
1584 "Uses diagonalisation of the exponent matrices to determine the "// &
1585 "matrix exponential exactly. Only implemented for GWBSE."))
1586 CALL section_add_keyword(section, keyword)
1587 CALL keyword_release(keyword)
1588
1589 CALL keyword_create(keyword, __location__, name="DENSITY_PROPAGATION", &
1590 description="The density matrix is propagated instead of the molecular orbitals. "// &
1591 "This can allow a linear scaling simulation. The density matrix is filtered with "// &
1592 "the threshold based on the EPS_FILTER keyword from the LS_SCF section", &
1593 usage="DENSITY_PROPAGATION .TRUE.", &
1594 default_l_val=.false., lone_keyword_l_val=.true.)
1595 CALL section_add_keyword(section, keyword)
1596 CALL keyword_release(keyword)
1597
1598 CALL keyword_create(keyword, __location__, name="SC_CHECK_START", &
1599 description="Speciefies how many iteration steps will be done without "// &
1600 "a check for self consistency. Can save some time in big calculations.", &
1601 usage="SC_CHECK_START 3", &
1602 default_i_val=0)
1603 CALL section_add_keyword(section, keyword)
1604 CALL keyword_release(keyword)
1605
1606 CALL keyword_create(keyword, __location__, name="EXP_ACCURACY", &
1607 description="Accuracy for the taylor and pade approximation. "// &
1608 "This is only an upper bound bound since the norm used for the guess "// &
1609 "is an upper bound for the needed one.", &
1610 usage="EXP_ACCURACY 1.0E-6", &
1611 default_r_val=1.0e-9_dp)
1612 CALL section_add_keyword(section, keyword)
1613 CALL keyword_release(keyword)
1614
1615 CALL keyword_create(keyword, __location__, name="PROPAGATOR", &
1616 description="Which propagator should be used for the orbitals", &
1617 usage="PROPAGATOR ETRS", default_i_val=do_etrs, &
1618 enum_c_vals=s2a("ETRS", "CN", "EM"), &
1619 enum_i_vals=[do_etrs, do_cn, do_em], &
1620 enum_desc=s2a("enforced time reversible symmetry", &
1621 "Crank Nicholson propagator", &
1622 "Exponential midpoint propagator"))
1623 CALL section_add_keyword(section, keyword)
1624 CALL keyword_release(keyword)
1625
1626 CALL keyword_create(keyword, __location__, name="INITIAL_WFN", &
1627 description="Controls the initial WFN used for propagation. "// &
1628 "Note that some energy contributions may not be "// &
1629 "initialized in the restart cases, for instance "// &
1630 "electronic entropy energy in the case of smearing.", &
1631 usage="INITIAL_WFN SCF_WFN", default_i_val=use_scf_wfn, &
1632 enum_c_vals=s2a("SCF_WFN", "RESTART_WFN", "RT_RESTART"), &
1633 enum_i_vals=[use_scf_wfn, use_restart_wfn, use_rt_restart], &
1634 enum_desc=s2a("An SCF run is performed to get the initial state.", &
1635 "A wavefunction from a previous SCF is propagated. Especially useful,"// &
1636 " if electronic constraints or restraints are used in the previous calculation, "// &
1637 "since these do not work in the rtp scheme.", &
1638 "use the wavefunction of a real time propagation/ehrenfest run"))
1639 CALL section_add_keyword(section, keyword)
1640 CALL keyword_release(keyword)
1641
1642 CALL keyword_create(keyword, __location__, name="APPLY_WFN_MIX_INIT_RESTART", &
1643 description="If set to True and in the case of INITIAL_WFN=RESTART_WFN, call the "// &
1644 "DFT%PRINT%WFN_MIX section to mix the read initial wfn. The starting wave-function of the "// &
1645 "RTP will be the mixed one. Setting this to True without a defined WFN_MIX section will "// &
1646 "not do anything as defining a WFN_MIX section without this keyword for RTP run with "// &
1647 "INITIAL_WFN=RESTART_WFN. Note that if INITIAL_WFN=SCF_WFN, this keyword is not needed to "// &
1648 "apply the mixing defined in the WFN_MIX section. Default is False.", &
1649 usage="APPLY_WFN_MIX_INIT_RESTART", &
1650 default_l_val=.false., lone_keyword_l_val=.true.)
1651 CALL section_add_keyword(section, keyword)
1652 CALL keyword_release(keyword)
1653
1654 CALL keyword_create(keyword, __location__, name="APPLY_DELTA_PULSE", &
1655 description="Applies a delta kick to the initial wfn (only RTP for now - the EMD"// &
1656 " case is not yet implemented). Only work for INITIAL_WFN=SCF_WFN", &
1657 usage="APPLY_DELTA_PULSE", &
1658 default_l_val=.false., lone_keyword_l_val=.true.)
1659 CALL section_add_keyword(section, keyword)
1660 CALL keyword_release(keyword)
1661
1662 CALL keyword_create(keyword, __location__, name="APPLY_DELTA_PULSE_MAG", &
1663 description="Applies a magnetic delta kick to the initial wfn (only RTP for now - the EMD"// &
1664 " case is not yet implemented). Only work for INITIAL_WFN=SCF_WFN", &
1665 usage="APPLY_DELTA_PULSE_MAG", &
1666 default_l_val=.false., lone_keyword_l_val=.true.)
1667 CALL section_add_keyword(section, keyword)
1668 CALL keyword_release(keyword)
1669
1670 CALL keyword_create(keyword, __location__, name="VELOCITY_GAUGE", &
1671 description="Perform propagation in the velocity gauge using the explicit vector potential"// &
1672 " only a constant vector potential as of now (corresonding to a delta-pulse)."// &
1673 " uses DELTA_PULSE_SCALE and DELTA_PULSE_DIRECTION to define the vector potential", &
1674 usage="VELOCITY_GAUGE T", &
1675 default_l_val=.false., lone_keyword_l_val=.true.)
1676 CALL section_add_keyword(section, keyword)
1677 CALL keyword_release(keyword)
1678
1679 CALL keyword_create(keyword, __location__, name="GAUGE_ORIG", &
1680 description="Define gauge origin for magnetic perturbation", &
1681 usage="GAUGE_ORIG COM", &
1682 enum_c_vals=s2a("COM", "COAC", "USER_DEFINED", "ZERO"), &
1683 enum_desc=s2a("Use Center of Mass", &
1684 "Use Center of Atomic Charges", &
1685 "Use User Defined Point (Keyword:REF_POINT)", &
1686 "Use Origin of Coordinate System"), &
1687 enum_i_vals=[use_mom_ref_com, &
1691 default_i_val=use_mom_ref_com)
1692 CALL section_add_keyword(section, keyword)
1693 CALL keyword_release(keyword)
1694
1695 CALL keyword_create(keyword, __location__, name="GAUGE_ORIG_MANUAL", &
1696 description="Manually defined gauge origin for magnetic perturbation [in Bohr!]", &
1697 usage="GAUGE_ORIG_MANUAL x y z", &
1698 repeats=.false., &
1699 n_var=3, default_r_vals=[0._dp, 0._dp, 0._dp], &
1700 type_of_var=real_t, &
1701 unit_str='bohr')
1702 CALL section_add_keyword(section, keyword)
1703 CALL keyword_release(keyword)
1704
1705 CALL keyword_create(keyword, __location__, name="VG_COM_NL", &
1706 description="apply gauge transformed non-local potential term"// &
1707 " only affects VELOCITY_GAUGE=.TRUE.", &
1708 usage="VG_COM_NL T", &
1709 default_l_val=.true., lone_keyword_l_val=.true.)
1710 CALL section_add_keyword(section, keyword)
1711 CALL keyword_release(keyword)
1712
1713 CALL keyword_create(keyword, __location__, name="COM_NL", &
1714 description="Include non-local commutator for periodic delta pulse."// &
1715 " only affects PERIODIC=.TRUE.", &
1716 usage="COM_NL", &
1717 default_l_val=.true., lone_keyword_l_val=.true.)
1718 CALL section_add_keyword(section, keyword)
1719 CALL keyword_release(keyword)
1720
1721 CALL keyword_create(keyword, __location__, name="LEN_REP", &
1722 description="Use length representation delta pulse (in conjunction with PERIODIC T)."// &
1723 " This corresponds to a 1st order perturbation in the length gauge."// &
1724 " Note that this is NOT compatible with a periodic calculation!"// &
1725 " Uses the reference point defined in DFT%PRINT%MOMENTS ", &
1726 usage="LEN_REP T", &
1727 default_l_val=.false., lone_keyword_l_val=.true.)
1728 CALL section_add_keyword(section, keyword)
1729 CALL keyword_release(keyword)
1730
1731 CALL keyword_create(keyword, __location__, name="PERIODIC", &
1732 description="Apply a delta-kick that is compatible with periodic boundary conditions"// &
1733 " for any value of DELTA_PULSE_SCALE. Uses perturbation theory for the preparation of"// &
1734 " the initial wfn with the velocity operator as perturbation."// &
1735 " If LEN_REP is .FALSE. this corresponds to a first order velocity gauge."// &
1736 " Note that the pulse is only applied when INITIAL_WFN is set to SCF_WFN,"// &
1737 " and not for restarts (RT_RESTART).", &
1738 usage="PERIODIC", &
1739 default_l_val=.true., lone_keyword_l_val=.true.)
1740 CALL section_add_keyword(section, keyword)
1741 CALL keyword_release(keyword)
1742
1743 CALL keyword_create(keyword, __location__, name="DELTA_PULSE_DIRECTION", &
1744 description="Direction of the applied electric field. The k vector is given as"// &
1745 " 2*Pi*[i,j,k]*inv(h_mat), which for PERIODIC .FALSE. yields exp(ikr) periodic with"// &
1746 " the unit cell, only if DELTA_PULSE_SCALE is set to unity. For an orthorhombic cell"// &
1747 " [1,0,0] yields [2*Pi/L_x,0,0]. For small cells, this results in a very large kick.", &
1748 usage="DELTA_PULSE_DIRECTION 1 1 1", n_var=3, default_i_vals=[1, 0, 0], &
1749 type_of_var=integer_t)
1750 CALL section_add_keyword(section, keyword)
1751 CALL keyword_release(keyword)
1752
1753 CALL keyword_create(keyword, __location__, name="DELTA_PULSE_SCALE", &
1754 description="Scale the k vector, which for PERIODIC .FALSE. results in exp(ikr) no"// &
1755 " longer being periodic with the unit cell. The norm of k is the strength of the"// &
1756 " applied electric field in atomic units.", &
1757 usage="DELTA_PULSE_SCALE 0.01 ", n_var=1, default_r_val=0.001_dp)
1758 CALL section_add_keyword(section, keyword)
1759 CALL keyword_release(keyword)
1760
1761 CALL keyword_create(keyword, __location__, name="HFX_BALANCE_IN_CORE", &
1762 description="If HFX is used, this keyword forces a redistribution/recalculation"// &
1763 " of the integrals, balanced with respect to the in core steps.", &
1764 usage="HFX_BALANCE_IN_CORE", &
1765 default_l_val=.false., lone_keyword_l_val=.true.)
1766 CALL section_add_keyword(section, keyword)
1767 CALL keyword_release(keyword)
1768
1769 CALL keyword_create(keyword, __location__, name="MCWEENY_MAX_ITER", &
1770 description="Determines the maximum amount of McWeeny steps used after each converged"// &
1771 " step in density propagation", &
1772 usage="MCWEENY_MAX_ITER 2", default_i_val=1)
1773 CALL section_add_keyword(section, keyword)
1774 CALL keyword_release(keyword)
1775
1776 CALL keyword_create( &
1777 keyword, __location__, name="ACCURACY_REFINEMENT", &
1778 description="If using density propagation some parts should be calculated with a higher accuracy than the rest"// &
1779 " to reduce numerical noise. This factor determines by how much the filtering threshold is"// &
1780 " reduced for these calculations.", &
1781 usage="ACCURACY_REFINEMENT", default_i_val=100)
1782 CALL section_add_keyword(section, keyword)
1783 CALL keyword_release(keyword)
1784
1785 CALL keyword_create(keyword, __location__, name="MCWEENY_EPS", &
1786 description="Threshold after which McWeeny is terminated", &
1787 usage="MCWEENY_EPS 0.00001", &
1788 default_r_val=0.0_dp)
1789 CALL section_add_keyword(section, keyword)
1790 CALL keyword_release(keyword)
1791
1792 NULLIFY (print_section)
1793 CALL section_create(print_section, __location__, name="PRINT", &
1794 description="Section of possible print options for an RTP runs", &
1795 repeats=.false.)
1796
1797 NULLIFY (print_key)
1798 CALL cp_print_key_section_create(print_key, __location__, "PROGRAM_RUN_INFO", &
1799 description="Controls the printing within real time propagation and Eherenfest dynamics", &
1800 print_level=low_print_level, filename="__STD_OUT__")
1801 CALL section_add_subsection(print_section, print_key)
1802 CALL section_release(print_key)
1803
1804 CALL cp_print_key_section_create(print_key, __location__, "RESTART", &
1805 description="Controls the dumping of the MO restart file during rtp. "// &
1806 "By default keeps a short history of three restarts. "// &
1807 "See also RESTART_HISTORY. In density propagation this controls the printing of "// &
1808 "density matrix.", &
1809 print_level=low_print_level, common_iter_levels=3, &
1810 each_iter_names=s2a("MD"), each_iter_values=[20], &
1811 add_last=add_last_numeric, filename="RESTART")
1812 CALL keyword_create(keyword, __location__, name="BACKUP_COPIES", &
1813 description="Specifies the maximum number of backup copies.", &
1814 usage="BACKUP_COPIES {int}", &
1815 default_i_val=1)
1816 CALL section_add_keyword(print_key, keyword)
1817 CALL keyword_release(keyword)
1818 CALL section_add_subsection(print_section, print_key)
1819 CALL section_release(print_key)
1820
1821 CALL cp_print_key_section_create(print_key, __location__, "RESTART_HISTORY", &
1822 description="Dumps unique MO restart files during the run keeping all of them. "// &
1823 "In density propagation it dumps the density matrix instead", &
1824 print_level=low_print_level, common_iter_levels=0, &
1825 each_iter_names=s2a("MD"), &
1826 each_iter_values=[500], &
1827 filename="RESTART")
1828 CALL keyword_create(keyword, __location__, name="BACKUP_COPIES", &
1829 description="Specifies the maximum number of backup copies.", &
1830 usage="BACKUP_COPIES {int}", &
1831 default_i_val=1)
1832 CALL section_add_keyword(print_key, keyword)
1833 CALL keyword_release(keyword)
1834 CALL section_add_subsection(print_section, print_key)
1835 CALL section_release(print_key)
1836
1837 CALL cp_print_key_section_create(print_key, __location__, "FIELD", &
1838 description="Print the time-dependent field applied during an EMD simulation in "// &
1839 "atomic unit.", &
1840 print_level=high_print_level, common_iter_levels=1, &
1841 each_iter_names=s2a("MD"), &
1842 each_iter_values=[1], &
1843 filename="FIELD")
1844 CALL section_add_subsection(print_section, print_key)
1845 CALL section_release(print_key)
1846
1847 CALL create_projection_rtp_section(print_key)
1848 CALL section_add_subsection(print_section, print_key)
1849 CALL section_release(print_key)
1850
1851 CALL cp_print_key_section_create(print_key, __location__, "CURRENT_INT", &
1852 description="Print the integral of the current density (only if the"// &
1853 " imaginary part of the density is NOT zero.", &
1854 print_level=high_print_level, common_iter_levels=1, &
1855 each_iter_names=s2a("MD"), &
1856 each_iter_values=[1], &
1857 filename="rtp_j_int")
1858 CALL section_add_subsection(print_section, print_key)
1859 CALL section_release(print_key)
1860
1861 CALL cp_print_key_section_create(print_key, __location__, "CURRENT", &
1862 description="Print the current during an EMD simulation to cube files.", &
1863 print_level=high_print_level, common_iter_levels=0, &
1864 each_iter_names=s2a("MD"), &
1865 each_iter_values=[20], &
1866 filename="current")
1867 CALL keyword_create(keyword, __location__, name="BACKUP_COPIES", &
1868 description="Specifies the maximum number of backup copies.", &
1869 usage="BACKUP_COPIES {int}", &
1870 default_i_val=1)
1871 CALL section_add_keyword(print_key, keyword)
1872 CALL keyword_release(keyword)
1873 CALL keyword_create(keyword, __location__, name="STRIDE", &
1874 description="The stride (X,Y,Z) used to write the cube file "// &
1875 "(larger values result in smaller cube files). You can provide 3 numbers (for X,Y,Z) or"// &
1876 " 1 number valid for all components.", &
1877 usage="STRIDE 2 2 2", n_var=-1, default_i_vals=[2, 2, 2], type_of_var=integer_t)
1878 CALL section_add_keyword(print_key, keyword)
1879 CALL keyword_release(keyword)
1880
1881 CALL section_add_subsection(print_section, print_key)
1882 CALL section_release(print_key)
1883
1884 ! Marek : Add print option for ASCII density files - DEVELPMENT ONLY?
1885 CALL cp_print_key_section_create(print_key, __location__, "DENSITY_MATRIX", &
1886 description="Prints the density matrix at iterations in clear text to a file", &
1887 print_level=high_print_level, common_iter_levels=0, &
1888 each_iter_names=s2a("MD"), &
1889 each_iter_values=[1], &
1890 filename="rho")
1891 CALL section_add_subsection(print_section, print_key)
1892 CALL section_release(print_key)
1893 ! Marek : Moments ASCII print
1894 CALL cp_print_key_section_create(print_key, __location__, "MOMENTS", &
1895 description="Prints the time-dependent electronic moments at "// &
1896 "iterations in clear text to a file.", &
1897 print_level=high_print_level, common_iter_levels=1, &
1898 each_iter_names=s2a("MD"), &
1899 each_iter_values=[1], &
1900 filename="__STD_OUT__")
1901 CALL keyword_create(keyword, __location__, name="REFERENCE", &
1902 variants=s2a("REF"), &
1903 description="Define the reference point for the calculation of the electrostatic moment.", &
1904 usage="REFERENCE COM", &
1905 enum_c_vals=s2a("COM", "COAC", "USER_DEFINED", "ZERO"), &
1906 enum_desc=s2a("Use Center of Mass", &
1907 "Use Center of Atomic Charges", &
1908 "Use User Defined Point (Keyword:REFERENCE_POINT)", &
1909 "Use Origin of Coordinate System"), &
1910 enum_i_vals=[use_mom_ref_com, &
1914 default_i_val=use_mom_ref_coac)
1915 CALL section_add_keyword(print_key, keyword)
1916 CALL keyword_release(keyword)
1917
1918 CALL keyword_create(keyword, __location__, name="REFERENCE_POINT", &
1919 variants=s2a("REF_POINT"), &
1920 description="Fixed reference point for the calculations of the electrostatic moment.", &
1921 usage="REFERENCE_POINT x y z", &
1922 repeats=.false., &
1923 n_var=3, default_r_vals=[0._dp, 0._dp, 0._dp], &
1924 type_of_var=real_t, &
1925 unit_str='angstrom')
1926 CALL section_add_keyword(print_key, keyword)
1927 CALL keyword_release(keyword)
1928 CALL section_add_subsection(print_section, print_key)
1929 CALL section_release(print_key)
1930 ! Marek : Fourier transform of MOMENTS ASCII print
1931 CALL cp_print_key_section_create(print_key, __location__, "MOMENTS_FT", &
1932 description="Prints the calculated Fourier transform of "// &
1933 "time-dependent moments. For calculations with real time pulse (not delta kick) "// &
1934 "can be supplied with starting time.", &
1935 print_level=medium_print_level, common_iter_levels=0, &
1936 each_iter_names=s2a("MD"), &
1937 each_iter_values=[1], &
1938 filename="MOMENTS_FT")
1939 CALL section_add_subsection(print_section, print_key)
1940 CALL section_release(print_key)
1941 ! Marek : Chosen element of (Fourier transformed) polarizability tensor (energy dependent) - text format
1942 CALL cp_print_key_section_create(print_key, __location__, "POLARIZABILITY", &
1943 description="Prints the chosen element of the energy dependent polarizability tensor "// &
1944 "to a specified file. The tensor is calculated as ratio of "// &
1945 "Fourier transform of the dipole "// &
1946 "moment trace and Fourier transform of the applied field "// &
1947 "(for delta kick, constant real field is applied.", &
1948 print_level=medium_print_level, common_iter_levels=0, &
1949 each_iter_names=s2a("MD"), &
1950 each_iter_values=[1], &
1951 filename="POLARIZABILITY")
1952 CALL keyword_create(keyword, __location__, "ELEMENT", &
1953 description="Specifies the element of polarizability which is to be printed out "// &
1954 "(indexing starts at 1). If not explicitly provided, RTBSE code tries to guess "// &
1955 "the optimal values - for applied electric field (both delta pulse and RT field) "// &
1956 "with only a single non-zero cartesian component, prints the 3 trivially available elements.", &
1957 type_of_var=integer_t, default_i_vals=[1, 1], n_var=2, usage="ELEMENT 1 1", repeats=.true.)
1958 CALL section_add_keyword(print_key, keyword)
1959 CALL keyword_release(keyword)
1960 CALL section_add_subsection(print_section, print_key)
1961 CALL section_release(print_key)
1962
1963 CALL cp_print_key_section_create(print_key, __location__, "E_CONSTITUENTS", &
1964 description="Print the energy constituents (relevant to RTP) which make up "// &
1965 "the Total Energy", &
1966 print_level=high_print_level, common_iter_levels=1, &
1967 each_iter_names=s2a("MD"), &
1968 each_iter_values=[1], &
1969 filename="rtp")
1970 CALL section_add_subsection(print_section, print_key)
1971 CALL section_release(print_key)
1972
1973 CALL section_add_subsection(section, print_section)
1974 CALL section_release(print_section)
1975
1976 ! RTBSE subsection
1977 NULLIFY (subsection)
1978 CALL create_rtbse_section(subsection)
1979 CALL section_add_subsection(section, subsection)
1980 CALL section_release(subsection)
1981 ! FT subsection
1982 CALL create_ft_section(subsection)
1983 CALL section_add_subsection(section, subsection)
1984 CALL section_release(subsection)
1985
1986 END SUBROUTINE create_rtp_section
1987! **************************************************************************************************
1988!> \brief Creates the subsection for specialized options of RTBSE code
1989!> \param section The created RTBSE section
1990!> \author Stepan Marek
1991! **************************************************************************************************
1992 SUBROUTINE create_rtbse_section(section)
1993 TYPE(section_type), POINTER :: section
1994
1995 TYPE(keyword_type), POINTER :: keyword
1996
1997 NULLIFY (keyword)
1998 cpassert(.NOT. ASSOCIATED(section))
1999
2000 CALL section_create(section, __location__, name="RTBSE", &
2001 description="Controls options for the real-time Bethe-Salpeter (RTBSE) propagation. "// &
2002 "Note that running RTBSE requires previous low-scaling "// &
2003 "[GW](#CP2K_INPUT.FORCE_EVAL.PROPERTIES.BANDSTRUCTURE.GW) calculation. Also note that "// &
2004 "designating this section as RTBSE run but choosing run type ENERGY leads to potential "// &
2005 "deallocation errors. More details (including description of output files) is available in "// &
2006 "the [methods](../../../../methods/properties/optical/rtbse) section of the documentation.", &
2007 repeats=.false., citations=[marek2025])
2008
2009 ! Marek : Controlling flow to RTBSE
2010 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
2011 description="Which method is used for the time propagation of electronic structure. "// &
2012 "By default, use the TDDFT method. Can also choose RT-BSE method, which propagates the lesser "// &
2013 "Green's function instead of density matrix/molecular orbitals.", &
2014 usage="&RTBSE TDDFT", &
2015 default_i_val=rtp_method_tddft, &
2016 lone_keyword_i_val=rtp_method_bse, &
2017 enum_c_vals=s2a("TDDFT", "RTBSE"), &
2018 enum_i_vals=[rtp_method_tddft, rtp_method_bse], &
2019 enum_desc=s2a("Use TDDFT for density matrix/MO propagation.", &
2020 "Use RT-BSE for Green's function propagation"))
2021 CALL section_add_keyword(section, keyword)
2022 CALL keyword_release(keyword)
2023
2024 ! Marek : Development option - run GWBSE starting from the KS Hamiltonian
2025 CALL keyword_create(keyword, __location__, name="RTBSE_HAMILTONIAN", &
2026 description="Which Hamiltonian to use as the single-particle Hamiltonian"// &
2027 " in the Green's propagator.", &
2028 usage="RTBSE_HAMILTONIAN G0W0", &
2029 default_i_val=rtp_bse_ham_g0w0, &
2030 enum_c_vals=s2a("KS", "G0W0"), &
2031 enum_i_vals=[rtp_bse_ham_ks, rtp_bse_ham_g0w0], &
2032 enum_desc=s2a("Use Kohn-Sham Hamiltonian for Green's propagation.", &
2033 "Use G0W0 Hamiltonian for Green's function propagation"))
2034 CALL section_add_keyword(section, keyword)
2035 CALL keyword_release(keyword)
2036
2037 END SUBROUTINE create_rtbse_section
2038! **************************************************************************************************
2039!> \brief Creates the subsection for Fourier transform options applicable to RTP output
2040!> \param ft_section The created FT section
2041!> \date 11.2025
2042!> \author Stepan Marek
2043! **************************************************************************************************
2044 SUBROUTINE create_ft_section(ft_section)
2045 TYPE(section_type), POINTER :: ft_section
2046
2047 TYPE(keyword_type), POINTER :: keyword
2048 TYPE(section_type), POINTER :: subsection
2049
2050 cpassert(.NOT. ASSOCIATED(ft_section))
2051
2052 ! Create the section itself
2053 CALL section_create(ft_section, __location__, &
2054 name="FT", &
2055 description="Define parameters for Fourier transforms used in RTP outputs.", &
2056 repeats=.false.)
2057
2058 ! Start time keyword
2059 NULLIFY (keyword)
2060 CALL keyword_create(keyword, __location__, "START_TIME", &
2061 description="The starting time from which damping is applied and from which on the trace is "// &
2062 "considered for the Fourier transform (Fourier transform is used for the calculation of "// &
2063 "MOMENTS_FT and POLARIZABILITY). Useful for real-time pulse - "// &
2064 "one can specify the center of the pulse as the starting point.", &
2065 type_of_var=real_t, &
2066 unit_str="fs", &
2067 default_r_val=0.0_dp)
2068 CALL section_add_keyword(ft_section, keyword)
2069 CALL keyword_release(keyword)
2070
2071 ! Damping keyword
2072 CALL keyword_create(keyword, __location__, "DAMPING", &
2073 description="Numerical Fourier transform (required for calculation of "// &
2074 "MOMENTS_FT and POLARIZABILITY) can oscillate "// &
2075 "when the final time trace values are far away from zero. "// &
2076 "This keyword controls the exponential damping added to the Fourier transform "// &
2077 "(Fourier transform is used for calculation of MOMENTS_FT and POLARIZABILITY). "// &
2078 "For negative values (the default), calculates the damping at the run time so that the last point "// &
2079 "in the time trace is reduced by factor e^(-4). When set manually, determines the time in which "// &
2080 "the moments trace is reduced by factor of e^(-1), except when set to zero, in which case "// &
2081 "the damping is not applied.", &
2082 type_of_var=real_t, &
2083 unit_str="fs", &
2084 default_r_val=-1.0_dp/femtoseconds)
2085 CALL section_add_keyword(ft_section, keyword)
2086 CALL keyword_release(keyword)
2087
2088 ! Create the Padé subsection
2089 NULLIFY (subsection)
2090 CALL section_create(subsection, __location__, name="PADE", &
2091 description=é"Defines the parameters for the Pad interpolation of the "// &
2092 "Fourier transforms used in the output of RTP. Only available with the GreenX library linked to CP2K.", &
2093 repeats=.false.)
2094
2095 ! Explicit presence of the section turns on the Padé interpolation
2096 CALL keyword_create(keyword, __location__, "_SECTION_PARAMETERS_", &
2097 description=é"Turns on the Pad interpolation", &
2098 type_of_var=logical_t, &
2099 default_l_val=.false., &
2100 lone_keyword_l_val=.true.)
2101 CALL section_add_keyword(subsection, keyword)
2102 CALL keyword_release(keyword)
2103
2104 ! Minimum interpolated energy
2105 CALL keyword_create(keyword, __location__, "E_MIN", &
2106 description=é"The minimum energy of the Pad interpolation output.", &
2107 type_of_var=real_t, &
2108 unit_str="eV", &
2109 default_r_val=0.0_dp)
2110 CALL section_add_keyword(subsection, keyword)
2111 CALL keyword_release(keyword)
2112
2113 ! Maximum interpolated energy
2114 CALL keyword_create(keyword, __location__, "E_MAX", &
2115 description=é"The maximum energy of the Pad interpolation output.", &
2116 type_of_var=real_t, &
2117 unit_str="eV", &
2118 default_r_val=100.0_dp)
2119 CALL section_add_keyword(subsection, keyword)
2120 CALL keyword_release(keyword)
2121
2122 ! Energy resolution
2123 CALL keyword_create(keyword, __location__, "E_STEP", &
2124 description=é"The energy resolution of the Pad interpolation output.", &
2125 type_of_var=real_t, &
2126 unit_str="eV", &
2127 default_r_val=0.02_dp/evolt)
2128 CALL section_add_keyword(subsection, keyword)
2129 CALL keyword_release(keyword)
2130
2131 ! Minimum fitting energy
2132 CALL keyword_create(keyword, __location__, "FIT_E_MIN", &
2133 description="The lower boundary in energy for the points "// &
2134 é"used in the fitting of Pad parameters. If negative, uses "// &
2135 "value of E_MIN (default).", &
2136 type_of_var=real_t, &
2137 unit_str="eV", &
2138 default_r_val=-1.0_dp)
2139 CALL section_add_keyword(subsection, keyword)
2140 CALL keyword_release(keyword)
2141
2142 ! Maximum fitting energy
2143 CALL keyword_create(keyword, __location__, "FIT_E_MAX", &
2144 description="The upper boundary in energy for the points "// &
2145 é"used in the fitting of Pad parameters. If negative, uses "// &
2146 "the value of E_MAX (default).", &
2147 type_of_var=real_t, &
2148 unit_str="eV", &
2149 default_r_val=-1.0_dp)
2150 CALL section_add_keyword(subsection, keyword)
2151 CALL keyword_release(keyword)
2152
2153 ! Add the Padé subsection
2154 CALL section_add_subsection(ft_section, subsection)
2155 CALL section_release(subsection)
2156 END SUBROUTINE create_ft_section
2157
2158! **************************************************************************************************
2159!> \brief Create CP2K input section for the SCCS model
2160!> \param section ...
2161!> \par History:
2162!> - Creation (10.10.2013,MK)
2163!> \author Matthias Krack (MK)
2164!> \version 1.0
2165! **************************************************************************************************
2166 SUBROUTINE create_sccs_section(section)
2167
2168 TYPE(section_type), POINTER :: section
2169
2170 TYPE(keyword_type), POINTER :: keyword
2171 TYPE(section_type), POINTER :: subsection
2172
2173 cpassert(.NOT. ASSOCIATED(section))
2174
2175 CALL section_create(section, __location__, &
2176 name="SCCS", &
2177 description="Define the parameters for self-consistent continuum solvation (SCCS) model", &
2178 citations=[fattebert2002, andreussi2012, yin2017], &
2179 n_keywords=8, &
2180 n_subsections=2, &
2181 repeats=.false.)
2182
2183 NULLIFY (keyword)
2184
2185 CALL keyword_create(keyword, __location__, &
2186 name="_SECTION_PARAMETERS_", &
2187 description="Controls the activation of the SCCS section", &
2188 usage="&SCCS ON", &
2189 default_l_val=.false., &
2190 lone_keyword_l_val=.true.)
2191 CALL section_add_keyword(section, keyword)
2192 CALL keyword_release(keyword)
2193
2194 CALL keyword_create(keyword, __location__, &
2195 name="ALPHA", &
2196 description="Solvent specific tunable parameter for the calculation of "// &
2197 "the repulsion term $G^\text{rep} = \alpha S$ "// &
2198 "where $S$ is the (quantum) surface of the cavity", &
2199 repeats=.false., &
2200 n_var=1, &
2201 type_of_var=real_t, &
2202 default_r_val=0.0_dp, &
2203 unit_str="mN/m")
2204 CALL section_add_keyword(section, keyword)
2205 CALL keyword_release(keyword)
2206
2207 CALL keyword_create(keyword, __location__, &
2208 name="BETA", &
2209 description="Solvent specific tunable parameter for the calculation of "// &
2210 "the dispersion term $G^\text{dis} = \beta V$ "// &
2211 "where $V$ is the (quantum) volume of the cavity", &
2212 repeats=.false., &
2213 n_var=1, &
2214 type_of_var=real_t, &
2215 default_r_val=0.0_dp, &
2216 unit_str="GPa")
2217 CALL section_add_keyword(section, keyword)
2218 CALL keyword_release(keyword)
2219
2220 CALL keyword_create(keyword, __location__, &
2221 name="DELTA_RHO", &
2222 description="Numerical increment for the calculation of the (quantum) "// &
2223 "surface of the solute cavity", &
2224 repeats=.false., &
2225 n_var=1, &
2226 type_of_var=real_t, &
2227 default_r_val=2.0e-5_dp)
2228 CALL section_add_keyword(section, keyword)
2229 CALL keyword_release(keyword)
2230
2231 CALL keyword_create(keyword, __location__, &
2232 name="DERIVATIVE_METHOD", &
2233 description="Method for the calculation of the numerical derivatives on the real-space grids", &
2234 usage="DERIVATIVE_METHOD cd5", &
2235 repeats=.false., &
2236 n_var=1, &
2237 default_i_val=sccs_derivative_fft, &
2238 enum_c_vals=s2a("FFT", "CD3", "CD5", "CD7"), &
2239 enum_i_vals=[sccs_derivative_fft, &
2243 enum_desc=s2a("Fast Fourier transformation", &
2244 "3-point stencil central differences", &
2245 "5-point stencil central differences", &
2246 "7-point stencil central differences"))
2247 CALL section_add_keyword(section, keyword)
2248 CALL keyword_release(keyword)
2249
2250 CALL keyword_create(keyword, __location__, &
2251 name="RELATIVE_PERMITTIVITY", &
2252 variants=s2a("DIELECTRIC_CONSTANT", "EPSILON_RELATIVE", "EPSILON_SOLVENT"), &
2253 description="Relative permittivity (dielectric constant) of the solvent (medium)", &
2254 repeats=.false., &
2255 n_var=1, &
2256 type_of_var=real_t, &
2257 default_r_val=80.0_dp, &
2258 usage="RELATIVE_PERMITTIVITY 78.36")
2259 CALL section_add_keyword(section, keyword)
2260 CALL keyword_release(keyword)
2261
2262 CALL keyword_create(keyword, __location__, &
2263 name="EPS_SCCS", &
2264 variants=s2a("EPS_ITER", "TAU_POL"), &
2265 description="Tolerance for the convergence of the polarisation density, "// &
2266 "i.e. requested accuracy for the SCCS iteration cycle", &
2267 repeats=.false., &
2268 n_var=1, &
2269 type_of_var=real_t, &
2270 default_r_val=1.0e-6_dp, &
2271 usage="EPS_ITER 1.0E-7")
2272 CALL section_add_keyword(section, keyword)
2273 CALL keyword_release(keyword)
2274
2275 CALL keyword_create(keyword, __location__, &
2276 name="EPS_SCF", &
2277 description="The SCCS iteration cycle is activated only if the SCF iteration cycle "// &
2278 "is converged to this threshold value", &
2279 repeats=.false., &
2280 n_var=1, &
2281 type_of_var=real_t, &
2282 default_r_val=0.5_dp, &
2283 usage="EPS_SCF 1.0E-2")
2284 CALL section_add_keyword(section, keyword)
2285 CALL keyword_release(keyword)
2286
2287 CALL keyword_create(keyword, __location__, &
2288 name="GAMMA", &
2289 variants=s2a("SURFACE_TENSION"), &
2290 description="Surface tension of the solvent used for the calculation of "// &
2291 "the cavitation term $G^\text{cav} = \gamma S$ "// &
2292 "where $S$ is the (quantum) surface of the cavity", &
2293 repeats=.false., &
2294 n_var=1, &
2295 type_of_var=real_t, &
2296 default_r_val=0.0_dp, &
2297 unit_str="mN/m")
2298 CALL section_add_keyword(section, keyword)
2299 CALL keyword_release(keyword)
2300
2301 CALL keyword_create(keyword, __location__, &
2302 name="MAX_ITER", &
2303 description="Maximum number of SCCS iteration steps performed to converge "// &
2304 "within the given tolerance", &
2305 repeats=.false., &
2306 n_var=1, &
2307 type_of_var=integer_t, &
2308 default_i_val=100, &
2309 usage="MAX_ITER 50")
2310 CALL section_add_keyword(section, keyword)
2311 CALL keyword_release(keyword)
2312
2313 CALL keyword_create(keyword, __location__, &
2314 name="METHOD", &
2315 description="Method used for the smoothing of the dielectric function", &
2316 usage="METHOD Fattebert-Gygi", &
2317 default_i_val=sccs_andreussi, &
2318 enum_c_vals=s2a("ANDREUSSI", "FATTEBERT-GYGI"), &
2319 enum_i_vals=[sccs_andreussi, sccs_fattebert_gygi], &
2320 enum_desc=s2a("Smoothing function proposed by Andreussi et al.", &
2321 "Smoothing function proposed by Fattebert and Gygi"))
2322 CALL section_add_keyword(section, keyword)
2323 CALL keyword_release(keyword)
2324
2325 CALL keyword_create(keyword, __location__, &
2326 name="MIXING", &
2327 variants=["ETA"], &
2328 description="Mixing parameter (Hartree damping) employed during the iteration procedure", &
2329 repeats=.false., &
2330 n_var=1, &
2331 type_of_var=real_t, &
2332 default_r_val=0.6_dp, &
2333 usage="MIXING 0.2")
2334 CALL section_add_keyword(section, keyword)
2335 CALL keyword_release(keyword)
2336
2337 NULLIFY (subsection)
2338
2339 CALL section_create(subsection, __location__, &
2340 name="ANDREUSSI", &
2341 description="Define the parameters of the dielectric smoothing function proposed by "// &
2342 "Andreussi et al.", &
2343 citations=[andreussi2012], &
2344 n_keywords=2, &
2345 n_subsections=0, &
2346 repeats=.false.)
2347
2348 CALL keyword_create(keyword, __location__, &
2349 name="RHO_MAX", &
2350 description="Maximum density value used for the smoothing of the dielectric function", &
2351 repeats=.false., &
2352 n_var=1, &
2353 type_of_var=real_t, &
2354 default_r_val=0.0035_dp, &
2355 usage="RHO_MAX 0.01")
2356 CALL section_add_keyword(subsection, keyword)
2357 CALL keyword_release(keyword)
2358
2359 CALL keyword_create(keyword, __location__, &
2360 name="RHO_MIN", &
2361 description="Minimum density value used for the smoothing of the dielectric function", &
2362 repeats=.false., &
2363 n_var=1, &
2364 type_of_var=real_t, &
2365 default_r_val=0.0001_dp, &
2366 usage="RHO_MIN 0.0003")
2367 CALL section_add_keyword(subsection, keyword)
2368 CALL keyword_release(keyword)
2369
2370 CALL section_add_subsection(section, subsection)
2371 CALL section_release(subsection)
2372
2373 CALL section_create(subsection, __location__, &
2374 name="FATTEBERT-GYGI", &
2375 description="Define the parameters of the dielectric smoothing function proposed by "// &
2376 "Fattebert and Gygi", &
2377 citations=[fattebert2002], &
2378 n_keywords=2, &
2379 n_subsections=0, &
2380 repeats=.false.)
2381
2382 CALL keyword_create(keyword, __location__, &
2383 name="BETA", &
2384 description="Parameter &beta; changes the width of the interface solute-solvent", &
2385 repeats=.false., &
2386 n_var=1, &
2387 type_of_var=real_t, &
2388 default_r_val=1.7_dp, &
2389 usage="BETA 1.3")
2390 CALL section_add_keyword(subsection, keyword)
2391 CALL keyword_release(keyword)
2392
2393 CALL keyword_create(keyword, __location__, &
2394 name="RHO_ZERO", &
2395 variants=["RHO0"], &
2396 description="Parameter $\rho_0$ defines the critical density in the middle "// &
2397 "of the interface solute-solvent", &
2398 repeats=.false., &
2399 n_var=1, &
2400 type_of_var=real_t, &
2401 default_r_val=0.0006_dp, &
2402 usage="RHO_ZERO 0.0004")
2403 CALL section_add_keyword(subsection, keyword)
2404 CALL keyword_release(keyword)
2405
2406 CALL section_add_subsection(section, subsection)
2407 CALL section_release(subsection)
2408
2409 END SUBROUTINE create_sccs_section
2410
2411END MODULE input_cp2k_dft
integer, parameter, public basis_sort_zet
integer, parameter, public basis_sort_default
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public vandevondele2005b
integer, save, public blochl1995
integer, save, public guidon2010
integer, save, public bengtsson1999
integer, save, public kunert2003
integer, save, public yin2017
integer, save, public avezac2005
integer, save, public andreussi2012
integer, save, public iannuzzi2006
integer, save, public fattebert2002
integer, save, public andermatt2016
integer, save, public merlot2014
integer, save, public perdew1981
integer, save, public marek2025
integer, save, public brelaz1979
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 high_print_level
integer, parameter, public add_last_numeric
integer, parameter, public silent_print_level
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
utils to manipulate splines on the regular grid of a pw
integer, parameter, public pw_interp
integer, parameter, public spline3_nopbc_interp
integer, parameter, public spline3_pbc_interp
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:1149
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public use_mom_ref_coac
integer, parameter, public sic_list_unpaired
integer, parameter, public sic_mauri_spz
integer, parameter, public do_bch
integer, parameter, public do_admm_purify_mo_no_diag
integer, parameter, public do_etrs
integer, parameter, public do_admm_aux_exch_func_opt_libxc
integer, parameter, public rel_zora_full
integer, parameter, public do_pade
integer, parameter, public do_admm_purify_none
integer, parameter, public use_mom_ref_user
integer, parameter, public rel_pot_full
integer, parameter, public do_admm_purify_none_dm
integer, parameter, public use_mom_ref_com
integer, parameter, public ehrenfest
integer, parameter, public use_restart_wfn
integer, parameter, public do_admm_purify_mcweeny
integer, parameter, public do_admm_blocking_purify_full
integer, parameter, public plus_u_lowdin
integer, parameter, public kg_tnadd_none
integer, parameter, public do_admm_aux_exch_func_sx_libxc
integer, parameter, public admm2_type
integer, parameter, public sic_list_all
integer, parameter, public do_cn
integer, parameter, public kg_tnadd_embed_ri
integer, parameter, public kg_tnadd_embed
integer, parameter, public sic_eo
integer, parameter, public sccs_derivative_cd5
integer, parameter, public do_admm_aux_exch_func_bee
integer, parameter, public rel_zora_mp
integer, parameter, public plus_u_mulliken_charges
integer, parameter, public use_scf_wfn
integer, parameter, public no_admm_type
integer, parameter, public rel_zora
integer, parameter, public do_admm_blocked_projection
integer, parameter, public kg_tnadd_atomic
integer, parameter, public do_admm_basis_projection
integer, parameter, public do_admm_aux_exch_func_default_libxc
integer, parameter, public do_admm_aux_exch_func_opt
integer, parameter, public do_admm_aux_exch_func_none
integer, parameter, public do_admm_purify_cauchy_subspace
integer, parameter, public plus_u_mulliken
integer, parameter, public kg_color_greedy
integer, parameter, public do_admm_aux_exch_func_bee_libxc
integer, parameter, public admm1_type
integer, parameter, public do_admm_aux_exch_func_pbex_libxc
integer, parameter, public kg_color_dsatur
integer, parameter, public do_admm_aux_exch_func_default
integer, parameter, public weight_type_unit
integer, parameter, public rtp_method_bse
integer, parameter, public admms_type
integer, parameter, public do_admm_charge_constrained_projection
integer, parameter, public rel_dkh
integer, parameter, public do_admm_purify_cauchy
integer, parameter, public sccs_fattebert_gygi
integer, parameter, public rel_trans_full
integer, parameter, public sccs_derivative_cd7
integer, parameter, public rel_trans_molecule
integer, parameter, public rel_trans_atom
integer, parameter, public weight_type_mass
integer, parameter, public use_rt_restart
integer, parameter, public do_exact
integer, parameter, public sccs_derivative_fft
integer, parameter, public use_mom_ref_zero
integer, parameter, public gaussian
integer, parameter, public rel_pot_erfc
integer, parameter, public rtp_method_tddft
integer, parameter, public rel_none
integer, parameter, public rtp_bse_ham_g0w0
integer, parameter, public do_taylor
integer, parameter, public do_admm_purify_mo_diag
integer, parameter, public do_em
integer, parameter, public sic_mauri_us
integer, parameter, public sic_none
integer, parameter, public rtp_bse_ham_ks
integer, parameter, public sccs_derivative_cd3
integer, parameter, public rel_sczora_mp
integer, parameter, public admmq_type
integer, parameter, public sccs_andreussi
integer, parameter, public sic_ad
integer, parameter, public do_admm_exch_scaling_none
integer, parameter, public do_arnoldi
integer, parameter, public admmp_type
integer, parameter, public do_admm_exch_scaling_merlot
integer, parameter, public real_time_propagation
integer, parameter, public numerical
integer, parameter, public do_admm_aux_exch_func_pbex
integer, parameter, public slater
input for the ALMO SCF section
subroutine, public create_almo_scf_section(section)
create the almo scf section
function that build the active space section of the input
subroutine, public create_active_space_section(section)
Create CP2K input section for the calculation of an active space Hamiltonian.
function that build the dft section of the input
subroutine, public create_bsse_section(section)
Create the BSSE section for counterpoise correction.
subroutine, public create_mgrid_section(section, create_subsections)
creates the multigrid
subroutine, public create_dft_section(section)
creates the dft section
subroutine, public create_interp_section(section)
creates the interpolation section
function that build the dft section of the input
subroutine, public create_ec_section(section)
creates the ENERGY CORRECTION section
Excited state input section.
subroutine, public create_exstate_section(section)
creates the EXCITED ENERGY section
function that build the input sections for external [potential, density VXC]
subroutine, public create_ext_pot_section(section)
Creates the section for applying an electrostatic external potential.
subroutine, public create_ext_vxc_section(section)
ZMP Creates the section for creating the external v_xc.
subroutine, public create_ext_den_section(section)
ZMP Creates the section for reading user supplied external density.
function that build the field section of the input
subroutine, public create_efield_section(section)
creates the section for time dependent nonperiodic fields
subroutine, public create_per_efield_section(section)
creates the section for static periodic fields
Harris input section.
subroutine, public create_harris_section(section)
creates the HARRIS_METHOD section
function that build the kpoints section of the input
subroutine, public create_kpoints_section(section)
Creates the Kpoints section SECTION: &kpoint... &end SCHEME [None, Gamma, Monkhorst-Pack,...
subroutine, public create_kpoint_set_section(section, section_name)
...
subroutine, public create_localize_section(section)
parameters fo the localization of wavefunctions
input for the linear scaling (LS) section
subroutine, public create_ls_scf_section(section)
creates the linear scaling scf section
function that build the poisson section of the input
subroutine, public create_poisson_section(section)
Creates the Poisson section.
function that build the print section of the dft input
subroutine, public create_print_dft_section(section)
Create the print dft section.
function that builds the projection of MO in RTP section of the input
subroutine, public create_projection_rtp_section(section)
creates the section for time dependent projection of the MOs
function that build the QS section of the input
subroutine, public create_qs_section(section)
creates the input section for the qs part
subroutine, public create_lrigpw_section(section)
input section for optional parameters for LRIGPW LRI: local resolution of identity
subroutine, public create_rsgrid_section(section)
...
function that build the scf section of the input
subroutine, public create_scf_section(section)
creates the structure of the section with the DFT SCF parameters
Functions that build SMEAGOL input section.
subroutine, public create_dft_smeagol_section(section)
Create SMEAGOL input section.
input section for NEGF based quantum transport calculations (integration with the quantum transport c...
subroutine, public create_transport_section(section)
creates the TRANSPORT section
function that build the XAS section of the input
subroutine, public create_xas_section(section)
makes the input section for core-level spectroscopy simulations
subroutine, public create_xas_tdp_section(section)
makes the input section for core-level spectroscopy simulations using linear response TDDFT
function that build the xc section of the input
subroutine, public create_xc_section(section)
creates the input section for the xc part
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, deprecation_notice)
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 lchar_t
integer, parameter, public logical_t
integer, parameter, public char_t
integer, parameter, public integer_t
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
Definition of physical constants:
Definition physcon.F:68
real(kind=dp), parameter, public femtoseconds
Definition physcon.F:153
real(kind=dp), parameter, public evolt
Definition physcon.F:183
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.
represent a keyword in the input
represent a section of the input file