(git:374b731)
Loading...
Searching...
No Matches
input_cp2k_ec.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief function that build the dft section of the input
10!> \par History
11!> 10.2005 moved out of input_cp2k [fawzi]
12!> \author fawzi
13! **************************************************************************************************
15 USE bibliography, ONLY: niklasson2003,&
21 USE input_constants, ONLY: &
41 USE input_val_types, ONLY: char_t,&
43 USE kinds, ONLY: dp
44 USE string_utilities, ONLY: s2a
45#include "./base/base_uses.f90"
46
47 IMPLICIT NONE
48 PRIVATE
49
50 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'input_cp2k_ec'
51
52 PUBLIC :: create_ec_section
53
54CONTAINS
55
56! **************************************************************************************************
57!> \brief creates the ENERGY CORRECTION section
58!> \param section ...
59!> \author JGH
60! **************************************************************************************************
61 SUBROUTINE create_ec_section(section)
62 TYPE(section_type), POINTER :: section
63
64 TYPE(keyword_type), POINTER :: keyword
65 TYPE(section_type), POINTER :: subsection
66
67 cpassert(.NOT. ASSOCIATED(section))
68
69 NULLIFY (keyword)
70 CALL section_create(section, __location__, name="ENERGY_CORRECTION", &
71 description="Sets the various options for the Energy Correction", &
72 n_keywords=0, n_subsections=2, repeats=.false.)
73
74 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
75 description="Controls the activation of the energy_correction", &
76 usage="&ENERGY_CORRECTION T", &
77 default_l_val=.false., &
78 lone_keyword_l_val=.true.)
79 CALL section_add_keyword(section, keyword)
80 CALL keyword_release(keyword)
81
82 ! add a special XC section
83 NULLIFY (subsection)
84 CALL create_xc_section(subsection)
85 CALL section_add_subsection(section, subsection)
86 CALL section_release(subsection)
87
88 ! add a section for solver keywords
89 NULLIFY (subsection)
90 CALL create_ec_solver_section(subsection)
91 CALL section_add_subsection(section, subsection)
92 CALL section_release(subsection)
93
94 ! add a print section for properties
95 NULLIFY (subsection)
96 CALL create_ec_print_section(subsection)
97 CALL section_add_subsection(section, subsection)
98 CALL section_release(subsection)
99
100 CALL keyword_create(keyword, __location__, name="ENERGY_FUNCTIONAL", &
101 description="Functional used in energy correction", &
102 usage="ENERGY_FUNCTIONAL HARRIS", &
103 default_i_val=ec_functional_harris, &
104 enum_c_vals=s2a("HARRIS", "DCDFT"), &
105 enum_desc=s2a("Harris functional", &
106 "Density-corrected DFT"), &
107 enum_i_vals=(/ec_functional_harris, ec_functional_dc/))
108 CALL section_add_keyword(section, keyword)
109 CALL keyword_release(keyword)
110
111 CALL keyword_create(keyword, __location__, name="HARRIS_BASIS", &
112 description="Specifies the type of basis to be used for the KG energy correction. "// &
113 "Options are: (1) the default orbital basis (ORBITAL); "// &
114 "(2) the primitive functions of the default orbital basis (PRIMITIVE); "// &
115 "(3) the basis set labeled in Kind section (HARRIS)", &
116 usage="HARRIS_BASIS ORBITAL", &
117 type_of_var=char_t, default_c_val="ORBITAL", n_var=-1)
118 CALL section_add_keyword(section, keyword)
119 CALL keyword_release(keyword)
120
121 CALL keyword_create(keyword, __location__, name="DEBUG_FORCES", &
122 description="Additional output to debug energy correction forces.", &
123 usage="DEBUG_FORCES T", default_l_val=.false., lone_keyword_l_val=.true.)
124 CALL section_add_keyword(section, keyword)
125 CALL keyword_release(keyword)
126 CALL keyword_create(keyword, __location__, name="DEBUG_STRESS", &
127 description="Additional output to debug energy correction forces.", &
128 usage="DEBUG_STRESS T", default_l_val=.false., lone_keyword_l_val=.true.)
129 CALL section_add_keyword(section, keyword)
130 CALL keyword_release(keyword)
131
132 CALL keyword_create(keyword, __location__, name="SKIP_EC", &
133 description="Skip EC calculation if ground-state calculation has not converged.", &
134 usage="SKIP_EC T", default_l_val=.false., lone_keyword_l_val=.true.)
135 CALL section_add_keyword(section, keyword)
136 CALL keyword_release(keyword)
137
138 CALL keyword_create(keyword, __location__, name="MAO", &
139 description="Use modified atomic orbitals (MAO) to solve Harris equation", &
140 usage="MAO T", default_l_val=.false., lone_keyword_l_val=.true.)
141 CALL section_add_keyword(section, keyword)
142 CALL keyword_release(keyword)
143
144 CALL keyword_create(keyword, __location__, name="MAO_MAX_ITER", &
145 description="Maximum iterations in MAO optimization. ", &
146 usage="MAO_MAX_ITER 100 ", default_i_val=0)
147 CALL section_add_keyword(section, keyword)
148 CALL keyword_release(keyword)
149
150 CALL keyword_create(keyword, __location__, name="MAO_EPS_GRAD", &
151 description="Threshold used for MAO iterations. ", &
152 usage="MAO_EPS_GRAD 1.0E-4 ", default_r_val=1.0e-5_dp)
153 CALL section_add_keyword(section, keyword)
154 CALL keyword_release(keyword)
155
156 CALL keyword_create(keyword, __location__, name="MAO_EPS1", &
157 description="Occupation threshold used to determine number of MAOs."// &
158 " KIND section MAO keyword sets the minimum.", &
159 usage="MAO_EPS1 0.001 ", default_r_val=1000.0_dp)
160 CALL section_add_keyword(section, keyword)
161 CALL keyword_release(keyword)
162
163 CALL keyword_create(keyword, __location__, name="MAO_IOLEVEL", &
164 description="Verbosity of MAO output: (0) no output ... (3) fully verbose", &
165 usage="MAO_IOLEVEL 0 ", default_i_val=1)
166 CALL section_add_keyword(section, keyword)
167 CALL keyword_release(keyword)
168
169 CALL keyword_create(keyword, __location__, name="ALGORITHM", &
170 description="Algorithm used to solve KS equation", &
171 usage="ALGORITHM DIAGONALIZATION", &
172 default_i_val=ec_diagonalization, &
173 enum_c_vals=s2a("DIAGONALIZATION", "MATRIX_SIGN", &
174 "TRS4", "TC2", "OT_DIAG"), &
175 enum_desc=s2a("Diagonalization of KS matrix.", &
176 "Matrix Sign algorithm", &
177 "Trace resetting trs4 algorithm", &
178 "Trace resetting tc2 algorithm", &
179 "OT diagonalization"), &
180 enum_i_vals=(/ec_diagonalization, ec_matrix_sign, &
182 CALL section_add_keyword(section, keyword)
183 CALL keyword_release(keyword)
184
185 CALL keyword_create(keyword, __location__, name="FACTORIZATION", &
186 description="Algorithm used to calculate factorization of overlap matrix", &
187 usage="FACTORIZATION CHOLESKY", &
188 default_i_val=kg_cholesky, &
189 enum_c_vals=s2a("CHOLESKY"), &
190 enum_desc=s2a("Cholesky factorization of overlap matrix"), &
191 enum_i_vals=(/kg_cholesky/))
192 CALL section_add_keyword(section, keyword)
193 CALL keyword_release(keyword)
194
195 CALL keyword_create(keyword, __location__, name="EPS_DEFAULT", &
196 description="Threshold used for accuracy estimates within energy correction. ", &
197 usage="EPS_DEFAULT 1.0E-7 ", default_r_val=1.0e-7_dp)
198 CALL section_add_keyword(section, keyword)
199 CALL keyword_release(keyword)
200
201 ! Keywords for LS solver of Harris functional
202 CALL keyword_create(keyword, __location__, name="EPS_FILTER", &
203 description="Threshold used for filtering matrix operations.", &
204 usage="EPS_FILTER 1.0E-12", default_r_val=1.0e-12_dp)
205 CALL section_add_keyword(section, keyword)
206 CALL keyword_release(keyword)
207
208 CALL keyword_create(keyword, __location__, name="EPS_LANCZOS", &
209 description="Threshold used for lanczos estimates.", &
210 usage="EPS_LANCZOS 1.0E-4", default_r_val=1.0e-3_dp)
211 CALL section_add_keyword(section, keyword)
212 CALL keyword_release(keyword)
213
214 CALL keyword_create(keyword, __location__, name="MAX_ITER_LANCZOS", &
215 description="Maximum number of lanczos iterations.", &
216 usage="MAX_ITER_LANCZOS ", default_i_val=128)
217 CALL section_add_keyword(section, keyword)
218 CALL keyword_release(keyword)
219
220 CALL keyword_create(keyword, __location__, name="MU", &
221 description="Value (or initial guess) for the chemical potential,"// &
222 " i.e. some suitable energy between HOMO and LUMO energy.", &
223 usage="MU 0.0", default_r_val=-0.1_dp)
224 CALL section_add_keyword(section, keyword)
225 CALL keyword_release(keyword)
226
227 CALL keyword_create(keyword, __location__, name="FIXED_MU", &
228 description="Should the calculation be performed at fixed chemical potential,"// &
229 " or should it be found fixing the number of electrons", &
230 usage="FIXED_MU .TRUE.", default_l_val=.false., lone_keyword_l_val=.true.)
231 CALL section_add_keyword(section, keyword)
232 CALL keyword_release(keyword)
233
234 CALL keyword_create(keyword, __location__, name="S_PRECONDITIONER", &
235 description="Preconditions S with some appropriate form.", &
236 usage="S_PRECONDITIONER MOLECULAR", &
237 default_i_val=ls_s_preconditioner_atomic, &
238 enum_c_vals=s2a("NONE", "ATOMIC", "MOLECULAR"), &
239 enum_desc=s2a("No preconditioner", &
240 "Using atomic blocks", &
241 "Using molecular sub-blocks. Recommended if molecules are defined and not too large."), &
243 CALL section_add_keyword(section, keyword)
244 CALL keyword_release(keyword)
245
246 CALL keyword_create(keyword, __location__, name="S_SQRT_METHOD", &
247 description="Method for the caclulation of the sqrt of S.", &
248 usage="S_SQRT_METHOD NEWTONSCHULZ", &
249 default_i_val=ls_s_sqrt_ns, &
250 enum_c_vals=s2a("NEWTONSCHULZ", "PROOT"), &
251 enum_desc=s2a("Using a Newton-Schulz-like iteration", &
252 "Using the p-th root method."), &
253 enum_i_vals=(/ls_s_sqrt_ns, ls_s_sqrt_proot/))
254 CALL section_add_keyword(section, keyword)
255 CALL keyword_release(keyword)
256
257 CALL keyword_create(keyword, __location__, name="S_SQRT_ORDER", &
258 variants=s2a("SIGN_SQRT_ORDER"), &
259 description="Order of the iteration method for the calculation of the sqrt of S.", &
260 usage="S_SQRT_ORDER 3", default_i_val=3)
261 CALL section_add_keyword(section, keyword)
262 CALL keyword_release(keyword)
263
264 CALL keyword_create(keyword, __location__, name="SIGN_METHOD", &
265 description="Method used for the computation of the sign matrix.", &
266 usage="SIGN_METHOD NEWTONSCHULZ", &
267 default_i_val=ls_scf_sign_ns, &
268 citations=(/vandevondele2012, niklasson2003/), &
269 enum_c_vals=s2a("NEWTONSCHULZ", "PROOT"), &
270 enum_desc=s2a("Newton-Schulz iteration.", &
271 "p-th order root iteration"), &
272 enum_i_vals=(/ls_scf_sign_ns, ls_scf_sign_proot/))
273 CALL section_add_keyword(section, keyword)
274 CALL keyword_release(keyword)
275
276 CALL keyword_create(keyword, __location__, name="SIGN_ORDER", &
277 description="Order of the method used for the computation of the sign matrix.", &
278 usage="SIGN_ORDER 2", &
279 default_i_val=2)
280 CALL section_add_keyword(section, keyword)
281 CALL keyword_release(keyword)
282
283 CALL keyword_create(keyword, __location__, name="DYNAMIC_THRESHOLD", &
284 description="Should the threshold for the purification be chosen dynamically", &
285 usage="DYNAMIC_THRESHOLD .TRUE.", default_l_val=.false., lone_keyword_l_val=.true.)
286 CALL section_add_keyword(section, keyword)
287 CALL keyword_release(keyword)
288
289 CALL keyword_create(keyword, __location__, name="NON_MONOTONIC", &
290 description="Should the purification be performed non-monotonically. Relevant for TC2 only.", &
291 usage="NON_MONOTONIC .TRUE.", default_l_val=.true., lone_keyword_l_val=.true.)
292 CALL section_add_keyword(section, keyword)
293 CALL keyword_release(keyword)
294
295 CALL keyword_create( &
296 keyword, __location__, name="MATRIX_CLUSTER_TYPE", &
297 description="Specify how atomic blocks should be clustered in the used matrices, in order to improve flop rate, "// &
298 "and possibly speedup the matrix multiply. Note that the atomic s_preconditioner can not be used. "// &
299 "Furthermore, since screening is on matrix blocks, "// &
300 "slightly more accurate results can be expected with molecular.", &
301 usage="MATRIX_CLUSTER_TYPE MOLECULAR", &
302 default_i_val=ls_cluster_atomic, &
303 enum_c_vals=s2a("ATOMIC", "MOLECULAR"), &
304 enum_desc=s2a("Using atomic blocks", &
305 "Using molecular blocks."), &
307 CALL section_add_keyword(section, keyword)
308 CALL keyword_release(keyword)
309
310 CALL keyword_create(keyword, __location__, name="S_INVERSION", &
311 description="Method used to compute the inverse of S.", &
312 usage="S_PRECONDITIONER MOLECULAR", &
313 default_i_val=ls_s_inversion_sign_sqrt, &
314 enum_c_vals=s2a("SIGN_SQRT", "HOTELLING"), &
315 enum_desc=s2a("Using the inverse sqrt as obtained from sign function iterations.", &
316 "Using the Hotellign iteration."), &
318 CALL section_add_keyword(section, keyword)
319 CALL keyword_release(keyword)
320
321 CALL keyword_create(keyword, __location__, name="REPORT_ALL_SPARSITIES", &
322 description="Run the sparsity report at the end of the SCF", &
323 usage="REPORT_ALL_SPARSITIES", default_l_val=.false., lone_keyword_l_val=.true.)
324 CALL section_add_keyword(section, keyword)
325 CALL keyword_release(keyword)
326
327 CALL keyword_create(keyword, __location__, name="CHECK_S_INV", &
328 description="Perform an accuracy check on the inverse/sqrt of the s matrix.", &
329 usage="CHECK_S_INV", default_l_val=.false., lone_keyword_l_val=.true.)
330 CALL section_add_keyword(section, keyword)
331 CALL keyword_release(keyword)
332
333 CALL keyword_create(keyword, __location__, name="OT_INITIAL_GUESS", &
334 description="Initial guess of density matrix used for OT Diagonalization", &
335 usage="OT_INITIAL_GUESS ATOMIC", &
336 default_i_val=ec_ot_atomic, &
337 enum_c_vals=s2a("ATOMIC", "GROUND_STATE"), &
338 enum_desc=s2a("Generate an atomic density using the atomic code", &
339 "Using the ground-state density."), &
340 enum_i_vals=(/ec_ot_atomic, ec_ot_gs/))
341 CALL section_add_keyword(section, keyword)
342 CALL keyword_release(keyword)
343
344 CALL keyword_create( &
345 keyword, __location__, &
346 name="ADMM", &
347 description="Decide whether to perform ADMM in the exact exchange calc. for DC-DFT. "// &
348 "The ADMM XC correction is governed by the AUXILIARY_DENSITY_MATRIX_METHOD section in &DFT. "// &
349 "In most cases, the Hartree-Fock exchange is not too expensive and there is no need for ADMM, "// &
350 "ADMM can however provide significant speedup and memory savings in case of diffuse basis sets. ", &
351 usage="ADMM", &
352 default_l_val=.false., &
353 lone_keyword_l_val=.true.)
354 CALL section_add_keyword(section, keyword)
355 CALL keyword_release(keyword)
356
357 END SUBROUTINE create_ec_section
358
359! **************************************************************************************************
360!> \brief creates the linear scaling solver section
361!> \param section ...
362!> \author Joost VandeVondele [2010-10], JGH [2019-12]
363! **************************************************************************************************
364 SUBROUTINE create_ec_solver_section(section)
365 TYPE(section_type), POINTER :: section
366
367 TYPE(keyword_type), POINTER :: keyword
368
369 cpassert(.NOT. ASSOCIATED(section))
370 CALL section_create(section, __location__, name="RESPONSE_SOLVER", &
371 description="Specifies the parameters of the linear scaling solver routines", &
372 n_keywords=24, n_subsections=3, repeats=.false., &
373 citations=(/vandevondele2012/))
374
375 NULLIFY (keyword)
376
377 CALL keyword_create(keyword, __location__, name="EPS", &
378 description="Target accuracy for the convergence of the conjugate gradient.", &
379 usage="EPS 1.e-10", default_r_val=1.e-12_dp)
380 CALL section_add_keyword(section, keyword)
381 CALL keyword_release(keyword)
382
383 CALL keyword_create(keyword, __location__, name="EPS_FILTER", &
384 description="Threshold used for filtering matrix operations.", &
385 usage="EPS_FILTER 1.0E-10", default_r_val=1.0e-10_dp)
386 CALL section_add_keyword(section, keyword)
387 CALL keyword_release(keyword)
388
389 CALL keyword_create(keyword, __location__, name="EPS_LANCZOS", &
390 description="Threshold used for lanczos estimates.", &
391 usage="EPS_LANCZOS 1.0E-4", default_r_val=1.0e-3_dp)
392 CALL section_add_keyword(section, keyword)
393 CALL keyword_release(keyword)
394
395 CALL keyword_create(keyword, __location__, name="MAX_ITER", &
396 description="Maximum number of conjugate gradient iteration "// &
397 "to be performed for one optimization.", &
398 usage="MAX_ITER 200", default_i_val=50)
399 CALL section_add_keyword(section, keyword)
400 CALL keyword_release(keyword)
401
402 CALL keyword_create(keyword, __location__, name="MAX_ITER_LANCZOS", &
403 description="Maximum number of lanczos iterations.", &
404 usage="MAX_ITER_LANCZOS 128", default_i_val=128)
405 CALL section_add_keyword(section, keyword)
406 CALL keyword_release(keyword)
407
408 CALL keyword_create(keyword, __location__, name="METHOD", &
409 description="Algorithm used to solve response equation. "// &
410 "Both solver are conjugate gradient based, but use either a vector (MO-coefficient) "// &
411 "or density matrix formalism in the orthonormal AO-basis to obtain response density", &
412 usage="METHOD SOLVER", &
413 default_i_val=ec_ls_solver, &
414 enum_c_vals=s2a("MO_SOLVER", "AO_ORTHO"), &
415 enum_desc=s2a("Solver based on MO (vector) formalism", &
416 "Solver based on density matrix formalism"), &
417 enum_i_vals=(/ec_mo_solver, ec_ls_solver/))
418 CALL section_add_keyword(section, keyword)
419 CALL keyword_release(keyword)
420
421 CALL keyword_create( &
422 keyword, __location__, name="PRECONDITIONER", &
423 description="Type of preconditioner to be used with MO conjugate gradient solver. "// &
424 "They differ in effectiveness, cost of construction, cost of application. "// &
425 "Properly preconditioned minimization can be orders of magnitude faster than doing nothing. "// &
426 "Only multi-level conjugate gradient preconditioner (MULTI_LEVEL) available for AO response solver (AO_ORTHO). ", &
427 usage="PRECONDITIONER FULL_ALL", &
428 default_i_val=precond_mlp, &
429 enum_c_vals=s2a("FULL_ALL", "FULL_SINGLE_INVERSE", "FULL_SINGLE", "FULL_KINETIC", "FULL_S_INVERSE", &
430 "MULTI_LEVEL", "NONE"), &
431 enum_desc=s2a("Most effective state selective preconditioner based on diagonalization, "// &
432 "requires the ENERGY_GAP parameter to be an underestimate of the HOMO-LUMO gap. "// &
433 "This preconditioner is recommended for almost all systems, except very large systems where "// &
434 "make_preconditioner would dominate the total computational cost.", &
435 "Based on H-eS cholesky inversion, similar to FULL_SINGLE in preconditioning efficiency "// &
436 "but cheaper to construct, "// &
437 "might be somewhat less robust. Recommended for large systems.", &
438 "Based on H-eS diagonalisation, not as good as FULL_ALL, but somewhat cheaper to apply. ", &
439 "Cholesky inversion of S and T, fast construction, robust, and relatively good, "// &
440 "use for very large systems.", &
441 "Cholesky inversion of S, not as good as FULL_KINETIC, yet equally expensive.", &
442 "Based on same CG as AO-solver itself, but uses cheaper linear transformation", &
443 "skip preconditioning"), &
446 CALL section_add_keyword(section, keyword)
447 CALL keyword_release(keyword)
448
449 CALL keyword_create(keyword, __location__, name="S_PRECONDITIONER", &
450 description="Preconditions S with some appropriate form.", &
451 usage="S_PRECONDITIONER MOLECULAR", &
452 default_i_val=ls_s_preconditioner_atomic, &
453 enum_c_vals=s2a("NONE", "ATOMIC", "MOLECULAR"), &
454 enum_desc=s2a("No preconditioner", &
455 "Using atomic blocks", &
456 "Using molecular sub-blocks. Recommended if molecules are defined and not too large."), &
458 CALL section_add_keyword(section, keyword)
459 CALL keyword_release(keyword)
460
461 CALL keyword_create(keyword, __location__, name="S_SQRT_METHOD", &
462 description="Method for the caclulation of the sqrt of S.", &
463 usage="S_SQRT_METHOD NEWTONSCHULZ", &
464 default_i_val=ls_s_sqrt_ns, &
465 enum_c_vals=s2a("NEWTONSCHULZ", "PROOT"), &
466 enum_desc=s2a("Using a Newton-Schulz-like iteration", &
467 "Using the p-th root method."), &
468 enum_i_vals=(/ls_s_sqrt_ns, ls_s_sqrt_proot/))
469 CALL section_add_keyword(section, keyword)
470 CALL keyword_release(keyword)
471
472 CALL keyword_create(keyword, __location__, name="S_SQRT_ORDER", &
473 variants=s2a("SIGN_SQRT_ORDER"), &
474 description="Order of the iteration method for the calculation of the sqrt of S.", &
475 usage="S_SQRT_ORDER 3", default_i_val=3)
476 CALL section_add_keyword(section, keyword)
477 CALL keyword_release(keyword)
478
479 CALL keyword_create( &
480 keyword, __location__, name="MATRIX_CLUSTER_TYPE", &
481 description="Specify how atomic blocks should be clustered in the used matrices, in order to improve flop rate, "// &
482 "and possibly speedup the matrix multiply. Note that the atomic s_preconditioner can not be used. "// &
483 "Furthermore, since screening is on matrix blocks, "// &
484 "slightly more accurate results can be expected with molecular.", &
485 usage="MATRIX_CLUSTER_TYPE MOLECULAR", &
486 default_i_val=ls_cluster_atomic, &
487 enum_c_vals=s2a("ATOMIC", "MOLECULAR"), &
488 enum_desc=s2a("Using atomic blocks", &
489 "Using molecular blocks."), &
491 CALL section_add_keyword(section, keyword)
492 CALL keyword_release(keyword)
493
494 CALL keyword_create(keyword, __location__, name="S_INVERSION", &
495 description="Method used to compute the inverse of S.", &
496 usage="S_PRECONDITIONER MOLECULAR", &
497 default_i_val=ls_s_inversion_sign_sqrt, &
498 enum_c_vals=s2a("SIGN_SQRT", "HOTELLING"), &
499 enum_desc=s2a("Using the inverse sqrt as obtained from sign function iterations.", &
500 "Using the Hotellign iteration."), &
502 CALL section_add_keyword(section, keyword)
503 CALL keyword_release(keyword)
504
505 CALL keyword_create(keyword, __location__, name="RESTART", &
506 description="Restart the response calculation if the restart file exists", &
507 usage="RESTART", &
508 default_l_val=.false., lone_keyword_l_val=.true.)
509 CALL section_add_keyword(section, keyword)
510 CALL keyword_release(keyword)
511
512 CALL keyword_create(keyword, __location__, name="RESTART_EVERY", &
513 description="Restart the conjugate gradient after the specified number of iterations.", &
514 usage="RESTART_EVERY 50", default_i_val=50)
515 CALL section_add_keyword(section, keyword)
516 CALL keyword_release(keyword)
517
518 END SUBROUTINE create_ec_solver_section
519
520! **************************************************************************************************
521!> \brief Create the print dft section
522!> \param section the section to create
523!> \author fbelle - from create_print_dft_section
524! **************************************************************************************************
525 SUBROUTINE create_ec_print_section(section)
526 TYPE(section_type), POINTER :: section
527
528 TYPE(keyword_type), POINTER :: keyword
529 TYPE(section_type), POINTER :: print_key
530
531 cpassert(.NOT. ASSOCIATED(section))
532 CALL section_create(section, __location__, name="PRINT", &
533 description="Section of possible print options in EC code.", &
534 n_keywords=0, n_subsections=1, repeats=.false.)
535
536 NULLIFY (print_key, keyword)
537
538 ! Output of BQB volumetric files
539 CALL cp_print_key_section_create(print_key, __location__, name="E_DENSITY_BQB", &
540 description="Controls the output of the electron density to the losslessly"// &
541 " compressed BQB file format, see [Brehm2018]"// &
542 " (via LibBQB see <https://brehm-research.de/bqb>)."// &
543 " Currently does not work with changing cell vector (NpT ensemble).", &
544 print_level=debug_print_level + 1, filename="", &
545 citations=(/brehm2018/))
546
547 CALL keyword_create(keyword, __location__, name="SKIP_FIRST", &
548 description="Skips the first step of a MD run (avoids duplicate step if restarted).", &
549 usage="SKIP_FIRST T", default_l_val=.false., lone_keyword_l_val=.true.)
550 CALL section_add_keyword(print_key, keyword)
551 CALL keyword_release(keyword)
552
553 CALL keyword_create(keyword, __location__, name="STORE_STEP_NUMBER", &
554 description="Stores the step number and simulation time in the comment line of each BQB"// &
555 " frame. Switch it off for binary compatibility with original CP2k CUBE files.", &
556 usage="STORE_STEP_NUMBER F", default_l_val=.true., lone_keyword_l_val=.true.)
557 CALL section_add_keyword(print_key, keyword)
558 CALL keyword_release(keyword)
559
560 CALL keyword_create(keyword, __location__, name="CHECK", &
561 description="Performs an on-the-fly decompression of each compressed BQB frame to check"// &
562 " whether the volumetric data exactly matches, and aborts the run if not so.", &
563 usage="CHECK T", default_l_val=.false., lone_keyword_l_val=.true.)
564 CALL section_add_keyword(print_key, keyword)
565 CALL keyword_release(keyword)
566
567 CALL keyword_create(keyword, __location__, name="OVERWRITE", &
568 description="Specify this keyword to overwrite the output BQB file if"// &
569 " it already exists. By default, the data is appended to an existing file.", &
570 usage="OVERWRITE T", default_l_val=.false., lone_keyword_l_val=.true.)
571 CALL section_add_keyword(print_key, keyword)
572 CALL keyword_release(keyword)
573
574 CALL keyword_create(keyword, __location__, name="HISTORY", &
575 description="Controls how many previous steps are taken into account for extrapolation in"// &
576 " compression. Use a value of 1 to compress the frames independently.", &
577 usage="HISTORY 10", n_var=1, default_i_val=10, type_of_var=integer_t)
578 CALL section_add_keyword(print_key, keyword)
579 CALL keyword_release(keyword)
580
581 CALL keyword_create(keyword, __location__, name="PARAMETER_KEY", &
582 description="Allows to supply previously optimized compression parameters via a"// &
583 " parameter key (alphanumeric character sequence starting with 'at')."// &
584 " Just leave away the 'at' sign here, because CP2k will otherwise"// &
585 " assume it is a variable name in the input", &
586 usage="PARAMETER_KEY <KEY>", n_var=1, default_c_val="", type_of_var=char_t)
587 CALL section_add_keyword(print_key, keyword)
588 CALL keyword_release(keyword)
589
590 CALL keyword_create(keyword, __location__, name="OPTIMIZE", &
591 description="Controls the time spent to optimize the parameters for compression efficiency.", &
592 usage="OPTIMIZE {OFF,QUICK,NORMAL,PATIENT,EXHAUSTIVE}", repeats=.false., n_var=1, &
593 default_i_val=bqb_opt_quick, &
594 enum_c_vals=s2a("OFF", "QUICK", "NORMAL", "PATIENT", "EXHAUSTIVE"), &
595 enum_desc=s2a("No optimization (use defaults)", "Quick optimization", &
596 "Standard optimization", "Precise optimization", "Exhaustive optimization"), &
598 CALL section_add_keyword(print_key, keyword)
599 CALL keyword_release(keyword)
600
601 CALL section_add_subsection(section, print_key)
602 CALL section_release(print_key)
603
604 ! Voronoi Integration via LibVori
605 NULLIFY (print_key)
606 CALL create_print_voronoi_section(print_key)
607 CALL section_add_subsection(section, print_key)
608 CALL section_release(print_key)
609
610 !Printing of Moments
611 CALL create_dipoles_section(print_key, "MOMENTS", high_print_level)
612 CALL keyword_create( &
613 keyword, __location__, &
614 name="MAX_MOMENT", &
615 description="Maximum moment to be calculated. Values higher than 1 not implemented under periodic boundaries.", &
616 usage="MAX_MOMENT {integer}", &
617 repeats=.false., &
618 n_var=1, &
619 type_of_var=integer_t, &
620 default_i_val=1)
621 CALL section_add_keyword(print_key, keyword)
622 CALL keyword_release(keyword)
623 CALL keyword_create(keyword, __location__, &
624 name="MAGNETIC", &
625 description="Calculate also magnetic moments, only implemented without periodic boundaries", &
626 usage="MAGNETIC yes", &
627 repeats=.false., &
628 n_var=1, &
629 default_l_val=.false., &
630 lone_keyword_l_val=.true.)
631 CALL section_add_keyword(print_key, keyword)
632 CALL keyword_release(keyword)
633 CALL section_add_subsection(section, print_key)
634 CALL section_release(print_key)
635
636 END SUBROUTINE create_ec_print_section
637
638END MODULE input_cp2k_ec
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public vandevondele2012
integer, save, public niklasson2003
integer, save, public brehm2018
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer, parameter, public debug_print_level
integer, parameter, public high_print_level
subroutine, public cp_print_key_section_create(print_key_section, location, name, description, print_level, each_iter_names, each_iter_values, add_last, filename, common_iter_levels, citations, unit_str)
creates a print_key section
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public ec_functional_harris
integer, parameter, public bqb_opt_quick
integer, parameter, public ec_functional_dc
integer, parameter, public ls_s_preconditioner_molecular
integer, parameter, public precond_mlp
integer, parameter, public ls_s_inversion_hotelling
integer, parameter, public ls_cluster_molecular
integer, parameter, public ec_ot_atomic
integer, parameter, public ls_s_preconditioner_atomic
integer, parameter, public ec_mo_solver
integer, parameter, public bqb_opt_normal
integer, parameter, public ot_precond_full_kinetic
integer, parameter, public bqb_opt_off
integer, parameter, public ec_ot_diag
integer, parameter, public ec_ot_gs
integer, parameter, public bqb_opt_exhaustive
integer, parameter, public ls_s_sqrt_proot
integer, parameter, public ot_precond_full_single
integer, parameter, public ls_s_sqrt_ns
integer, parameter, public ls_s_preconditioner_none
integer, parameter, public ls_scf_sign_proot
integer, parameter, public ot_precond_none
integer, parameter, public ot_precond_full_single_inverse
integer, parameter, public ec_matrix_tc2
integer, parameter, public ec_diagonalization
integer, parameter, public ec_matrix_trs4
integer, parameter, public ls_scf_sign_ns
integer, parameter, public ec_matrix_sign
integer, parameter, public ls_s_inversion_sign_sqrt
integer, parameter, public bqb_opt_patient
integer, parameter, public ot_precond_s_inverse
integer, parameter, public kg_cholesky
integer, parameter, public ot_precond_full_all
integer, parameter, public ls_cluster_atomic
integer, parameter, public ec_ls_solver
function that build the dft section of the input
subroutine, public create_ec_section(section)
creates the ENERGY CORRECTION section
creates the mm section of the input
subroutine, public create_dipoles_section(print_key, label, print_level)
creates the input section for the qs part
function that build the dft section of the input
subroutine, public create_print_voronoi_section(print_key)
Create the print voronoi section.
function that build the xc section of the input
subroutine, public create_xc_section(section)
creates the input section for the xc part
represents keywords in an input
subroutine, public keyword_release(keyword)
releases the given keyword (see doc/ReferenceCounting.html)
subroutine, public keyword_create(keyword, location, name, description, usage, type_of_var, n_var, repeats, variants, default_val, default_l_val, default_r_val, default_lc_val, default_c_val, default_i_val, default_l_vals, default_r_vals, default_c_vals, default_i_vals, lone_keyword_val, lone_keyword_l_val, lone_keyword_r_val, lone_keyword_c_val, lone_keyword_i_val, lone_keyword_l_vals, lone_keyword_r_vals, lone_keyword_c_vals, lone_keyword_i_vals, enum_c_vals, enum_i_vals, enum, enum_strict, enum_desc, unit_str, citations, deprecation_notice, removed)
creates a keyword object
objects that represent the structure of input sections and the data contained in an input section
subroutine, public section_create(section, location, name, description, n_keywords, n_subsections, repeats, citations)
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 char_t
integer, parameter, public integer_t
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
Utilities for string manipulations.
represent a keyword in the input
represent a section of the input file