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