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