(git:419edc0)
Loading...
Searching...
No Matches
input_cp2k_scf.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2025 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 scf section of the input
10!> \par History
11!> 10.2005 moved out of input_cp2k [fawzi]
12!> 07.2024 moved out of input_cp2k_dft [JGH]
13!> \author fawzi
14! **************************************************************************************************
16 USE bibliography, ONLY: becke1988b,&
28 USE cp_units, ONLY: cp_unit_to_cp2k
29 USE input_constants, ONLY: &
60 USE input_val_types, ONLY: integer_t,&
61 real_t
62 USE kinds, ONLY: dp
66 USE string_utilities, ONLY: newline,&
67 s2a
68#include "./base/base_uses.f90"
69
70 IMPLICIT NONE
71 PRIVATE
72
73 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'input_cp2k_scf'
74
76
77CONTAINS
78
79! **************************************************************************************************
80!> \brief creates the structure of the section with the DFT SCF parameters
81!> \param section will contain the SCF section
82!> \author fawzi
83! **************************************************************************************************
84 SUBROUTINE create_scf_section(section)
85 TYPE(section_type), POINTER :: section
86
87 TYPE(keyword_type), POINTER :: keyword
88 TYPE(section_type), POINTER :: print_key, subsection
89
90 NULLIFY (print_key)
91
92 cpassert(.NOT. ASSOCIATED(section))
93 CALL section_create(section, __location__, name="scf", &
94 description="Parameters needed to perform an SCF run.", &
95 n_keywords=18, n_subsections=7, repeats=.false.)
96
97 NULLIFY (subsection)
98
99 CALL create_ot_section(subsection)
100 CALL section_add_subsection(section, subsection)
101 CALL section_release(subsection)
102
103 CALL create_diagonalization_section(subsection)
104 CALL section_add_subsection(section, subsection)
105 CALL section_release(subsection)
106
107 CALL create_outer_scf_section(subsection)
108 CALL section_add_subsection(section, subsection)
109 CALL section_release(subsection)
110
111 CALL create_smear_section(subsection)
112 CALL section_add_subsection(section, subsection)
113 CALL section_release(subsection)
114
115 CALL create_mixing_section(subsection)
116 CALL section_add_subsection(section, subsection)
117 CALL section_release(subsection)
118
119 CALL create_mom_section(subsection)
120 CALL section_add_subsection(section, subsection)
121 CALL section_release(subsection)
122
123 NULLIFY (keyword)
124
125 CALL keyword_create(keyword, __location__, name="MAX_ITER_LUMO", &
126 variants=(/"MAX_ITER_LUMOS"/), &
127 description="Maximum number of iterations for the calculation of the LUMO energies "// &
128 "with the OT eigensolver.", &
129 usage="MAX_ITER_LUMO 100", default_i_val=299)
130 CALL section_add_keyword(section, keyword)
131 CALL keyword_release(keyword)
132
133 CALL keyword_create(keyword, __location__, name="EPS_LUMO", &
134 variants=(/"EPS_LUMOS"/), &
135 description="Target accuracy for the calculation of the LUMO energies with the OT eigensolver.", &
136 usage="EPS_LUMO 1.0E-6", default_r_val=1.0e-5_dp)
137 CALL section_add_keyword(section, keyword)
138 CALL keyword_release(keyword)
139
140 CALL keyword_create(keyword, __location__, name="MAX_SCF", &
141 description="Maximum number of SCF iteration to be performed for one optimization", &
142 usage="MAX_SCF 200", default_i_val=50)
143 CALL section_add_keyword(section, keyword)
144 CALL keyword_release(keyword)
145
146 CALL keyword_create(keyword, __location__, name="MAX_SCF_HISTORY", variants=(/"MAX_SCF_HIST"/), &
147 description="Maximum number of SCF iterations after the history pipeline is filled", &
148 usage="MAX_SCF_HISTORY 1", default_i_val=0, lone_keyword_i_val=1)
149 CALL section_add_keyword(section, keyword)
150 CALL keyword_release(keyword)
151
152 CALL keyword_create(keyword, __location__, name="MAX_DIIS", &
153 variants=(/"MAX_DIIS_BUFFER_SIZE"/), &
154 description="Maximum number of DIIS vectors to be used", &
155 usage="MAX_DIIS 3", default_i_val=4)
156 CALL section_add_keyword(section, keyword)
157 CALL keyword_release(keyword)
158
159 CALL keyword_create(keyword, __location__, name="LEVEL_SHIFT", &
160 variants=(/"LSHIFT"/), &
161 description="Use level shifting to improve convergence", &
162 unit_str="au_e", &
163 usage="LEVEL_SHIFT 0.1", &
164 default_r_val=0.0_dp)
165 CALL section_add_keyword(section, keyword)
166 CALL keyword_release(keyword)
167
168 CALL keyword_create(keyword, __location__, name="EPS_SCF", &
169 description="Target accuracy for the SCF convergence.", &
170 usage="EPS_SCF 1.e-6", default_r_val=1.e-5_dp)
171 CALL section_add_keyword(section, keyword)
172 CALL keyword_release(keyword)
173
174 CALL keyword_create(keyword, __location__, name="EPS_SCF_HISTORY", variants=(/"EPS_SCF_HIST"/), &
175 description="Target accuracy for the SCF convergence after the history pipeline is filled.", &
176 usage="EPS_SCF_HISTORY 1.e-5", default_r_val=0.0_dp, lone_keyword_r_val=1.0e-5_dp)
177 CALL section_add_keyword(section, keyword)
178 CALL keyword_release(keyword)
179
180 CALL keyword_create(keyword, __location__, name="CHOLESKY", &
181 description="If the cholesky method should be used for computing "// &
182 "the inverse of S, and in this case calling which Lapack routines", &
183 usage="CHOLESKY REDUCE", default_i_val=cholesky_restore, &
184 enum_c_vals=s2a("OFF", "REDUCE", "RESTORE", "INVERSE", "INVERSE_DBCSR"), &
185 enum_desc=s2a("The cholesky algorithm is not used", "Reduce is called", &
186 "Reduce is replaced by two restore", &
187 "Restore uses operator multiply by inverse of the triangular matrix", &
188 "Like inverse, but matrix stored as dbcsr, sparce matrix algebra used when possible"), &
190 CALL section_add_keyword(section, keyword)
191 CALL keyword_release(keyword)
192
193 CALL keyword_create(keyword, __location__, name="EPS_EIGVAL", &
194 description="Throw away linear combinations of basis functions with a small eigenvalue in S", &
195 usage="EPS_EIGVAL 1.0", default_r_val=1.0e-5_dp)
196 CALL section_add_keyword(section, keyword)
197 CALL keyword_release(keyword)
198
199 CALL keyword_create(keyword, __location__, name="EPS_DIIS", &
200 description="Threshold on the convergence to start using DIAG/DIIS or OT/DIIS."// &
201 " Default for OT/DIIS is never to switch.", &
202 usage="EPS_DIIS 5.0e-2", default_r_val=0.1_dp)
203 CALL section_add_keyword(section, keyword)
204 CALL keyword_release(keyword)
205
206 CALL keyword_create( &
207 keyword, __location__, name="SCF_GUESS", &
208 description="Change the initial guess for the wavefunction.", &
209 usage="SCF_GUESS RESTART", default_i_val=atomic_guess, &
210 enum_c_vals=s2a("ATOMIC", "RESTART", "RANDOM", "CORE", &
211 "HISTORY_RESTART", "MOPAC", "EHT", "SPARSE", "NONE"), &
212 enum_desc=s2a("Generate an atomic density using the atomic code", &
213 "Use the RESTART file as an initial guess (and ATOMIC if not present).", &
214 "Use random wavefunction coefficients.", &
215 "Diagonalize the core hamiltonian for an initial guess.", &
216 "Extrapolated from previous RESTART files.", &
217 "Use same guess as MOPAC for semi-empirical methods or a simple diagonal density matrix for other methods", &
218 "Use the EHT (gfn0-xTB) code to generate an initial wavefunction.", &
219 "Generate a sparse wavefunction using the atomic code (for OT based methods)", &
220 "Skip initial guess (only for non-self consistent methods)."), &
223 CALL section_add_keyword(section, keyword)
224 CALL keyword_release(keyword)
225
226 CALL keyword_create(keyword, __location__, name="NROW_BLOCK", &
227 description="sets the number of rows in a scalapack block", &
228 usage="NROW_BLOCK 31", default_i_val=32)
229 CALL section_add_keyword(section, keyword)
230 CALL keyword_release(keyword)
231
232 CALL keyword_create(keyword, __location__, name="NCOL_BLOCK", &
233 description="Sets the number of columns in a scalapack block", &
234 usage="NCOL_BLOCK 31", default_i_val=32)
235 CALL section_add_keyword(section, keyword)
236 CALL keyword_release(keyword)
237
238 CALL keyword_create(keyword, __location__, name="ADDED_MOS", &
239 description="Number of additional MOS added for each spin. Use -1 to add all available. "// &
240 "alpha/beta spin can be specified independently (if spin-polarized calculation requested).", &
241 usage="ADDED_MOS", default_i_val=0, n_var=-1)
242 CALL section_add_keyword(section, keyword)
243 CALL keyword_release(keyword)
244
245 CALL keyword_create(keyword, __location__, &
246 name="ROKS_SCHEME", &
247 description="Selects the ROKS scheme when ROKS is applied.", &
248 usage="ROKS_SCHEME HIGH-SPIN", &
249 repeats=.false., &
250 n_var=1, &
251 enum_c_vals=s2a("GENERAL", "HIGH-SPIN"), &
252 enum_i_vals=(/general_roks, high_spin_roks/), &
253 default_i_val=high_spin_roks)
254 CALL section_add_keyword(section, keyword)
255 CALL keyword_release(keyword)
256
257 CALL keyword_create(keyword, __location__, &
258 name="ROKS_F", &
259 variants=(/"F_ROKS"/), &
260 description="Allows to define the parameter f for the "// &
261 "general ROKS scheme.", &
262 usage="ROKS_F 1/2", &
263 repeats=.false., &
264 n_var=1, &
265 type_of_var=real_t, &
266 default_r_val=0.5_dp)
267 CALL section_add_keyword(section, keyword)
268 CALL keyword_release(keyword)
269
270 CALL keyword_create(keyword, __location__, &
271 name="ROKS_PARAMETERS", &
272 variants=(/"ROKS_PARAMETER"/), &
273 description="Allows to define all parameters for the high-spin "// &
274 "ROKS scheme explicitly. "// &
275 "The full set of 6 parameters has to be specified "// &
276 "in the order acc, bcc, aoo, boo, avv, bvv", &
277 usage="ROKS_PARAMETERS 1/2 1/2 1/2 1/2 1/2 1/2", &
278 repeats=.false., &
279 n_var=6, &
280 type_of_var=real_t, &
281 default_r_vals=(/-0.5_dp, 1.5_dp, 0.5_dp, 0.5_dp, 1.5_dp, -0.5_dp/))
282 CALL section_add_keyword(section, keyword)
283 CALL keyword_release(keyword)
284
285 CALL keyword_create(keyword, __location__, name="IGNORE_CONVERGENCE_FAILURE", &
286 description="If true, only a warning is issued if an SCF "// &
287 "iteration has not converged. By default, a run is aborted "// &
288 "if the required convergence criteria have not been achieved.", &
289 usage="IGNORE_CONVERGENCE_FAILURE logical_value", &
290 default_l_val=.false., &
291 lone_keyword_l_val=.true.)
292 CALL section_add_keyword(section, keyword)
293 CALL keyword_release(keyword)
294
295 CALL keyword_create(keyword, __location__, name="FORCE_SCF_CALCULATION", &
296 description="Request a SCF type solution even for nonSCF methods. ", &
297 usage="FORCE_SCF_CALCULATION logical_value", &
298 default_l_val=.false., &
299 lone_keyword_l_val=.true.)
300 CALL section_add_keyword(section, keyword)
301 CALL keyword_release(keyword)
302
303 CALL section_create(subsection, __location__, name="PRINT", &
304 description="Printing of information during the SCF.", repeats=.false.)
305
306 CALL cp_print_key_section_create(print_key, __location__, "RESTART", &
307 description="Controls the dumping of the MO restart file during SCF. "// &
308 "By default keeps a short history of three restarts. "// &
309 "See also RESTART_HISTORY", &
310 print_level=low_print_level, common_iter_levels=3, &
311 each_iter_names=s2a("QS_SCF"), each_iter_values=(/20/), &
312 add_last=add_last_numeric, filename="RESTART")
313 CALL keyword_create(keyword, __location__, name="BACKUP_COPIES", &
314 description="Specifies the maximum number of backup copies.", &
315 usage="BACKUP_COPIES {int}", &
316 default_i_val=1)
317 CALL section_add_keyword(print_key, keyword)
318 CALL keyword_release(keyword)
319 CALL section_add_subsection(subsection, print_key)
320 CALL section_release(print_key)
321
323 print_key, __location__, "RESTART_HISTORY", &
324 description="Dumps unique MO restart files during the run keeping all of them.", &
325 print_level=low_print_level, common_iter_levels=0, &
326 each_iter_names=s2a("__ROOT__", "MD", "GEO_OPT", "ROT_OPT", "NEB", "METADYNAMICS", "QS_SCF"), &
327 each_iter_values=(/500, 500, 500, 500, 500, 500, 500/), &
328 filename="RESTART")
329 CALL keyword_create(keyword, __location__, name="BACKUP_COPIES", &
330 description="Specifies the maximum number of backup copies.", &
331 usage="BACKUP_COPIES {int}", &
332 default_i_val=1)
333 CALL section_add_keyword(print_key, keyword)
334 CALL keyword_release(keyword)
335 CALL section_add_subsection(subsection, print_key)
336 CALL section_release(print_key)
337
338 CALL cp_print_key_section_create(print_key, __location__, "iteration_info", &
339 description="Controls the printing of basic iteration information during the SCF.", &
340 print_level=low_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
341 CALL keyword_create(keyword, __location__, name="time_cumul", &
342 description="If the printkey is activated switches the printing of timings"// &
343 " to cumulative (over the SCF).", &
344 default_l_val=.false., lone_keyword_l_val=.true.)
345 CALL section_add_keyword(print_key, keyword)
346 CALL keyword_release(keyword)
347 CALL section_add_subsection(subsection, print_key)
348 CALL section_release(print_key)
349
350 CALL cp_print_key_section_create(print_key, __location__, "program_run_info", &
351 description="Controls the printing of basic information during the SCF.", &
352 print_level=low_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
353 CALL section_add_subsection(subsection, print_key)
354 CALL section_release(print_key)
355
356 CALL cp_print_key_section_create(print_key, __location__, "MO_ORTHONORMALITY", &
357 description="Controls the printing relative to the orthonormality of MOs (CT S C).", &
358 print_level=high_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
359 CALL section_add_subsection(subsection, print_key)
360 CALL section_release(print_key)
361
362 CALL cp_print_key_section_create(print_key, __location__, "MO_MAGNITUDE", &
363 description="Prints the min/max eigenvalues of the overlap of the MOs without S (CT C).", &
364 print_level=high_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
365 CALL section_add_subsection(subsection, print_key)
366 CALL section_release(print_key)
367
368 CALL cp_print_key_section_create(print_key, __location__, "detailed_energy", &
369 description="Controls the printing of detailed energy information.", &
370 print_level=high_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
371 CALL section_add_subsection(subsection, print_key)
372 CALL section_release(print_key)
373
374 CALL cp_print_key_section_create(print_key, __location__, "diis_info", &
375 description="Controls the printing of DIIS information.", &
376 print_level=high_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
377 CALL section_add_subsection(subsection, print_key)
378 CALL section_release(print_key)
379
380 CALL cp_print_key_section_create(print_key, __location__, "total_densities", &
381 description="Controls the printing of total densities.", &
382 print_level=high_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
383 CALL section_add_subsection(subsection, print_key)
384 CALL section_release(print_key)
385
386 CALL cp_print_key_section_create(print_key, __location__, "Lanczos", &
387 description="Controls the printing of information on Lanczos refinement iterations.", &
388 print_level=high_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
389 CALL section_add_subsection(subsection, print_key)
390 CALL section_release(print_key)
391
393 print_key, __location__, "DIAG_SUB_SCF", &
394 description="Controls the printing of information on subspace diagonalization internal loop. ", &
395 print_level=high_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
396 CALL section_add_subsection(subsection, print_key)
397 CALL section_release(print_key)
398
399 CALL cp_print_key_section_create(print_key, __location__, "Davidson", &
400 description="Controls the printing of information on Davidson iterations.", &
401 print_level=high_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
402 CALL section_add_subsection(subsection, print_key)
403 CALL section_release(print_key)
404
405 CALL cp_print_key_section_create(print_key, __location__, "FILTER_MATRIX", &
406 description="Controls the printing of information on Filter Matrix method.", &
407 print_level=high_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
408 CALL section_add_subsection(subsection, print_key)
409 CALL section_release(print_key)
410
411 CALL keyword_create(keyword, __location__, name="DM_RESTART_WRITE", &
412 description="Write the density matrix into a binary file at the end of the SCF.", &
413 usage="DM_RESTART_WRITE", default_l_val=.false., lone_keyword_l_val=.true.)
414 CALL section_add_keyword(subsection, keyword)
415 CALL keyword_release(keyword)
416
417 CALL section_add_subsection(section, subsection)
418 CALL section_release(subsection)
419
420 END SUBROUTINE create_scf_section
421
422! **************************************************************************************************
423!> \brief creates the structure of the section with SCF parameters
424!> controlling an other loop
425!> \param section will contain the SCF section
426!> \author Joost VandeVondele [2006.03]
427! **************************************************************************************************
428 SUBROUTINE create_outer_scf_section(section)
429 TYPE(section_type), POINTER :: section
430
431 TYPE(keyword_type), POINTER :: keyword
432 TYPE(section_type), POINTER :: subsection
433
434 cpassert(.NOT. ASSOCIATED(section))
435 CALL section_create(section, __location__, name="OUTER_SCF", &
436 description="parameters controlling the outer SCF loop", &
437 n_keywords=13, n_subsections=1, repeats=.false.)
438
439 NULLIFY (keyword)
440
441 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
442 description="controls the activation of the outer SCF loop", &
443 usage="&OUTER_SCF ON", default_l_val=.false., lone_keyword_l_val=.true.)
444 CALL section_add_keyword(section, keyword)
445 CALL keyword_release(keyword)
446
447 ! add CDFT_OPT section
448 NULLIFY (subsection)
449 CALL create_cdft_opt_section(subsection)
450 CALL section_add_subsection(section, subsection)
451 CALL section_release(subsection)
452
453 CALL keyword_create(keyword, __location__, name="TYPE", &
454 description="Specifies which kind of outer SCF should be employed", &
455 usage="TYPE DDAPC_CONSTRAINT ", &
456 default_i_val=outer_scf_none, &
457 enum_c_vals=s2a("DDAPC_CONSTRAINT", "S2_CONSTRAINT", &
458 "BASIS_CENTER_OPT", "CDFT_CONSTRAINT", "NONE"), &
459 enum_desc=s2a("Enforce a constraint on the DDAPC, requires the corresponding section", &
460 "Enforce a constraint on the S2, requires the corresponding section", &
461 "Optimize positions of basis functions, if atom types FLOATING_BASIS_CENTER "// &
462 "are defined", &
463 "Enforce a constraint on a generic CDFT weight population. "// &
464 "Requires the corresponding section QS&CDFT"// &
465 " which determines the type of weight used.", &
466 "Do nothing in the outer loop, useful for resetting the inner loop,"), &
469 CALL section_add_keyword(section, keyword)
470 CALL keyword_release(keyword)
471
472 CALL keyword_create(keyword, __location__, name="OPTIMIZER", &
473 description="Method used to bring the outer loop to a stationary point", &
474 usage="OPTIMIZER SD", &
475 default_i_val=outer_scf_optimizer_none, &
476 enum_c_vals=s2a("SD", "DIIS", "NONE", "BISECT", "BROYDEN", "NEWTON", "SECANT", "NEWTON_LS"), &
477 enum_desc=s2a("Takes steps in the direction of the gradient, multiplied by step_size", &
478 "Uses a Direct Inversion in the Iterative Subspace method", &
479 "Do nothing, useful only with the none type", &
480 "Bisection of the gradient, useful for difficult one dimensional cases", &
481 "Broyden's method. Variant defined in BROYDEN_TYPE.", &
482 "Newton's method. Only compatible with CDFT constraints.", &
483 "Secant method. Only for one dimensional cases. See Broyden for "// &
484 "multidimensional cases.", &
485 "Newton's method with backtracking line search to find the optimal step size. "// &
486 "Only compatible with CDFT constraints. Starts from the regular Newton solution "// &
487 "and successively reduces the step size until the L2 norm of the CDFT gradient "// &
488 "decreases or MAX_LS steps is reached. Potentially very expensive because "// &
489 "each iteration performs a full SCF calculation."), &
494 CALL section_add_keyword(section, keyword)
495 CALL keyword_release(keyword)
496
497 CALL keyword_create(keyword, __location__, name="BISECT_TRUST_COUNT", &
498 description="Maximum number of times the same point will be used in bisection,"// &
499 " a small number guards against the effect of wrongly converged states.", &
500 usage="BISECT_TRUST_COUNT 5", default_i_val=10)
501 CALL section_add_keyword(section, keyword)
502 CALL keyword_release(keyword)
503
504 CALL keyword_create(keyword, __location__, name="EPS_SCF", &
505 description="The target gradient of the outer SCF variables. "// &
506 "Notice that the EPS_SCF of the inner loop also determines "// &
507 "the value that can be reached in the outer loop, "// &
508 "typically EPS_SCF of the outer loop must be smaller "// &
509 "than or equal to EPS_SCF of the inner loop.", &
510 usage="EPS_SCF 1.0E-6 ", default_r_val=1.0e-5_dp)
511 CALL section_add_keyword(section, keyword)
512 CALL keyword_release(keyword)
513
514 CALL keyword_create(keyword, __location__, name="DIIS_BUFFER_LENGTH", &
515 description="Maximum number of DIIS vectors used ", &
516 usage="DIIS_BUFFER_LENGTH 5", default_i_val=3)
517 CALL section_add_keyword(section, keyword)
518 CALL keyword_release(keyword)
519
520 CALL keyword_create(keyword, __location__, name="EXTRAPOLATION_ORDER", &
521 description="Number of past states used in the extrapolation of the variables during e.g. MD", &
522 usage="EXTRAPOLATION_ORDER 5", default_i_val=3)
523 CALL section_add_keyword(section, keyword)
524 CALL keyword_release(keyword)
525
526 CALL keyword_create(keyword, __location__, name="MAX_SCF", &
527 description="The maximum number of outer loops ", &
528 usage="MAX_SCF 20", default_i_val=50)
529 CALL section_add_keyword(section, keyword)
530 CALL keyword_release(keyword)
531
532 CALL keyword_create(keyword, __location__, name="STEP_SIZE", &
533 description="The initial step_size used in the optimizer (currently steepest descent). "// &
534 "Note that in cases where a sadle point is sought for (constrained DFT),"// &
535 " this can be negative. For Newton and Broyden optimizers, use a value less/higher than "// &
536 "the default 1.0 (in absolute value, the sign is not significant) to active an under/overrelaxed "// &
537 "optimizer.", &
538 usage="STEP_SIZE -1.0", default_r_val=0.5_dp)
539 CALL section_add_keyword(section, keyword)
540 CALL keyword_release(keyword)
541
542 END SUBROUTINE create_outer_scf_section
543
544! **************************************************************************************************
545!> \brief makes the orbital transformation section
546!> \param section ...
547!> \par History
548!> 11.2004 created [Joost VandeVondele]
549! **************************************************************************************************
550 SUBROUTINE create_ot_section(section)
551 TYPE(section_type), POINTER :: section
552
553 TYPE(keyword_type), POINTER :: keyword
554
555 cpassert(.NOT. ASSOCIATED(section))
556 CALL section_create(section, __location__, name="OT", &
557 description="Sets the various options for the orbital transformation (OT) method. "// &
558 "Default settings already provide an efficient, yet robust method. "// &
559 "Most systems benefit from using the FULL_ALL preconditioner "// &
560 "combined with a small value (0.001) of ENERGY_GAP. "// &
561 "Well-behaved systems might benefit from using a DIIS minimizer. "//newline//newline// &
562 "**Advantages:** "// &
563 "It's fast, because no expensive diagonalisation is performed. "// &
564 "If preconditioned correctly, method guaranteed to find minimum. "//newline//newline// &
565 "**Disadvantages:** "// &
566 "Sensitive to preconditioning. A good preconditioner can be expensive. "// &
567 "No smearing, or advanced SCF mixing possible: POOR convergence for metallic systems.", &
568 n_keywords=27, n_subsections=0, repeats=.false., &
569 citations=(/vandevondele2003, weber2008/))
570
571 NULLIFY (keyword)
572
573 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
574 description="controls the activation of the ot method", &
575 usage="&OT T", &
576 default_l_val=.false., &
577 lone_keyword_l_val=.true.)
578 CALL section_add_keyword(section, keyword)
579 CALL keyword_release(keyword)
580
581 CALL keyword_create(keyword, __location__, name="ALGORITHM", &
582 description="Algorithm to be used for OT", &
583 usage="ALGORITHM STRICT", &
584 default_i_val=ot_algo_taylor_or_diag, &
585 enum_c_vals=s2a("STRICT", "IRAC"), &
586 enum_desc=s2a("Strict orthogonality: Taylor or diagonalization based algorithm.", &
587 "Orbital Transformation based Iterative Refinement "// &
588 "of the Approximative Congruence transformation (OT/IR)."), &
589 enum_i_vals=(/ot_algo_taylor_or_diag, ot_algo_irac/), &
591 CALL section_add_keyword(section, keyword)
592 CALL keyword_release(keyword)
593
594 CALL keyword_create(keyword, __location__, name="IRAC_DEGREE", &
595 description="The refinement polynomial degree (2, 3 or 4).", &
596 usage="IRAC_DEGREE 4", &
597 default_i_val=4)
598 CALL section_add_keyword(section, keyword)
599 CALL keyword_release(keyword)
600
601 CALL keyword_create(keyword, __location__, name="MAX_IRAC", &
602 description="Maximum allowed refinement iteration.", &
603 usage="MAX_IRAC 5", &
604 default_i_val=50)
605 CALL section_add_keyword(section, keyword)
606 CALL keyword_release(keyword)
607
608 CALL keyword_create(keyword, __location__, name="ORTHO_IRAC", &
609 description="The orthogonality method.", &
610 usage="ORTHO_IRAC POLY", &
611 default_i_val=ot_chol_irac, &
612 enum_c_vals=s2a("CHOL", "POLY", "LWDN"), &
613 enum_desc=s2a("Cholesky.", "Polynomial.", "Loewdin."), &
614 enum_i_vals=(/ot_chol_irac, ot_poly_irac, ot_lwdn_irac/))
615 CALL section_add_keyword(section, keyword)
616 CALL keyword_release(keyword)
617
618 CALL keyword_create(keyword, __location__, name="EPS_IRAC_FILTER_MATRIX", &
619 description="Sets the threshold for filtering the matrices.", &
620 usage="EPS_IRAC_FILTER_MATRIX 1.0E-5", &
621 default_r_val=0.0_dp)
622 CALL section_add_keyword(section, keyword)
623 CALL keyword_release(keyword)
624
625 CALL keyword_create(keyword, __location__, name="EPS_IRAC", &
626 description="Targeted accuracy during the refinement iteration.", &
627 usage="EPS_IRAC 1.0E-5", &
628 default_r_val=1.0e-10_dp)
629 CALL section_add_keyword(section, keyword)
630 CALL keyword_release(keyword)
631
632 CALL keyword_create(keyword, __location__, name="EPS_IRAC_QUICK_EXIT", &
633 description="Only one extra refinement iteration is "// &
634 "done when the norm is below this value.", &
635 usage="EPS_IRAC_QUICK_EXIT 1.0E-2", &
636 default_r_val=1.0e-5_dp)
637 CALL section_add_keyword(section, keyword)
638 CALL keyword_release(keyword)
639
640 CALL keyword_create(keyword, __location__, name="EPS_IRAC_SWITCH", &
641 description="The algorithm switches to the polynomial "// &
642 "refinement when the norm is below this value.", &
643 usage="EPS_IRAC_SWITCH 1.0E-3", &
644 default_r_val=1.0e-2_dp)
645 CALL section_add_keyword(section, keyword)
646 CALL keyword_release(keyword)
647
648 CALL keyword_create(keyword, __location__, name="ON_THE_FLY_LOC", &
649 description="On the fly localization of the molecular orbitals. "// &
650 "Can only be used with OT/IRAC.", &
651 usage="ON_THE_FLY_LOC T", &
652 default_l_val=.false.)
653 CALL section_add_keyword(section, keyword)
654 CALL keyword_release(keyword)
655
656 CALL keyword_create( &
657 keyword, __location__, name="MINIMIZER", &
658 description="Minimizer to be used with the OT method", &
659 usage="MINIMIZER DIIS", &
660 default_i_val=ot_mini_cg, &
661 enum_c_vals=s2a("SD", "CG", "DIIS", "BROYDEN"), &
662 enum_desc=s2a("Steepest descent: not recommended", "Conjugate Gradients: most reliable, use for difficult systems."// &
663 " The total energy should decrease at every OT CG step if the line search is appropriate.", &
664 "Direct inversion in the iterative subspace: less reliable than CG, but sometimes about 50% faster", &
665 "Broyden mixing approximating the inverse Hessian"), &
667 CALL section_add_keyword(section, keyword)
668 CALL keyword_release(keyword)
669
670 CALL keyword_create(keyword, __location__, name="SAFE_DIIS", &
671 variants=(/"SAFER_DIIS"/), &
672 description="Reject DIIS steps if they point away from the"// &
673 " minimum, do SD in that case.", &
674 usage="SAFE_DIIS ON", default_l_val=.true.)
675 CALL section_add_keyword(section, keyword)
676 CALL keyword_release(keyword)
677
678 CALL keyword_create(keyword, __location__, name="MAX_SCF_DIIS", &
679 description="Maximum DIIS SCF inner loop cycles. This can be used to extend"// &
680 " SCF cycles after a switch to DIIS (see eps_diis).", &
681 usage="MAX_SCF_DIIS 20", &
682 default_i_val=0)
683 CALL section_add_keyword(section, keyword)
684 CALL keyword_release(keyword)
685
686 CALL keyword_create(keyword, __location__, name="N_HISTORY_VEC", &
687 variants=s2a("NDIIS", "N_DIIS", "N_BROYDEN"), &
688 description="Number of history vectors to be used with DIIS or BROYDEN", &
689 usage="N_DIIS 4", &
690 default_i_val=7)
691 CALL section_add_keyword(section, keyword)
692 CALL keyword_release(keyword)
693
694 CALL keyword_create(keyword, __location__, name="BROYDEN_BETA", &
695 description="Underrelaxation for the broyden mixer", &
696 usage="BROYDEN_BETA 0.9", &
697 default_r_val=0.9_dp)
698 CALL section_add_keyword(section, keyword)
699 CALL keyword_release(keyword)
700
701 CALL keyword_create(keyword, __location__, name="BROYDEN_GAMMA", &
702 description="Backtracking parameter", &
703 usage="BROYDEN_GAMMA 0.5", &
704 default_r_val=0.5_dp)
705 CALL section_add_keyword(section, keyword)
706 CALL keyword_release(keyword)
707
708 CALL keyword_create(keyword, __location__, name="BROYDEN_SIGMA", &
709 description="Curvature of energy functional.", &
710 usage="BROYDEN_SIGMA 0.25", &
711 default_r_val=0.25_dp)
712 CALL section_add_keyword(section, keyword)
713 CALL keyword_release(keyword)
714
715 CALL keyword_create(keyword, __location__, name="BROYDEN_ETA", &
716 description="Dampening of estimated energy curvature.", &
717 usage="BROYDEN_ETA 0.7", &
718 default_r_val=0.7_dp)
719 CALL section_add_keyword(section, keyword)
720 CALL keyword_release(keyword)
721
722 CALL keyword_create(keyword, __location__, name="BROYDEN_OMEGA", &
723 description="Growth limit of curvature.", &
724 usage="BROYDEN_OMEGA 1.1", &
725 default_r_val=1.1_dp)
726 CALL section_add_keyword(section, keyword)
727 CALL keyword_release(keyword)
728
729 CALL keyword_create(keyword, __location__, name="BROYDEN_SIGMA_DECREASE", &
730 description="Reduction of curvature on bad approximation.", &
731 usage="BROYDEN_SIGMA_DECREASE 0.7", &
732 default_r_val=0.7_dp)
733 CALL section_add_keyword(section, keyword)
734 CALL keyword_release(keyword)
735
736 CALL keyword_create(keyword, __location__, name="BROYDEN_SIGMA_MIN", &
737 description="Minimum adaptive curvature.", &
738 usage="BROYDEN_SIGMA_MIN 0.05", &
739 default_r_val=0.05_dp)
740 CALL section_add_keyword(section, keyword)
741 CALL keyword_release(keyword)
742
743 CALL keyword_create(keyword, __location__, name="BROYDEN_FORGET_HISTORY", &
744 description="Forget history on bad approximation", &
745 usage="BROYDEN_FORGET_HISTORY OFF", default_l_val=.false., &
746 lone_keyword_l_val=.true.)
747 CALL section_add_keyword(section, keyword)
748 CALL keyword_release(keyword)
749
750 CALL keyword_create(keyword, __location__, name="BROYDEN_ADAPTIVE_SIGMA", &
751 description="Enable adaptive curvature estimation", &
752 usage="BROYDEN_ADAPTIVE_SIGMA ON", default_l_val=.true., &
753 lone_keyword_l_val=.true.)
754 CALL section_add_keyword(section, keyword)
755 CALL keyword_release(keyword)
756
757 CALL keyword_create(keyword, __location__, name="BROYDEN_ENABLE_FLIP", &
758 description="Ensure positive definite update", &
759 usage="BROYDEN_ENABLE_FLIP ON", default_l_val=.true., &
760 lone_keyword_l_val=.true.)
761 CALL section_add_keyword(section, keyword)
762 CALL keyword_release(keyword)
763
764 CALL keyword_create(keyword, __location__, name="LINESEARCH", &
765 variants=(/"LINE_SEARCH"/), &
766 description="1D line search algorithm to be used with the OT minimizer,"// &
767 " in increasing order of robustness and cost. MINIMIZER CG combined with"// &
768 " LINESEARCH GOLD should always find an electronic minimum."// &
769 " Whereas the 2PNT minimizer is almost always OK, 3PNT might be needed for systems"// &
770 " in which successive OT CG steps do not decrease the total energy.", &
771 usage="LINESEARCH GOLD", &
772 default_i_val=ls_2pnt, &
773 enum_c_vals=s2a("ADAPT", "NONE", "2PNT", "3PNT", "GOLD"), &
774 enum_desc=s2a("extrapolates usually based on 3 points,"// &
775 " uses additional points on demand, very robust.", &
776 "always take steps of fixed length", &
777 "extrapolate based on 3 points", &
778 "extrapolate based on 2 points", &
779 "perform 1D golden section search of the minimum (very expensive)"), &
780 enum_i_vals=(/ls_adapt, ls_none, ls_2pnt, ls_3pnt, ls_gold/))
781 CALL section_add_keyword(section, keyword)
782 CALL keyword_release(keyword)
783
784 CALL keyword_create( &
785 keyword, __location__, name="STEPSIZE", &
786 description="Initial stepsize used for the line search, sometimes this parameter can be reduced to stabilize DIIS"// &
787 " or to improve the CG behavior in the first few steps."// &
788 " The optimal value depends on the quality of the preconditioner."// &
789 " A negative values leaves the choice to CP2K depending on the preconditioner.", &
790 usage="STEPSIZE 0.4", &
791 default_r_val=-1.0_dp)
792 CALL section_add_keyword(section, keyword)
793 CALL keyword_release(keyword)
794
795 CALL keyword_create(keyword, __location__, name="GOLD_TARGET", &
796 description="Target relative uncertainty in the location of the minimum for LINESEARCH GOLD", &
797 usage="GOLD_TARGET 0.1", &
798 default_r_val=0.01_dp)
799 CALL section_add_keyword(section, keyword)
800 CALL keyword_release(keyword)
801
802 CALL keyword_create( &
803 keyword, __location__, name="PRECONDITIONER", &
804 description="Type of preconditioner to be used with all minimization schemes. "// &
805 "They differ in effectiveness, cost of construction, cost of application. "// &
806 "Properly preconditioned minimization can be orders of magnitude faster than doing nothing.", &
807 usage="PRECONDITIONER FULL_ALL", &
808 default_i_val=ot_precond_full_kinetic, &
809 enum_c_vals=s2a("FULL_ALL", "FULL_SINGLE_INVERSE", "FULL_SINGLE", "FULL_KINETIC", "FULL_S_INVERSE", &
810 "NONE"), &
811 enum_desc=s2a("Most effective state selective preconditioner based on diagonalization, "// &
812 "requires the ENERGY_GAP parameter to be an underestimate of the HOMO-LUMO gap. "// &
813 "This preconditioner is recommended for almost all systems, except very large systems where "// &
814 "make_preconditioner would dominate the total computational cost.", &
815 "Based on H-eS cholesky inversion, similar to FULL_SINGLE in preconditioning efficiency "// &
816 "but cheaper to construct, "// &
817 "might be somewhat less robust. Recommended for large systems.", &
818 "Based on H-eS diagonalisation, not as good as FULL_ALL, but somewhat cheaper to apply. ", &
819 "Cholesky inversion of S and T, fast construction, robust, and relatively good, "// &
820 "use for very large systems.", &
821 "Cholesky inversion of S, not as good as FULL_KINETIC, yet equally expensive.", &
822 "skip preconditioning"), &
826 CALL section_add_keyword(section, keyword)
827 CALL keyword_release(keyword)
828
829 CALL keyword_create(keyword, __location__, name="CHOLESKY", &
830 description="If FULL_ALL the cholesky decomposition of the S matrix is used. "// &
831 "Options on the algorithm to be used.", &
832 usage="CHOLESKY REDUCE", default_i_val=cholesky_reduce, &
833 enum_c_vals=s2a("OFF", "REDUCE", "RESTORE", "INVERSE", "INVERSE_DBCSR"), &
834 enum_desc=s2a("The cholesky algorithm is not used", "Reduce is called", &
835 "Reduce is replaced by two restore", &
836 "Restore uses operator multiply by inverse of the triangular matrix", &
837 "Like inverse, but matrix stored as dbcsr, sparce matrix algebra used when possible"), &
839 CALL section_add_keyword(section, keyword)
840 CALL keyword_release(keyword)
841
842 CALL keyword_create( &
843 keyword, __location__, name="PRECOND_SOLVER", &
844 description="How the preconditioner is applied to the residual.", &
845 usage="PRECOND_SOLVER DIRECT", &
846 default_i_val=ot_precond_solver_default, &
847 enum_c_vals=s2a("DEFAULT", "DIRECT", "INVERSE_CHOLESKY", "INVERSE_UPDATE"), &
848 enum_desc=s2a("the default", "Cholesky decomposition followed by triangular solve "// &
849 "(works for FULL_KINETIC/SINGLE_INVERSE/S_INVERSE)", &
850 "Cholesky decomposition followed by explicit inversion "// &
851 "(works for FULL_KINETIC/SINGLE_INVERSE/S_INVERSE)", &
852 "Performs a Hotelling update of the inverse if a previous preconditioner is present. "// &
853 "Mainly useful for GPU accelerated systems (works for FULL_KINETIC/SINGLE_INVERSE/S_INVERSE)"), &
854 enum_i_vals=(/ot_precond_solver_default, &
858 CALL section_add_keyword(section, keyword)
859 CALL keyword_release(keyword)
860
861 CALL keyword_create( &
862 keyword, __location__, name="ENERGY_GAP", &
863 description="Should be an estimate for the energy gap [a.u.] (HOMO-LUMO) and is used in preconditioning, "// &
864 "especially effective with the FULL_ALL preconditioner, in which case it should be an underestimate "// &
865 "of the gap (can be a small number, e.g. 0.002)."// &
866 " FULL_SINGLE_INVERSE takes it as lower bound (values below 0.05 can cause stability issues)."// &
867 " In general, higher values will tame the preconditioner in case of poor initial guesses."// &
868 " A negative value will leave the choice to CP2K depending on type of preconditioner.", &
869 usage="ENERGY_GAP 0.001", &
870 default_r_val=-1.0_dp)
871 CALL section_add_keyword(section, keyword)
872 CALL keyword_release(keyword)
873
874 CALL keyword_create( &
875 keyword, __location__, name="EPS_TAYLOR", &
876 variants=(/"EPSTAYLOR"/), &
877 description="Target accuracy of the taylor expansion for the matrix functions, should normally be kept as is.", &
878 usage="EPS_TAYLOR 1.0E-15", &
879 default_r_val=1.0e-16_dp)
880 CALL section_add_keyword(section, keyword)
881 CALL keyword_release(keyword)
882
883 CALL keyword_create( &
884 keyword, __location__, name="MAX_TAYLOR", &
885 description="Maximum order of the Taylor expansion before diagonalisation is preferred, for large parallel runs"// &
886 " a slightly higher order could sometimes result in a small speedup.", &
887 usage="MAX_TAYLOR 5", &
888 default_i_val=4)
889 CALL section_add_keyword(section, keyword)
890 CALL keyword_release(keyword)
891
892 CALL keyword_create(keyword, __location__, name="ROTATION", &
893 description="Introduce additional variables so that rotations of the occupied"// &
894 " subspace are allowed as well, only needed for cases where the energy is not invariant under"// &
895 " a rotation of the occupied subspace such as non-singlet restricted calculations"// &
896 " or fractional occupations.", &
897 usage="ROTATION", lone_keyword_l_val=.true., &
898 default_l_val=.false.)
899 CALL section_add_keyword(section, keyword)
900 CALL keyword_release(keyword)
901
902 CALL keyword_create(keyword, __location__, name="ENERGIES", &
903 description="Optimize orbital energies for use in Fermi-Dirac smearing "// &
904 "(requires ROTATION and FD smearing to be active).", &
905 usage="ENERGIES", lone_keyword_l_val=.true., &
906 default_l_val=.false.)
907 CALL section_add_keyword(section, keyword)
908 CALL keyword_release(keyword)
909
910 CALL keyword_create(keyword, __location__, name="OCCUPATION_PRECONDITIONER", &
911 description="Preconditioner with the occupation numbers (FD smearing)", &
912 usage="OCCUPATION_PRECONDITIONER", lone_keyword_l_val=.true., &
913 default_l_val=.false.)
914 CALL section_add_keyword(section, keyword)
915 CALL keyword_release(keyword)
916
917 CALL keyword_create(keyword, __location__, name="NONDIAG_ENERGY", &
918 description="Add a non-diagonal energy penalty (FD smearing)", &
919 usage="NONDIAG_ENERGY", lone_keyword_l_val=.true., &
920 default_l_val=.false.)
921 CALL section_add_keyword(section, keyword)
922 CALL keyword_release(keyword)
923
924 CALL keyword_create(keyword, __location__, name="NONDIAG_ENERGY_STRENGTH", &
925 description="The prefactor for the non-diagonal energy penalty (FD smearing)", &
926 usage="NONDIAG_ENERGY_STRENGTH", default_r_val=1.0_dp)
927 CALL section_add_keyword(section, keyword)
928 CALL keyword_release(keyword)
929
930 END SUBROUTINE create_ot_section
931
932! **************************************************************************************************
933!> \brief creates the diagonalization section
934!> \param section ...
935!> \par History
936!> 10.2008 created [JGH]
937! **************************************************************************************************
938 SUBROUTINE create_diagonalization_section(section)
939 TYPE(section_type), POINTER :: section
940
941 TYPE(keyword_type), POINTER :: keyword
942 TYPE(section_type), POINTER :: subsection
943
944 cpassert(.NOT. ASSOCIATED(section))
945 CALL section_create(section, __location__, name="DIAGONALIZATION", &
946 description="Set up type and parameters for Kohn-Sham matrix diagonalization.", &
947 n_keywords=0, n_subsections=1, repeats=.false.)
948
949 NULLIFY (keyword)
950
951 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
952 description="controls the activation of the diagonalization method", &
953 usage="&DIAGONALIZATION T", &
954 default_l_val=.false., &
955 lone_keyword_l_val=.true.)
956 CALL section_add_keyword(section, keyword)
957 CALL keyword_release(keyword)
958
959 CALL keyword_create(keyword, __location__, name="ALGORITHM", &
960 description="Algorithm to be used for diagonalization", &
961 usage="ALGORITHM STANDARD", &
962 default_i_val=diag_standard, &
963 enum_c_vals=s2a("STANDARD", "OT", "LANCZOS", "DAVIDSON", "FILTER_MATRIX"), &
964 enum_desc=s2a("Standard diagonalization: LAPACK methods or Jacobi.", &
965 "Iterative diagonalization using OT method", &
966 "Block Krylov-space approach to self-consistent diagonalisation", &
967 "Preconditioned blocked Davidson", &
968 "Filter matrix diagonalization"), &
971 CALL section_add_keyword(section, keyword)
972 CALL keyword_release(keyword)
973
974 CALL keyword_create(keyword, __location__, name="JACOBI_THRESHOLD", &
975 description="Controls the accuracy of the pseudo-diagonalization method using Jacobi rotations", &
976 usage="JACOBI_THRESHOLD 1.0E-6", &
977 default_r_val=1.0e-7_dp, &
978 citations=(/stewart1982/))
979 CALL section_add_keyword(section, keyword)
980 CALL keyword_release(keyword)
981
982 CALL keyword_create(keyword, __location__, name="EPS_JACOBI", &
983 description="Below this threshold value for the SCF convergence the pseudo-diagonalization "// &
984 "method using Jacobi rotations is activated. This method is much faster than a "// &
985 "real diagonalization and it is even speeding up while achieving full convergence. "// &
986 "However, it needs a pre-converged wavefunction obtained by at least one real "// &
987 "diagonalization which is further optimized while keeping the original eigenvalue "// &
988 "spectrum. The MO eigenvalues are NOT updated. The method might be useful to speed "// &
989 "up calculations for large systems e.g. using a semi-empirical method.", &
990 usage="EPS_JACOBI 1.0E-5", &
991 default_r_val=0.0_dp, &
992 citations=(/stewart1982/))
993 CALL section_add_keyword(section, keyword)
994 CALL keyword_release(keyword)
995
996 CALL keyword_create(keyword, __location__, name="EPS_ADAPT", &
997 description="Required accuracy in iterative diagonalization as compared to current SCF convergence", &
998 usage="EPS_ADAPT 0.01", &
999 default_r_val=0._dp)
1000 CALL section_add_keyword(section, keyword)
1001 CALL keyword_release(keyword)
1002
1003 CALL keyword_create(keyword, __location__, name="MAX_ITER", &
1004 description="Maximum number of iterations in iterative diagonalization", &
1005 usage="MAX_ITER 20", &
1006 default_i_val=2)
1007 CALL section_add_keyword(section, keyword)
1008 CALL keyword_release(keyword)
1009
1010 CALL keyword_create(keyword, __location__, name="EPS_ITER", &
1011 description="Required accuracy in iterative diagonalization", &
1012 usage="EPS_ITER 1.e-8", &
1013 default_r_val=1.e-8_dp)
1014 CALL section_add_keyword(section, keyword)
1015 CALL keyword_release(keyword)
1016
1017 NULLIFY (subsection)
1018 CALL create_ot_section(subsection)
1019 CALL section_add_subsection(section, subsection)
1020 CALL section_release(subsection)
1021
1022 NULLIFY (subsection)
1023 CALL create_krylov_section(subsection)
1024 CALL section_add_subsection(section, subsection)
1025 CALL section_release(subsection)
1026
1027 NULLIFY (subsection)
1028 CALL create_diag_subspace_section(subsection)
1029 CALL section_add_subsection(section, subsection)
1030 CALL section_release(subsection)
1031
1032 NULLIFY (subsection)
1033 CALL create_davidson_section(subsection)
1034 CALL section_add_subsection(section, subsection)
1035 CALL section_release(subsection)
1036
1037 NULLIFY (subsection)
1038 CALL create_filtermatrix_section(subsection)
1039 CALL section_add_subsection(section, subsection)
1040 CALL section_release(subsection)
1041
1042 END SUBROUTINE create_diagonalization_section
1043
1044! **************************************************************************************************
1045!> \brief ...
1046!> \param section ...
1047! **************************************************************************************************
1048 SUBROUTINE create_davidson_section(section)
1049 TYPE(section_type), POINTER :: section
1050
1051 TYPE(keyword_type), POINTER :: keyword
1052
1053 cpassert(.NOT. ASSOCIATED(section))
1054 CALL section_create(section, __location__, name="DAVIDSON", &
1055 description=" ", &
1056 n_keywords=2, n_subsections=0, repeats=.false.)
1057
1058 NULLIFY (keyword)
1059
1060 CALL keyword_create( &
1061 keyword, __location__, name="PRECONDITIONER", &
1062 description="Type of preconditioner to be used with all minimization schemes. ", &
1063 usage="PRECONDITIONER FULL_ALL", &
1064 default_i_val=ot_precond_full_all, &
1065 enum_c_vals=s2a("FULL_ALL", "FULL_SINGLE_INVERSE", "NONE"), &
1066 enum_desc=s2a("Most effective state selective preconditioner based on diagonalization ", &
1067 "Based on H-eS cholesky inversion, similar to FULL_SINGLE in preconditioning efficiency "// &
1068 "but cheaper to construct, might be somewhat less robust. Recommended for large systems.", &
1069 "skip preconditioning"), &
1071 citations=(/vandevondele2003/))
1072 CALL section_add_keyword(section, keyword)
1073 CALL keyword_release(keyword)
1074
1075 CALL keyword_create(keyword, __location__, name="PRECOND_SOLVER", &
1076 description="How the preconditioner is applied to the residual.", &
1077 usage="PRECOND_SOLVER DIRECT", &
1078 default_i_val=ot_precond_solver_default, &
1079 enum_c_vals=s2a("DEFAULT", "DIRECT", "INVERSE_CHOLESKY"), &
1080 enum_desc=s2a("the default", "Cholesky decomposition followed by triangular solve "// &
1081 "(works for FULL_KINETIC/SINGLE_INVERSE/S_INVERSE)", &
1082 "Cholesky decomposition followed by explicit inversion "// &
1083 "(works for FULL_KINETIC/SINGLE_INVERSE/S_INVERSE)"), &
1084 enum_i_vals=(/ot_precond_solver_default, &
1087 CALL section_add_keyword(section, keyword)
1088 CALL keyword_release(keyword)
1089
1090 CALL keyword_create( &
1091 keyword, __location__, name="ENERGY_GAP", &
1092 description="Should be an estimate for the energy gap [a.u.] (HOMO-LUMO) and is used in preconditioning, "// &
1093 "especially effective with the FULL_ALL preconditioner, in which case it should be an underestimate "// &
1094 "of the gap (0.001 doing normally fine). For the other preconditioners, making this value larger (0.2)"// &
1095 " will tame the preconditioner in case of poor initial guesses.", &
1096 usage="ENERGY_GAP 0.001", &
1097 default_r_val=0.2_dp)
1098 CALL section_add_keyword(section, keyword)
1099 CALL keyword_release(keyword)
1100
1101 CALL keyword_create(keyword, __location__, name="NEW_PREC_EACH", &
1102 description="Number of SCF iterations after which a new Preconditioner is computed", &
1103 usage="NEW_PREC_EACH 10", default_i_val=20)
1104 CALL section_add_keyword(section, keyword)
1105 CALL keyword_release(keyword)
1106
1107 CALL keyword_create(keyword, __location__, name="FIRST_PREC", &
1108 description="First SCF iteration at which a Preconditioner is employed", &
1109 usage="FIRST_PREC 1", default_i_val=1)
1110 CALL section_add_keyword(section, keyword)
1111 CALL keyword_release(keyword)
1112
1113 CALL keyword_create(keyword, __location__, name="CONV_MOS_PERCENT", &
1114 description="Minimal percent of MOS that have to converge within the Davidson loop"// &
1115 " before the SCF iteration is completed and a new Hamiltonian is computed", &
1116 usage="CONV_MOS_PERCENT 0.8", default_r_val=0.5_dp)
1117 CALL section_add_keyword(section, keyword)
1118 CALL keyword_release(keyword)
1119
1120 CALL keyword_create(keyword, __location__, name="SPARSE_MOS", &
1121 description="Use MOS as sparse matrix and avoid as much as possible multiplications with full matrices", &
1122 usage="SPARSE_MOS", default_l_val=.true., &
1123 lone_keyword_l_val=.true.)
1124 CALL section_add_keyword(section, keyword)
1125 CALL keyword_release(keyword)
1126
1127 END SUBROUTINE create_davidson_section
1128
1129! **************************************************************************************************
1130!> \brief ...
1131!> \param section ...
1132! **************************************************************************************************
1133 SUBROUTINE create_krylov_section(section)
1134 TYPE(section_type), POINTER :: section
1135
1136 TYPE(keyword_type), POINTER :: keyword
1137
1138 cpassert(.NOT. ASSOCIATED(section))
1139 CALL section_create(section, __location__, name="KRYLOV", &
1140 description=" ", &
1141 n_keywords=2, n_subsections=0, repeats=.false.)
1142
1143 NULLIFY (keyword)
1144
1145 CALL keyword_create(keyword, __location__, name="NKRYLOV", &
1146 description="Dimension of the Krylov space used for the Lanczos refinement", &
1147 usage="NKRYLOV 20", &
1148 default_i_val=4)
1149 CALL section_add_keyword(section, keyword)
1150 CALL keyword_release(keyword)
1151
1152 CALL keyword_create(keyword, __location__, name="NBLOCK", &
1153 description="Size of the block of vectors refined simultaneously by the Lanczos procedure", &
1154 usage="NBLOCK 1", &
1155 default_i_val=32)
1156 CALL section_add_keyword(section, keyword)
1157 CALL keyword_release(keyword)
1158
1159 CALL keyword_create(keyword, __location__, name="EPS_KRYLOV", &
1160 description="Convergence criterion for the MOs", &
1161 usage="EPS_KRYLOV 0.00001", &
1162 default_r_val=0.0000001_dp)
1163 CALL section_add_keyword(section, keyword)
1164 CALL keyword_release(keyword)
1165
1166 CALL keyword_create(keyword, __location__, name="EPS_STD_DIAG", &
1167 description="Level of convergence to be reached before starting the Lanczos procedure."// &
1168 " Above this threshold a standard diagonalization method is used."// &
1169 " If negative Lanczos is started at the first iteration", &
1170 usage="EPS_STD_DIAG 0.001", &
1171 default_r_val=-1.0_dp)
1172 CALL section_add_keyword(section, keyword)
1173 CALL keyword_release(keyword)
1174
1175 CALL keyword_create(keyword, __location__, name="CHECK_MOS_CONV", &
1176 description="This requires to check the convergence of MOS also when standard "// &
1177 "diagonalization steps are performed, if the block krylov approach is active.", &
1178 usage="CHECK_MOS_CONV T", &
1179 default_l_val=.false., &
1180 lone_keyword_l_val=.true.)
1181 CALL section_add_keyword(section, keyword)
1182 CALL keyword_release(keyword)
1183
1184 END SUBROUTINE create_krylov_section
1185
1186! **************************************************************************************************
1187!> \brief ...
1188!> \param section ...
1189! **************************************************************************************************
1190 SUBROUTINE create_diag_subspace_section(section)
1191 TYPE(section_type), POINTER :: section
1192
1193 TYPE(keyword_type), POINTER :: keyword
1194 TYPE(section_type), POINTER :: subsection
1195
1196 cpassert(.NOT. ASSOCIATED(section))
1197 CALL section_create(section, __location__, name="DIAG_SUB_SCF", &
1198 description="Activation of self-consistenf subspace refinement by diagonalization "// &
1199 "of H by adjusting the occupation but keeping the MOS unchanged.", &
1200 n_keywords=2, n_subsections=1, repeats=.false.)
1201
1202 NULLIFY (keyword, subsection)
1203
1204 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
1205 description="controls the activation of inner SCF loop to refine occupations in MOS subspace", &
1206 usage="&DIAG_SUB_SCF T", &
1207 default_l_val=.false., &
1208 lone_keyword_l_val=.true.)
1209 CALL section_add_keyword(section, keyword)
1210 CALL keyword_release(keyword)
1211
1212 CALL keyword_create(keyword, __location__, name="MAX_ITER", &
1213 description="Maximum number of iterations for the SCF inner loop", &
1214 usage="MAX_ITER 20", &
1215 default_i_val=2)
1216 CALL section_add_keyword(section, keyword)
1217 CALL keyword_release(keyword)
1218
1219 CALL keyword_create(keyword, __location__, name="EPS_ENE", &
1220 description="Required energy accuracy for convergence of subspace diagonalization", &
1221 usage="EPS_ENE 1.e-8", &
1222 default_r_val=1.e-4_dp)
1223 CALL section_add_keyword(section, keyword)
1224 CALL keyword_release(keyword)
1225
1226 CALL keyword_create(keyword, __location__, name="EPS_ADAPT_SCF", &
1227 description="Required density matrix accuracy as compared to current SCF convergence", &
1228 usage="EPS_ADAPT_SCF 1.e-1", &
1229 default_r_val=1._dp)
1230 CALL section_add_keyword(section, keyword)
1231 CALL keyword_release(keyword)
1232
1233 CALL keyword_create( &
1234 keyword, __location__, name="EPS_SKIP_SUB_DIAG", &
1235 description="Level of convergence to be reached before starting the internal loop of subspace rotations."// &
1236 " Above this threshold only the outer diagonalization method is used."// &
1237 " If negative the subspace rotation is started at the first iteration", &
1238 usage="EPS_SKIP_SUB_DIAG 0.001", &
1239 default_r_val=-1.0_dp)
1240 CALL section_add_keyword(section, keyword)
1241 CALL keyword_release(keyword)
1242
1243 CALL create_mixing_section(subsection)
1244 CALL section_add_subsection(section, subsection)
1245 CALL section_release(subsection)
1246 END SUBROUTINE create_diag_subspace_section
1247
1248! **************************************************************************************************
1249!> \brief Create CP2K input section for the smearing of occupation numbers
1250!> \param section ...
1251!> \date 27.08.2008
1252!> \author Matthias Krack (MK)
1253!> \version 1.0
1254! **************************************************************************************************
1255 SUBROUTINE create_smear_section(section)
1256
1257 TYPE(section_type), POINTER :: section
1258
1259 TYPE(keyword_type), POINTER :: keyword
1260
1261 cpassert(.NOT. ASSOCIATED(section))
1262
1263 CALL section_create(section, __location__, &
1264 name="SMEAR", &
1265 description="Define the smearing of the MO occupation numbers", &
1266 n_keywords=6, &
1267 n_subsections=0, &
1268 repeats=.false.)
1269
1270 NULLIFY (keyword)
1271
1272 CALL keyword_create(keyword, __location__, &
1273 name="_SECTION_PARAMETERS_", &
1274 description="Controls the activation of smearing", &
1275 usage="&SMEAR ON", &
1276 default_l_val=.false., &
1277 lone_keyword_l_val=.true.)
1278 CALL section_add_keyword(section, keyword)
1279 CALL keyword_release(keyword)
1280
1281 CALL keyword_create(keyword, __location__, &
1282 name="METHOD", &
1283 description="Smearing method to be applied", &
1284 usage="METHOD Fermi_Dirac", &
1285 default_i_val=smear_energy_window, &
1286 enum_c_vals=s2a("FERMI_DIRAC", "ENERGY_WINDOW", "LIST"), &
1288 enum_desc=s2a("Fermi-Dirac distribution defined by the keyword ELECTRONIC_TEMPERATURE", &
1289 "Energy window defined by the keyword WINDOW_SIZE", &
1290 "Use a fixed list of occupations"))
1291 CALL section_add_keyword(section, keyword)
1292 CALL keyword_release(keyword)
1293
1294 CALL keyword_create(keyword, __location__, &
1295 name="LIST", &
1296 description="A list of fractional occupations to use. Must match the number of states "// &
1297 "and sum up to the correct number of electrons", &
1298 repeats=.false., &
1299 n_var=-1, &
1300 type_of_var=real_t, &
1301 usage="LIST 2.0 0.6666 0.6666 0.66666 0.0 0.0")
1302 CALL section_add_keyword(section, keyword)
1303 CALL keyword_release(keyword)
1304
1305 CALL keyword_create(keyword, __location__, &
1306 name="ELECTRONIC_TEMPERATURE", &
1307 variants=s2a("ELEC_TEMP", "TELEC"), &
1308 description="Electronic temperature in the case of Fermi-Dirac smearing", &
1309 repeats=.false., &
1310 n_var=1, &
1311 type_of_var=real_t, &
1312 default_r_val=cp_unit_to_cp2k(value=300.0_dp, unit_str="K"), &
1313 unit_str="K", &
1314 usage="ELECTRONIC_TEMPERATURE [K] 300")
1315 CALL section_add_keyword(section, keyword)
1316 CALL keyword_release(keyword)
1317
1318 CALL keyword_create(keyword, __location__, &
1319 name="EPS_FERMI_DIRAC", &
1320 description="Accuracy checks on occupation numbers use this as a tolerance", &
1321 repeats=.false., &
1322 n_var=1, &
1323 type_of_var=real_t, &
1324 default_r_val=1.0e-10_dp, &
1325 usage="EPS_FERMI_DIRAC 1.0E-6")
1326 CALL section_add_keyword(section, keyword)
1327 CALL keyword_release(keyword)
1328
1329 CALL keyword_create(keyword, __location__, &
1330 name="WINDOW_SIZE", &
1331 description="Size of the energy window centred at the Fermi level", &
1332 repeats=.false., &
1333 n_var=1, &
1334 type_of_var=real_t, &
1335 default_r_val=0.0_dp, &
1336 unit_str="au_e", &
1337 usage="WINDOW_SIZE [eV] 0.3")
1338 CALL section_add_keyword(section, keyword)
1339 CALL keyword_release(keyword)
1340
1341 CALL keyword_create(keyword, __location__, name="FIXED_MAGNETIC_MOMENT", &
1342 description="Imposed difference between the numbers of electrons of spin up "// &
1343 "and spin down: m = n(up) - n(down). A negative value (default) allows "// &
1344 "for a change of the magnetic moment. -1 specifically keeps an integer "// &
1345 "number of spin up and spin down electrons.", &
1346 repeats=.false., &
1347 n_var=1, &
1348 type_of_var=real_t, &
1349 default_r_val=-100.0_dp, &
1350 usage="FIXED_MAGNETIC_MOMENT 1.5")
1351 CALL section_add_keyword(section, keyword)
1352 CALL keyword_release(keyword)
1353
1354 END SUBROUTINE create_smear_section
1355
1356! **************************************************************************************************
1357!> \brief Creates the input section for defining CDFT constraints.
1358!> \param section the section to create
1359! **************************************************************************************************
1361 TYPE(section_type), POINTER :: section
1362
1363 TYPE(keyword_type), POINTER :: keyword
1364 TYPE(section_type), POINTER :: group_section, print_key, subsection
1365
1366 NULLIFY (keyword, subsection, group_section, print_key)
1367
1368 cpassert(.NOT. ASSOCIATED(section))
1369 CALL section_create(section, __location__, name="CDFT", &
1370 description="Parameters needed to set up a constrained DFT calculation."// &
1371 " Each repetition of the ATOM_GROUP section defines a new constraint."// &
1372 " The constraint(s) is (are) converged in a separate external SCF loop with settings"// &
1373 " read from the OUTER_SCF section. Supported constraints: Becke and Gaussian"// &
1374 " Hirshfeld (partial).", n_keywords=8, n_subsections=2, &
1375 repeats=.false., citations=(/holmberg2017, holmberg2018/))
1376
1377 NULLIFY (subsection, keyword)
1378 CALL create_outer_scf_section(subsection)
1379 CALL section_add_subsection(section, subsection)
1380 CALL section_release(subsection)
1381
1382 CALL create_becke_constraint_section(subsection)
1383 CALL section_add_subsection(section, subsection)
1384 CALL section_release(subsection)
1385
1386 CALL create_hirshfeld_constraint_section(subsection)
1387 CALL section_add_subsection(section, subsection)
1388 CALL section_release(subsection)
1389
1390 CALL keyword_create(keyword, __location__, name="TYPE_OF_CONSTRAINT", &
1391 description="Specifies the type of constraint used.", &
1392 usage="TYPE_OF_CONSTRAINT (NONE|HIRSHFELD|BECKE)", &
1393 enum_c_vals=s2a("NONE", "HIRSHFELD", "BECKE"), &
1396 enum_desc=s2a("No constraint (disables section).", &
1397 "Gaussian Hirshfeld constraint. Partial implementation: no forces. "// &
1398 "Requires corresponding section. Not as extensively tested.", &
1399 "Becke constraint. Requires corresponding section."), &
1400 citations=(/becke1988b/), &
1401 default_i_val=outer_scf_none)
1402 CALL section_add_keyword(section, keyword)
1403 CALL keyword_release(keyword)
1404
1405 CALL keyword_create(keyword, __location__, name="STRENGTH", &
1406 description="Constraint force constants (Lagrange multipliers). "// &
1407 "Give one value per constraint group.", &
1408 type_of_var=real_t, n_var=-1, &
1409 default_r_val=0.0_dp)
1410 CALL section_add_keyword(section, keyword)
1411 CALL keyword_release(keyword)
1412
1413 CALL keyword_create(keyword, __location__, name="TARGET", &
1414 description="Constraint target values. Give one value per constraint group. "// &
1415 "The target value is the desired number of valence electrons, spin moment, or the number of "// &
1416 "alpha or beta electrons on the atoms that define the constraint, suitably multiplied by "// &
1417 "atomic coefficients in case a relative constraint between two sets of atoms is employed. "// &
1418 "Note that core charges are not subtracted from the target value.", &
1419 usage="TARGET {real}", repeats=.false., &
1420 type_of_var=real_t, n_var=-1, &
1421 default_r_val=0.0_dp)
1422 CALL section_add_keyword(section, keyword)
1423 CALL keyword_release(keyword)
1424
1425 CALL keyword_create(keyword, __location__, name="ATOMIC_CHARGES", &
1426 description="Calculate atomic CDFT charges with selected weight function"// &
1427 " (Z = Z_core - Z_CDFT). With fragment based constraints, charges are"// &
1428 " relative to the fragment reference state i.e. Z = Z_CDFT -"// &
1429 " Z_frag_reference. Note: if the number of atoms is greater than the"// &
1430 " default pw_pool max cache, calculation of atomic CDFT charges"// &
1431 " will prompt harmless warnings during deallocation of atomic grids.", &
1432 usage="ATOMIC_CHARGES", &
1433 default_l_val=.false., lone_keyword_l_val=.true.)
1434 CALL section_add_keyword(section, keyword)
1435 CALL keyword_release(keyword)
1436
1437 CALL keyword_create(keyword, __location__, name="FRAGMENT_A_FILE_NAME", variants=(/"FRAGMENT_A_FILE"/), &
1438 description="Name of the reference total electron density cube file for fragment A."// &
1439 " May include a path. The reference electron density needs to be outputted"// &
1440 " on the same grid as the full system (same cutoff and cell, output stride 1).", &
1441 usage="FRAGMENT_A_FILE_NAME <FILENAME>", &
1442 default_lc_val="fragment_a.cube")
1443 CALL section_add_keyword(section, keyword)
1444 CALL keyword_release(keyword)
1445
1446 CALL keyword_create(keyword, __location__, name="FRAGMENT_B_FILE_NAME", variants=(/"FRAGMENT_B_FILE"/), &
1447 description="Name of the reference total electron density cube file for fragment B."// &
1448 " May include a path. The reference electron density needs to be outputted"// &
1449 " on the same grid as the full system (same cutoff and cell, output stride 1).", &
1450 usage="FRAGMENT_B_FILE_NAME <FILENAME>", &
1451 default_lc_val="fragment_b.cube")
1452 CALL section_add_keyword(section, keyword)
1453 CALL keyword_release(keyword)
1454
1455 CALL keyword_create(keyword, __location__, name="FRAGMENT_A_SPIN_FILE", &
1456 variants=(/"FRAGMENT_A_SPIN_FILE_NAME"/), &
1457 description="Name of the reference spin density cube file for fragment A."// &
1458 " May include a path. The reference spin density needs to be outputted"// &
1459 " on the same grid as the full system (same cutoff and cell, output stride 1).", &
1460 usage="FRAGMENT_A_SPIN_FILE <FILENAME>", &
1461 default_lc_val="fragment_a_spin.cube")
1462 CALL section_add_keyword(section, keyword)
1463 CALL keyword_release(keyword)
1464
1465 CALL keyword_create(keyword, __location__, name="FRAGMENT_B_SPIN_FILE", &
1466 variants=(/"FRAGMENT_B_SPIN_FILE_NAME"/), &
1467 description="Name of the reference spin density cube file for fragment B."// &
1468 " May include a path. The reference spin density needs to be outputted"// &
1469 " on the same grid as the full system (same cutoff and cell, output stride 1).", &
1470 usage="FRAGMENT_B_SPIN_FILE <FILENAME>", &
1471 default_lc_val="fragment_b_spin.cube")
1472 CALL section_add_keyword(section, keyword)
1473 CALL keyword_release(keyword)
1474
1475 CALL keyword_create(keyword, __location__, name="FLIP_FRAGMENT_A", &
1476 description="Logical which determines if the reference spin difference density "// &
1477 "(rho_alpha-rho_beta) for fragment A should be flipped. With default (off) "// &
1478 "value, the fragment is constrained to have more alpha than beta electrons "// &
1479 "if the isolated fragment has unpaired electrons. Useful in conjunction with "// &
1480 "FLIP_FRAGMENT_B.", &
1481 usage="FLIP_FRAGMENT_A", &
1482 default_l_val=.false., lone_keyword_l_val=.true.)
1483 CALL section_add_keyword(section, keyword)
1484 CALL keyword_release(keyword)
1485
1486 CALL keyword_create(keyword, __location__, name="FLIP_FRAGMENT_B", &
1487 description="Logical which determines if the reference spin difference density "// &
1488 "(rho_alpha-rho_beta) for fragment B should be flipped. With default (off) "// &
1489 "value, the fragment is constrained to have more alpha than beta electrons "// &
1490 "if the isolated fragment has unpaired electrons. Useful in conjunction with "// &
1491 "FLIP_FRAGMENT_A.", &
1492 usage="FLIP_FRAGMENT_B", &
1493 default_l_val=.false., lone_keyword_l_val=.true.)
1494 CALL section_add_keyword(section, keyword)
1495 CALL keyword_release(keyword)
1496
1497 CALL cp_print_key_section_create(print_key, __location__, "PROGRAM_RUN_INFO", &
1498 description="Controls the printing of basic info about the method.", &
1499 print_level=low_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
1500
1501 CALL section_create(subsection, __location__, name="WEIGHT_FUNCTION", &
1502 description="Controls the printing of cube files with "// &
1503 "the CDFT weight function(s). Intended for single-point testing. "// &
1504 "In multistep simulations, generated cube files are overwritten each step.", &
1505 n_keywords=1, n_subsections=0, repeats=.false.)
1506
1507 CALL keyword_create(keyword, __location__, name="STRIDE", &
1508 description="The stride (X,Y,Z) used to write the cube file "// &
1509 "(larger values result in smaller cube files). You can provide 3 numbers (for X,Y,Z) or"// &
1510 " 1 number valid for all components.", &
1511 usage="STRIDE 2 2 2", n_var=-1, default_i_vals=(/2, 2, 2/), type_of_var=integer_t)
1512 CALL section_add_keyword(subsection, keyword)
1513 CALL keyword_release(keyword)
1514
1515 CALL section_add_subsection(print_key, subsection)
1516 CALL section_release(subsection)
1517
1518 CALL section_add_subsection(section, print_key)
1519 CALL section_release(print_key)
1520
1521 CALL section_create(group_section, __location__, name="ATOM_GROUP", &
1522 description="Define a group of atoms for use in a CDFT constraint. Each repetition of "// &
1523 "this section creates a new constraint.", &
1524 n_keywords=4, n_subsections=0, repeats=.true.)
1525
1526 CALL keyword_create(keyword, __location__, name="ATOMS", &
1527 description="Specifies the list of atoms that are included in the constraint group.", &
1528 usage="ATOMS {integer} {integer} .. {integer}", &
1529 n_var=-1, type_of_var=integer_t)
1530 CALL section_add_keyword(group_section, keyword)
1531 CALL keyword_release(keyword)
1532
1533 CALL keyword_create(keyword, __location__, name="COEFF", &
1534 description="Defines coefficients for the atoms in the list of atoms. Accepts values +/-1.0.", &
1535 usage="COEFF 1.0 -1.0", repeats=.true., &
1536 type_of_var=real_t, n_var=-1)
1537 CALL section_add_keyword(group_section, keyword)
1538 CALL keyword_release(keyword)
1539
1540 CALL keyword_create(keyword, __location__, name="CONSTRAINT_TYPE ", &
1541 description="Determines what type of constraint to apply. ", &
1542 usage="CONSTRAINT_TYPE (CHARGE|MAGNETIZATION|ALPHA|BETA)", &
1543 enum_c_vals=s2a("CHARGE", "MAGNETIZATION", "ALPHA", "BETA"), &
1546 enum_desc=s2a("Total charge density constraint (rho_alpha + rho_beta).", &
1547 "Magnetization density constraint (rho_alpha - rho_beta).", &
1548 "Alpha spin density constraint.", &
1549 "Beta spin density constraint."), &
1550 default_i_val=cdft_charge_constraint)
1551 CALL section_add_keyword(group_section, keyword)
1552 CALL keyword_release(keyword)
1553
1554 CALL keyword_create(keyword, __location__, name="FRAGMENT_CONSTRAINT", &
1555 description="Use a fragment based constraint. "// &
1556 "Takes as input the electron densities of two isolated fragments in the "// &
1557 "same geometry that they have in the full system. "// &
1558 "The isolated fragment densities are read from cube files defined in FRAGMENT_{A,B}_FILE. "// &
1559 "For magnetization density constraints, additional files containing the spin difference "// &
1560 "densities must be defined with the keywords FRAGMENT_{A,B}_SPIN_FILE. "// &
1561 "With this keyword active, the target value of the constraint is calculated from the "// &
1562 "the superposition of the isolated fragment densities. Supports only static calculations.", &
1563 usage="FRAGMENT_CONSTRAINT", &
1564 default_l_val=.false., lone_keyword_l_val=.true.)
1565 CALL section_add_keyword(group_section, keyword)
1566 CALL keyword_release(keyword)
1567
1568 CALL section_add_subsection(section, group_section)
1569 CALL section_release(group_section)
1570
1571 CALL section_create(group_section, __location__, name="DUMMY_ATOMS", &
1572 description="Define an extra group of atoms for which only atomic CDFT charges "// &
1573 "should be computed. The section cannot contain any constraint "// &
1574 "atoms that were included in section ATOM_GROUP.", &
1575 n_keywords=1, n_subsections=0, repeats=.true.)
1576
1577 CALL keyword_create(keyword, __location__, name="ATOMS", &
1578 description="Specifies the list of atoms that are included in the DUMMY_ATOMS group.", &
1579 usage="ATOMS {integer} {integer} .. {integer}", &
1580 n_var=-1, type_of_var=integer_t)
1581 CALL section_add_keyword(group_section, keyword)
1582 CALL keyword_release(keyword)
1583
1584 CALL section_add_subsection(section, group_section)
1585 CALL section_release(group_section)
1586
1587 CALL keyword_create(keyword, __location__, name="REUSE_PRECOND", &
1588 description="Reuse a previously built OT preconditioner between subsequent CDFT SCF iterations "// &
1589 "if the inner OT SCF loop converged in PRECOND_FREQ steps or less. Intended mainly for MD "// &
1590 "simulations with the FULL_ALL preconditioner to speed up the final iterations of the CDFT SCF loop.", &
1591 usage="REUSE_PRECOND yes", repeats=.false., n_var=1, &
1592 default_l_val=.false., lone_keyword_l_val=.true.)
1593 CALL section_add_keyword(section, keyword)
1594 CALL keyword_release(keyword)
1595
1596 CALL keyword_create(keyword, __location__, name="PRECOND_FREQ", &
1597 description="See REUSE_PRECOND.", &
1598 usage="PRECOND_FREQ {int}", default_i_val=0)
1599 CALL section_add_keyword(section, keyword)
1600 CALL keyword_release(keyword)
1601
1602 CALL keyword_create(keyword, __location__, name="MAX_REUSE", &
1603 description="Determines how many times a previously built preconditioner can be reused.", &
1604 usage="MAX_REUSE {int}", default_i_val=0)
1605 CALL section_add_keyword(section, keyword)
1606 CALL keyword_release(keyword)
1607
1608 CALL keyword_create(keyword, __location__, name="PURGE_HISTORY", &
1609 description="Purge wavefunction and constraint history to improve SCF convergence during MD."// &
1610 " Counts how often the convergence of the first CDFT SCF iteration takes 2 or more outer SCF"// &
1611 " iterations and purges the history if the counter exceeds PURGE_FREQ, and PURGE_OFFSET"// &
1612 " MD steps have passed since the last purge."// &
1613 " The counter is zeroed after each purge.", &
1614 usage="PURGE_HISTORY yes", repeats=.false., n_var=1, &
1615 default_l_val=.false., lone_keyword_l_val=.true.)
1616 CALL section_add_keyword(section, keyword)
1617 CALL keyword_release(keyword)
1618
1619 CALL keyword_create(keyword, __location__, name="PURGE_FREQ", &
1620 description="See PURGE_HISTORY.", &
1621 usage="PURGE_FREQ {int} ", default_i_val=1)
1622 CALL section_add_keyword(section, keyword)
1623 CALL keyword_release(keyword)
1624
1625 CALL keyword_create(keyword, __location__, name="PURGE_OFFSET", &
1626 description="See PURGE_HISTORY.", &
1627 usage="PURGE_OFFSET {int} ", default_i_val=1)
1628 CALL section_add_keyword(section, keyword)
1629 CALL keyword_release(keyword)
1630
1631 CALL keyword_create(keyword, __location__, name="COUNTER", &
1632 description="A counter to track the total number of energy evaluations. Needed by"// &
1633 " some optimizers to print information. Useful mainly for restarts.", &
1634 usage="COUNTER {int} ", default_i_val=0)
1635 CALL section_add_keyword(section, keyword)
1636 CALL keyword_release(keyword)
1637
1638 CALL keyword_create(keyword, __location__, name="IN_MEMORY", &
1639 description="Precompute gradients due to constraint during"// &
1640 " initial formation of constraint and store them in memory. Does"// &
1641 " nothing if forces are not calculated.", &
1642 usage="IN_MEMORY", &
1643 default_l_val=.false., lone_keyword_l_val=.true.)
1644 CALL section_add_keyword(section, keyword)
1645 CALL keyword_release(keyword)
1646
1647 END SUBROUTINE create_cdft_control_section
1648
1649! **************************************************************************************************
1650!> \brief Creates the input section for defining Gaussian Hirshfeld CDFT constraints.
1651!> \param section the section to create
1652! **************************************************************************************************
1653 SUBROUTINE create_hirshfeld_constraint_section(section)
1654 TYPE(section_type), POINTER :: section
1655
1656 TYPE(keyword_type), POINTER :: keyword
1657
1658 NULLIFY (keyword)
1659
1660 cpassert(.NOT. ASSOCIATED(section))
1661 CALL section_create(section, __location__, name="HIRSHFELD_CONSTRAINT", &
1662 description="Parameters for CDFT with a Gaussian Hirshfeld constraint.", &
1663 n_keywords=11, n_subsections=0, repeats=.false.)
1664
1665 CALL keyword_create(keyword, __location__, name="SHAPE_FUNCTION", &
1666 description="Type of shape function used for Hirshfeld partitioning.", &
1667 usage="SHAPE_FUNCTION {Gaussian,Density}", repeats=.false., n_var=1, &
1668 default_i_val=shape_function_gaussian, &
1669 enum_c_vals=s2a("GAUSSIAN", "DENSITY"), &
1670 enum_desc=s2a("One Gaussian per atom with radius determined by the keyword GAUSSIAN_SHAPE.", &
1671 "Atomic density expanded in terms of multiple Gaussians."), &
1673 CALL section_add_keyword(section, keyword)
1674 CALL keyword_release(keyword)
1675
1676 CALL keyword_create(keyword, __location__, name="GAUSSIAN_SHAPE", &
1677 description="Specifies the type of Gaussian used for SHAPE_FUNCTION GAUSSIAN.", &
1678 usage="GAUSSIAN_SHAPE (SINGLE|VDW|COVALENT|USER)", &
1679 enum_c_vals=s2a("DEFAULT", "SINGLE", "VDW", "COVALENT", "USER"), &
1681 enum_desc=s2a("Use covalent radii (in angstrom) to construct Gaussians, but fixed"// &
1682 " 1.0_dp radius for elements with a radius larger than this value.", &
1683 "Single Gaussian for all atom types with radius given by GAUSSIAN_RADIUS.", &
1684 "Use van der Waals radii to construct Gaussians.", &
1685 "Use covalent radii to construct Gaussians.", &
1686 "Use user defined radii (keyword ATOMIC_RADII) to construct Gaussians."), &
1687 default_i_val=radius_default)
1688 CALL section_add_keyword(section, keyword)
1689 CALL keyword_release(keyword)
1690
1691 CALL keyword_create(keyword, __location__, name="GAUSSIAN_RADIUS", &
1692 description="Radius parameter controlling the creation of Gaussians.", &
1693 usage="GAUSSIAN_RADIUS <REAL>", &
1694 unit_str="angstrom", &
1695 default_r_val=cp_unit_to_cp2k(3.0_dp, "angstrom"), &
1696 type_of_var=real_t, n_var=1)
1697 CALL section_add_keyword(section, keyword)
1698 CALL keyword_release(keyword)
1699
1700 CALL keyword_create(keyword, __location__, name="ATOMIC_RADII", &
1701 description="Defines custom radii to setup the spherical Gaussians. "// &
1702 "Give one value per element in the same order as they "// &
1703 "appear in the input coordinates.", &
1704 usage="ATOMIC_RADII {real} {real} {real}", repeats=.false., &
1705 unit_str="angstrom", &
1706 type_of_var=real_t, n_var=-1)
1707 CALL section_add_keyword(section, keyword)
1708 CALL keyword_release(keyword)
1709
1710 CALL keyword_create(keyword, __location__, name="USE_BOHR", &
1711 description="Convert the Gaussian radius from angstrom to bohr. This results in a larger "// &
1712 "Gaussian than without unit conversion.", &
1713 usage="USE_BOHR .TRUE.", &
1714 default_l_val=.false., lone_keyword_l_val=.true.)
1715 CALL section_add_keyword(section, keyword)
1716 CALL keyword_release(keyword)
1717
1718 CALL keyword_create(keyword, __location__, name="PRINT_DENSITY", &
1719 description="Logical to control printing of Hirshfeld densities to .cube file.", &
1720 usage="PRINT_DENSITY TRUE", &
1721 default_l_val=.false., lone_keyword_l_val=.true.)
1722 CALL section_add_keyword(section, keyword)
1723 CALL keyword_release(keyword)
1724
1725 CALL keyword_create(keyword, __location__, name="ATOMS_MEMORY", &
1726 description="Number of atomic gradients to store in memory.", &
1727 usage="ATOMS_MEMORY", &
1728 n_var=1, type_of_var=integer_t, &
1729 default_i_val=80)
1730 CALL section_add_keyword(section, keyword)
1731 CALL keyword_release(keyword)
1732
1733 CALL keyword_create(keyword, __location__, name="USE_ATOMIC_CUTOFF", &
1734 description="Logical to control use of ATOMIC_CUTOFF.", &
1735 usage="USE_ATOMIC_CUTOFF TRUE", &
1736 default_l_val=.true., lone_keyword_l_val=.true.)
1737 CALL section_add_keyword(section, keyword)
1738 CALL keyword_release(keyword)
1739
1740 CALL keyword_create(keyword, __location__, name="EPS_CUTOFF", &
1741 description="Numerical cutoff for calculation of weight function.", &
1742 usage="EPS_CUTOFF {real} ", default_r_val=1.0e-12_dp)
1743 CALL section_add_keyword(section, keyword)
1744 CALL keyword_release(keyword)
1745
1746 CALL keyword_create(keyword, __location__, name="ATOMIC_CUTOFF", &
1747 description="Numerical cutoff for calculation of Hirshfeld densities.", &
1748 usage="ATOMIC_CUTOFF {real} ", default_r_val=1.0e-12_dp)
1749 CALL section_add_keyword(section, keyword)
1750 CALL keyword_release(keyword)
1751
1752 END SUBROUTINE create_hirshfeld_constraint_section
1753
1754! **************************************************************************************************
1755!> \brief Create input section to define CDFT constraint settings specific to Becke weight function.
1756!> \param section the section to create
1757! **************************************************************************************************
1758 SUBROUTINE create_becke_constraint_section(section)
1759 TYPE(section_type), POINTER :: section
1760
1761 TYPE(keyword_type), POINTER :: keyword
1762
1763 NULLIFY (keyword)
1764 cpassert(.NOT. ASSOCIATED(section))
1765 CALL section_create(section, __location__, name="BECKE_CONSTRAINT", &
1766 description="Define settings influencing the construction of the Becke weight function.", &
1767 n_keywords=13, repeats=.false., citations=(/becke1988b/))
1768
1769 CALL keyword_create(keyword, __location__, name="ADJUST_SIZE", &
1770 description="Adjust Becke cell boundaries with atomic"// &
1771 " radii to generate a heteronuclear cutoff profile. These"// &
1772 " radii are defined with the keyword ATOMIC_RADII.", &
1773 usage="ADJUST_SIZE", &
1774 default_l_val=.false., lone_keyword_l_val=.true.)
1775 CALL section_add_keyword(section, keyword)
1776 CALL keyword_release(keyword)
1777
1778 CALL keyword_create(keyword, __location__, name="ATOMIC_RADII", &
1779 description="Defines atomic radii to generate a heteronuclear cutoff profile."// &
1780 " Give one value per element in the same order as they"// &
1781 " appear in the input coordinates.", &
1782 usage="ATOMIC_RADII {real} {real} {real}", repeats=.false., &
1783 unit_str="angstrom", &
1784 type_of_var=real_t, n_var=-1)
1785 CALL section_add_keyword(section, keyword)
1786 CALL keyword_release(keyword)
1787
1788 CALL keyword_create(keyword, __location__, name="SHOULD_SKIP", &
1789 description="If grid point is farther than GLOBAL_CUTOFF from all constraint atoms, "// &
1790 "move directly to next grid point, thus saving computational resources.", &
1791 usage="SHOULD_SKIP", &
1792 default_l_val=.false., lone_keyword_l_val=.true.)
1793 CALL section_add_keyword(section, keyword)
1794 CALL keyword_release(keyword)
1795
1796 CALL keyword_create(keyword, __location__, name="CAVITY_CONFINE", &
1797 description="Activates Gaussian cavity confinement. The constraint is evaluated only inside "// &
1798 "the cavity. The cavity is formed by summing spherical Gaussians centered on the constraint atoms.", &
1799 usage="CAVITY_CONFINE", &
1800 default_l_val=.false., lone_keyword_l_val=.true.)
1801 CALL section_add_keyword(section, keyword)
1802 CALL keyword_release(keyword)
1803
1804 CALL keyword_create(keyword, __location__, name="CAVITY_SHAPE", &
1805 description="Specifies the type of Gaussian cavity used.", &
1806 usage="CAVITY_SHAPE (SINGLE|VDW|COVALENT|USER)", &
1807 enum_c_vals=s2a("DEFAULT", "SINGLE", "VDW", "COVALENT", "USER"), &
1809 enum_desc=s2a("Use covalent radii (in angstrom) to construct Gaussians, but fixed"// &
1810 " 1.0_dp radius for elements with a radius larger than this value.", &
1811 "Single Gaussian for all atom types with radius given by CAVITY_RADIUS.", &
1812 "Use van der Waals radii to construct Gaussians.", &
1813 "Use covalent radii to construct Gaussians.", &
1814 "Use user defined radii (keyword ATOMIC_RADII) to construct Gaussians."), &
1815 default_i_val=radius_default)
1816 CALL section_add_keyword(section, keyword)
1817 CALL keyword_release(keyword)
1818
1819 CALL keyword_create(keyword, __location__, name="CAVITY_USE_BOHR", &
1820 description="Convert the cavity radius from angstrom to bohr. This results in a larger"// &
1821 " confinement cavity than without unit conversion.", &
1822 usage="CAVITY_USE_BOHR TRUE", &
1823 default_l_val=.false., lone_keyword_l_val=.true.)
1824 CALL section_add_keyword(section, keyword)
1825 CALL keyword_release(keyword)
1826
1827 CALL keyword_create(keyword, __location__, name="CAVITY_PRINT", &
1828 description="Print cavity in Gaussian cube file format. Currently, printing options"// &
1829 " are hardcoded.", &
1830 usage="CAVITY_PRINT", &
1831 default_l_val=.false., lone_keyword_l_val=.true.)
1832 CALL section_add_keyword(section, keyword)
1833 CALL keyword_release(keyword)
1834
1835 CALL keyword_create(keyword, __location__, name="CAVITY_RADIUS", &
1836 description="Radius parameter controlling the creation of Gaussian cavity confinement.", &
1837 usage="CAVITY_RADIUS <REAL>", &
1838 unit_str="angstrom", &
1839 default_r_val=cp_unit_to_cp2k(3.0_dp, "angstrom"), &
1840 type_of_var=real_t, n_var=1)
1841 CALL section_add_keyword(section, keyword)
1842 CALL keyword_release(keyword)
1843
1844 CALL keyword_create(keyword, __location__, name="EPS_CAVITY", &
1845 description="Density threshold for cavity creation. Grid points where the Gaussian"// &
1846 " density falls below the threshold are ignored.", &
1847 usage="EPS_CAVITY {real} ", default_r_val=1.0e-6_dp)
1848 CALL section_add_keyword(section, keyword)
1849 CALL keyword_release(keyword)
1850
1851 CALL keyword_create(keyword, __location__, name="CUTOFF_TYPE", &
1852 description="Specifies the type of cutoff used when building the Becke weight function.", &
1853 usage="CUTOFF_TYPE (GLOBAL|ELEMENT)", &
1854 enum_c_vals=s2a("GLOBAL", "ELEMENT"), &
1855 enum_i_vals=(/becke_cutoff_global, becke_cutoff_element/), &
1856 enum_desc=s2a("Use a single value for all elements. Read from GLOBAL_CUTOFF.", &
1857 "Use a different value for all elements. Values read from ELEMENT_CUTOFF."), &
1858 default_i_val=becke_cutoff_global)
1859 CALL section_add_keyword(section, keyword)
1860 CALL keyword_release(keyword)
1861
1862 CALL keyword_create(keyword, __location__, name="GLOBAL_CUTOFF", &
1863 description="Parameter used to select which atoms contribute to the"// &
1864 " weight function at each real space grid point.", &
1865 usage="GLOBAL_CUTOFF <REAL>", &
1866 unit_str="angstrom", &
1867 default_r_val=cp_unit_to_cp2k(3.1750632515_dp, "angstrom"), &
1868 type_of_var=real_t, n_var=1)
1869 CALL section_add_keyword(section, keyword)
1870 CALL keyword_release(keyword)
1871
1872 CALL keyword_create(keyword, __location__, name="ELEMENT_CUTOFF", &
1873 description="Defines element specific cutoffs to decide which atoms contribute to the"// &
1874 " weight function at each real space grid point. Give one value per element in the same"// &
1875 " order as they appear in the coordinates.", &
1876 usage="ELEMENT_CUTOFF {real} {real} {real}", repeats=.false., &
1877 unit_str="angstrom", &
1878 type_of_var=real_t, n_var=-1)
1879 CALL section_add_keyword(section, keyword)
1880 CALL keyword_release(keyword)
1881
1882 CALL keyword_create(keyword, __location__, name="IN_MEMORY", &
1883 description="Precompute gradients due to Becke constraint during"// &
1884 " initial formation of constraint and store them in memory. Useful"// &
1885 " in combination with confinement, memory intensive otherwise. Does"// &
1886 " nothing if forces are not calculated.", &
1887 usage="IN_MEMORY", &
1888 default_l_val=.false., lone_keyword_l_val=.true.)
1889 CALL section_add_keyword(section, keyword)
1890 CALL keyword_release(keyword)
1891
1892 END SUBROUTINE create_becke_constraint_section
1893
1894! **************************************************************************************************
1895!> \brief creates the input section for parameters related to CDFT specific optimizers
1896!> \param section the section to be created
1897!> \par History
1898!> 03.2018 separated from create_outer_scf_section [Nico Holmberg]
1899!> \author Nico Holmberg
1900! **************************************************************************************************
1901 SUBROUTINE create_cdft_opt_section(section)
1902 TYPE(section_type), POINTER :: section
1903
1904 TYPE(keyword_type), POINTER :: keyword
1905
1906 cpassert(.NOT. ASSOCIATED(section))
1907 CALL section_create(section, __location__, name="CDFT_OPT", &
1908 description="Parameters controlling optimization methods that are compatible "// &
1909 "only with CDFT based constraints (i.e. CDFT SCF is active). Specifically, "// &
1910 "the control parameters for the Broyden and Newton optimizers are defined in this "// &
1911 "section.", &
1912 n_keywords=10, n_subsections=0, repeats=.false.)
1913
1914 NULLIFY (keyword)
1915
1916 CALL keyword_create(keyword, __location__, name="BROYDEN_TYPE", &
1917 description="Specifies the Broyden optimizer variant to use.", &
1918 usage="BROYDEN_TYPE BT1", &
1919 default_i_val=broyden_type_1, &
1920 enum_c_vals=s2a("BT1", "BT1_EXPLICIT", "BT2", "BT2_EXPLICIT", &
1921 "BT1_LS", "BT1_EXPLICIT_LS", "BT2_LS", "BT2_EXPLICIT_LS"), &
1922 enum_desc=s2a("Broyden's first method, also known as the good method. The initial Jacobian"// &
1923 " is built from MD history if available. Otherwise switches to SD for one"// &
1924 " SCF iteration until a Jacobian can be built from the SCF history.", &
1925 "Same as BT1, but computes the explicit Jacobian with finite differences. "// &
1926 "Requires a CDFT SCF procedure to be active.", &
1927 "Same as BT1, but uses Broyden's second method, also known as the bad method.", &
1928 "Same as BT1_EXPLICIT, but using Broyden's second method.", &
1929 "Same as BT1, but uses backtracking line search for optimizing the step size "// &
1930 "(see optimizer NEWTON_LS).", &
1931 "Same as BT1_EXPLICIT, but uses backtracking line search for optimizing the step size.", &
1932 "Same as BT2, but uses backtracking line search for optimizing the step size.", &
1933 "Same as BT2_EXPLICIT, but uses backtracking line search for optimizing the step size."), &
1937 CALL section_add_keyword(section, keyword)
1938 CALL keyword_release(keyword)
1939
1940 CALL keyword_create(keyword, __location__, name="JACOBIAN_TYPE", &
1941 description="Finite difference method used to calculate the inverse Jacobian "// &
1942 "needed by some optimizers. Compatible only with CDFT constraints.", &
1943 usage="JACOBIAN_TYPE FD1", &
1944 default_i_val=jacobian_fd1, &
1945 enum_c_vals=s2a("FD1", "FD1_BACKWARD", "FD2", "FD2_BACKWARD", "FD1_CENTRAL"), &
1946 enum_desc=s2a("First order forward difference (one extra energy evaluation per constraint).", &
1947 "First order backward difference (one extra energy evaluation per constraint).", &
1948 "Second order forward difference (two extra energy evaluations per constraint).", &
1949 "Second order backward difference (two extra energy evaluations per constraint).", &
1950 "First order central difference (two extra energy evaluations per constraint)."), &
1953 CALL section_add_keyword(section, keyword)
1954 CALL keyword_release(keyword)
1955
1956 CALL keyword_create(keyword, __location__, name="JACOBIAN_STEP", &
1957 description="Step size to use in the calculation of the inverse Jacobian with finite differences. "// &
1958 "Expects one value for all constraints, or one value per constraint.", &
1959 usage="JACOBIAN_STEP 5.0E-3 ", n_var=-1, default_r_val=5.0e-3_dp)
1960 CALL section_add_keyword(section, keyword)
1961 CALL keyword_release(keyword)
1962
1963 CALL keyword_create(keyword, __location__, name="JACOBIAN_FREQ", &
1964 description="Defines parameters that control how often the explicit Jacobian is built,"// &
1965 " which is needed by some optimizers. Expects two values. The first value"// &
1966 " determines how many consecutive CDFT SCF iterations should skip a rebuild,"// &
1967 " whereas the latter how many MD steps. The values can be zero (meaning never"// &
1968 " rebuild) or positive. Both values cannot be zero.", &
1969 usage="JACOBIAN_FREQ 1 1", n_var=2, &
1970 default_i_vals=(/1, 1/), type_of_var=integer_t)
1971 CALL section_add_keyword(section, keyword)
1972 CALL keyword_release(keyword)
1973
1974 CALL keyword_create(keyword, __location__, name="JACOBIAN_RESTART", &
1975 description="Restart the inverse Jacobian using the vector defined with keyword JACOBIAN_VECTOR.", &
1976 usage="JACOBIAN_RESTART TRUE", &
1977 default_l_val=.false., lone_keyword_l_val=.true.)
1978 CALL section_add_keyword(section, keyword)
1979 CALL keyword_release(keyword)
1980
1981 CALL keyword_create(keyword, __location__, name="JACOBIAN_VECTOR", &
1982 description="Defines the inverse Jacobian matrix. Useful for restarting calculations. "// &
1983 "Expects n^2 values where n is the total number of constraints. "// &
1984 "The matrix should be given in row major order.", &
1985 usage="JACOBIAN_VECTOR 1.0 0.0", n_var=-1, type_of_var=real_t)
1986 CALL section_add_keyword(section, keyword)
1987 CALL keyword_release(keyword)
1988
1989 CALL keyword_create(keyword, __location__, name="MAX_LS", &
1990 description="The maximum number of backtracking line search steps to perform.", &
1991 usage="MAX_LS 5", default_i_val=5)
1992 CALL section_add_keyword(section, keyword)
1993 CALL keyword_release(keyword)
1994
1995 CALL keyword_create(keyword, __location__, name="FACTOR_LS", &
1996 description="Control parameter for backtracking line search. The step size is reduced by "// &
1997 "this factor on every line search iteration. Value must be between 0 and 1 (exclusive).", &
1998 usage="FACTOR_LS 0.5", default_r_val=0.5_dp)
1999 CALL section_add_keyword(section, keyword)
2000 CALL keyword_release(keyword)
2001
2002 CALL keyword_create(keyword, __location__, name="CONTINUE_LS", &
2003 description="Continue backtracking line search until MAX_LS steps are reached or the "// &
2004 "norm of the CDFT gradient no longer decreases. Default (false) behavior exits the "// &
2005 "line search procedure on the first step that the gradient decreases.", &
2006 usage="CONTINUE_LS TRUE", &
2007 default_l_val=.false., lone_keyword_l_val=.true.)
2008 CALL section_add_keyword(section, keyword)
2009 CALL keyword_release(keyword)
2010
2011 END SUBROUTINE create_cdft_opt_section
2012
2013END MODULE input_cp2k_scf
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public vandevondele2003
integer, save, public holmberg2017
integer, save, public vandevondele2005a
integer, save, public schiffmann2015
integer, save, public weber2008
integer, save, public holmberg2018
integer, save, public becke1988b
integer, save, public stewart1982
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 high_print_level
integer, parameter, public add_last_numeric
subroutine, public cp_print_key_section_create(print_key_section, location, name, description, print_level, each_iter_names, each_iter_values, add_last, filename, common_iter_levels, citations, unit_str)
creates a print_key section
unit conversion facility
Definition cp_units.F:30
real(kind=dp) function, public cp_unit_to_cp2k(value, unit_str, defaults, power)
converts to the internal cp2k units to the given unit
Definition cp_units.F:1150
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public jacobian_fd1_central
integer, parameter, public smear_fermi_dirac
integer, parameter, public core_guess
integer, parameter, public mopac_guess
integer, parameter, public ls_adapt
integer, parameter, public radius_vdw
integer, parameter, public outer_scf_optimizer_sd
integer, parameter, public cholesky_restore
integer, parameter, public ot_chol_irac
integer, parameter, public cdft_beta_constraint
integer, parameter, public broyden_type_2_explicit_ls
integer, parameter, public cdft_magnetization_constraint
integer, parameter, public outer_scf_optimizer_bisect
integer, parameter, public outer_scf_optimizer_secant
integer, parameter, public smear_energy_window
integer, parameter, public ls_3pnt
integer, parameter, public diag_block_krylov
integer, parameter, public becke_cutoff_element
integer, parameter, public outer_scf_cdft_constraint
integer, parameter, public broyden_type_1_explicit
integer, parameter, public no_guess
integer, parameter, public broyden_type_2_ls
integer, parameter, public jacobian_fd2
integer, parameter, public broyden_type_1
integer, parameter, public radius_default
integer, parameter, public outer_scf_optimizer_broyden
integer, parameter, public atomic_guess
integer, parameter, public broyden_type_1_explicit_ls
integer, parameter, public ot_algo_irac
integer, parameter, public outer_scf_basis_center_opt
integer, parameter, public cholesky_dbcsr
integer, parameter, public broyden_type_2_explicit
integer, parameter, public outer_scf_s2_constraint
integer, parameter, public ot_algo_taylor_or_diag
integer, parameter, public ot_poly_irac
integer, parameter, public history_guess
integer, parameter, public cholesky_off
integer, parameter, public cdft_charge_constraint
integer, parameter, public smear_list
integer, parameter, public jacobian_fd1
integer, parameter, public broyden_type_2
integer, parameter, public ot_mini_cg
integer, parameter, public ot_precond_full_kinetic
integer, parameter, public cholesky_reduce
integer, parameter, public high_spin_roks
integer, parameter, public cholesky_inverse
integer, parameter, public ot_mini_diis
integer, parameter, public diag_ot
integer, parameter, public outer_scf_ddapc_constraint
integer, parameter, public ot_precond_solver_default
integer, parameter, public outer_scf_becke_constraint
integer, parameter, public radius_user
integer, parameter, public random_guess
integer, parameter, public ot_precond_full_single
integer, parameter, public ot_precond_solver_inv_chol
integer, parameter, public shape_function_density
integer, parameter, public outer_scf_hirshfeld_constraint
integer, parameter, public radius_covalent
integer, parameter, public jacobian_fd1_backward
integer, parameter, public ot_precond_none
integer, parameter, public jacobian_fd2_backward
integer, parameter, public ls_2pnt
integer, parameter, public ls_none
integer, parameter, public ot_precond_full_single_inverse
integer, parameter, public ot_lwdn_irac
integer, parameter, public diag_filter_matrix
integer, parameter, public ls_gold
integer, parameter, public sparse_guess
integer, parameter, public diag_block_davidson
integer, parameter, public shape_function_gaussian
integer, parameter, public radius_single
integer, parameter, public outer_scf_optimizer_newton_ls
integer, parameter, public gaussian
integer, parameter, public outer_scf_optimizer_none
integer, parameter, public outer_scf_optimizer_newton
integer, parameter, public general_roks
integer, parameter, public outer_scf_none
integer, parameter, public diag_standard
integer, parameter, public becke_cutoff_global
integer, parameter, public cdft_alpha_constraint
integer, parameter, public broyden_type_1_ls
integer, parameter, public restart_guess
integer, parameter, public ot_precond_s_inverse
integer, parameter, public ot_mini_broyden
integer, parameter, public ot_precond_solver_update
integer, parameter, public outer_scf_optimizer_diis
integer, parameter, public ot_mini_sd
integer, parameter, public numerical
integer, parameter, public ot_precond_full_all
integer, parameter, public ot_precond_solver_direct
integer, parameter, public eht_guess
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
subroutine, public create_cdft_control_section(section)
Creates the input section for defining CDFT constraints.
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 integer_t
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
module that contains the definitions of the scf types
subroutine, public create_mixing_section(section, ls_scf)
Create CP2K input section for the mixing of the density matrix to be used only with diagonalization m...
subroutine, public create_filtermatrix_section(section)
Input section for filter matrix diagonalisation method.
Definition qs_fb_input.F:37
manage control variables for the maximum overlap method
subroutine, public create_mom_section(section)
Create CP2K input section for variable occupancy using the Maximum Overlap Method....
Utilities for string manipulations.
character(len=1), parameter, public newline
represent a keyword in the input
represent a section of the input file