(git:5f3bc36)
Loading...
Searching...
No Matches
input_cp2k_properties_dft.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 dft section of the input
10!> \par History
11!> 01.2013 moved out of input_cp2k_dft [MI]
12!> \author MI
13! **************************************************************************************************
15 USE bibliography, ONLY: futera2017, &
17 kondov2007, &
18 luber2014, &
22 weber2009, &
31 USE cp_units, ONLY: cp_unit_to_cp2k
32 USE input_constants, ONLY: &
60 USE input_val_types, ONLY: char_t, &
61 integer_t, &
62 lchar_t, &
63 logical_t, &
64 real_t
66 USE kinds, ONLY: dp
67 USE string_utilities, ONLY: s2a
68#include "./base/base_uses.f90"
69
70 IMPLICIT NONE
71 PRIVATE
72
73 LOGICAL, PRIVATE, PARAMETER :: debug_this_module = .true.
74 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'input_cp2k_properties_dft'
75
77
78CONTAINS
79
80! **************************************************************************************************
81!> \brief Create the PROPERTIES section
82!> \param section the section to create
83!> \author teo
84! **************************************************************************************************
85 SUBROUTINE create_properties_section(section)
86 TYPE(section_type), POINTER :: section
87
88 TYPE(keyword_type), POINTER :: keyword
89 TYPE(section_type), POINTER :: subsection
90
91 cpassert(.NOT. ASSOCIATED(section))
92 CALL section_create(section, __location__, name="PROPERTIES", &
93 description="This section is used to set up the PROPERTIES calculation.", &
94 n_keywords=0, n_subsections=6, repeats=.false.)
95
96 NULLIFY (subsection, keyword)
97
98 CALL create_linres_section(subsection, create_subsections=.true.)
99 CALL section_add_subsection(section, subsection)
100 CALL section_release(subsection)
101
102 CALL create_et_coupling_section(subsection)
103 CALL section_add_subsection(section, subsection)
104 CALL section_release(subsection)
105
106 CALL create_resp_section(subsection)
107 CALL section_add_subsection(section, subsection)
108 CALL section_release(subsection)
109
110 CALL create_atprop_section(subsection)
111 CALL section_add_subsection(section, subsection)
112 CALL section_release(subsection)
113
114 CALL cp_print_key_section_create(subsection, __location__, name="FIT_CHARGE", &
115 description="This section is used to print the density derived atomic point charges. "// &
116 "The fit of the charges is controlled through the DENSITY_FITTING section", &
117 print_level=high_print_level, filename="__STD_OUT__")
118 CALL keyword_create(keyword, __location__, name="TYPE_OF_DENSITY", &
119 description="Specifies the type of density used for the fitting", &
120 usage="TYPE_OF_DENSITY (FULL|SPIN)", &
121 enum_c_vals=s2a("FULL", "SPIN"), &
122 enum_i_vals=(/do_full_density, do_spin_density/), &
123 enum_desc=s2a("Full density", "Spin density"), &
124 default_i_val=do_full_density)
125 CALL section_add_keyword(subsection, keyword)
126 CALL keyword_release(keyword)
127 CALL section_add_subsection(section, subsection)
128 CALL section_release(subsection)
129
130 CALL create_tddfpt2_section(subsection)
131 CALL section_add_subsection(section, subsection)
132 CALL section_release(subsection)
133
134 CALL create_rixs_section(subsection)
135 CALL section_add_subsection(section, subsection)
136 CALL section_release(subsection)
137
138 CALL create_bandstructure_section(subsection)
139 CALL section_add_subsection(section, subsection)
140 CALL section_release(subsection)
141
142 CALL create_tipscan_section(subsection)
143 CALL section_add_subsection(section, subsection)
144 CALL section_release(subsection)
145
146 END SUBROUTINE create_properties_section
147
148! **************************************************************************************************
149!> \brief creates the input structure used to activate
150!> a resonant inelastic xray scattering (RIXS) calculation
151! **************************************************************************************************
152
153 SUBROUTINE create_rixs_section(section)
154 TYPE(section_type), POINTER :: section
155 TYPE(section_type), POINTER :: subsection, print_key
156 TYPE(keyword_type), POINTER :: keyword
157
158 cpassert(.NOT. ASSOCIATED(section))
159
160 NULLIFY (keyword, subsection, print_key)
161
162 CALL section_create(section, __location__, name="RIXS", &
163 description="Resonant Inelastic Xray Scattering using XAS_TDP and TDDFPT.", &
164 n_keywords=1, n_subsections=3, repeats=.false., &
165 citations=(/vazdacruz2021/))
166
167 CALL create_tddfpt2_section(subsection)
168 CALL section_add_subsection(section, subsection)
169 CALL section_release(subsection)
170
171 CALL create_xas_tdp_section(subsection)
172 CALL section_add_subsection(section, subsection)
173 CALL section_release(subsection)
174
175 CALL section_create(subsection, __location__, "PRINT", "Controls the printing of information "// &
176 "during RIXS calculations", repeats=.false.)
177
178 CALL cp_print_key_section_create(print_key, __location__, name="SPECTRUM", &
179 description="Controles the printing of the RIXS spectrum "// &
180 "in output files", &
181 print_level=low_print_level, filename="", &
182 common_iter_levels=3)
183 CALL section_add_subsection(subsection, print_key)
184 CALL section_release(print_key)
185
186 CALL section_add_subsection(section, subsection)
187 CALL section_release(subsection)
188
189 END SUBROUTINE create_rixs_section
190
191! **************************************************************************************************
192!> \brief creates the input structure used to activate
193!> a linear response calculation
194!> Available properties : none
195!> \param section the section to create
196!> \param create_subsections indicates whether or not subsections should be created
197!> \param default_set_tdlr default parameters to be used if called from TDDFPT
198!> \author MI
199! **************************************************************************************************
200 SUBROUTINE create_linres_section(section, create_subsections, default_set_tdlr)
201 TYPE(section_type), POINTER :: section
202 LOGICAL, INTENT(in) :: create_subsections
203 LOGICAL, INTENT(IN), OPTIONAL :: default_set_tdlr
204
205 INTEGER :: def_max_iter, def_precond
206 REAL(kind=dp) :: def_egap, def_eps, def_eps_filter
207 TYPE(keyword_type), POINTER :: keyword
208 TYPE(section_type), POINTER :: print_key, subsection
209
210 CHARACTER(len=256) :: desc
211
212 NULLIFY (keyword, print_key)
213
214 IF (PRESENT(default_set_tdlr)) THEN
215 def_egap = 0.02_dp
216 def_eps = 1.0e-10_dp
217 def_eps_filter = 1.0e-15_dp
218 def_max_iter = 100
220 desc = "Controls the parameters of the LINRES force calculations for excited states."
221 ELSE
222 def_egap = 0.2_dp
223 def_eps = 1.e-6_dp
224 def_eps_filter = 0.0_dp
225 def_max_iter = 50
226 def_precond = ot_precond_none
227 desc = "The linear response is used to calculate one of the following properties: nmr, epr, raman, ..."
228 END IF
229
230 cpassert(.NOT. ASSOCIATED(section))
231 CALL section_create(section, __location__, name="linres", &
232 description=desc, n_keywords=5, n_subsections=2, repeats=.false., &
233 citations=(/putrino2000/))
234
235 CALL keyword_create(keyword, __location__, name="EPS", &
236 description="target accuracy for the convergence of the conjugate gradient.", &
237 usage="EPS 1.e-6", default_r_val=def_eps)
238 CALL section_add_keyword(section, keyword)
239 CALL keyword_release(keyword)
240
241 CALL keyword_create(keyword, __location__, name="EPS_FILTER", &
242 description="Filter threshold for response density matrix.", &
243 usage="EPS_FILTER 1.e-8", default_r_val=def_eps_filter)
244 CALL section_add_keyword(section, keyword)
245 CALL keyword_release(keyword)
246
247 CALL keyword_create(keyword, __location__, name="MAX_ITER", &
248 description="Maximum number of conjugate gradient iteration to be performed for one optimization.", &
249 usage="MAX_ITER 200", default_i_val=def_max_iter)
250 CALL section_add_keyword(section, keyword)
251 CALL keyword_release(keyword)
252
253 CALL keyword_create(keyword, __location__, name="RESTART_EVERY", &
254 description="Restart the conjugate gradient after the specified number of iterations.", &
255 usage="RESTART_EVERY 200", default_i_val=50)
256 CALL section_add_keyword(section, keyword)
257 CALL keyword_release(keyword)
258
259 CALL keyword_create( &
260 keyword, __location__, name="PRECONDITIONER", &
261 description="Type of preconditioner to be used with all minimization schemes. "// &
262 "They differ in effectiveness, cost of construction, cost of application. "// &
263 "Properly preconditioned minimization can be orders of magnitude faster than doing nothing.", &
264 usage="PRECONDITIONER FULL_ALL", &
265 default_i_val=def_precond, &
266 enum_c_vals=s2a("FULL_ALL", "FULL_SINGLE_INVERSE", "FULL_SINGLE", "FULL_KINETIC", "FULL_S_INVERSE", &
267 "NONE"), &
268 enum_desc=s2a("Most effective state selective preconditioner based on diagonalization, "// &
269 "requires the ENERGY_GAP parameter to be an underestimate of the HOMO-LUMO gap. "// &
270 "This preconditioner is recommended for almost all systems, except very large systems where "// &
271 "make_preconditioner would dominate the total computational cost.", &
272 "Based on H-eS cholesky inversion, similar to FULL_SINGLE in preconditioning efficiency "// &
273 "but cheaper to construct, "// &
274 "might be somewhat less robust. Recommended for large systems.", &
275 "Based on H-eS diagonalisation, not as good as FULL_ALL, but somewhat cheaper to apply. ", &
276 "Cholesky inversion of S and T, fast construction, robust, and relatively good, "// &
277 "use for very large systems.", &
278 "Cholesky inversion of S, not as good as FULL_KINETIC, yet equally expensive.", &
279 "skip preconditioning"), &
282 CALL section_add_keyword(section, keyword)
283 CALL keyword_release(keyword)
284
285 CALL keyword_create(keyword, __location__, name="ENERGY_GAP", &
286 description="Energy gap estimate [a.u.] for preconditioning", &
287 usage="ENERGY_GAP 0.1", &
288 default_r_val=def_egap)
289 CALL section_add_keyword(section, keyword)
290 CALL keyword_release(keyword)
291
292 CALL keyword_create(keyword, __location__, name="EVERY_N_STEP", &
293 description="Perform a linear response calculation every N-th step for MD run", &
294 usage="EVERY_N_STEP 50", default_i_val=1)
295 CALL section_add_keyword(section, keyword)
296 CALL keyword_release(keyword)
297
298 CALL keyword_create(keyword, __location__, name="RESTART", &
299 description="Restart the response calculation if the restart file exists", &
300 usage="RESTART", &
301 default_l_val=.false., lone_keyword_l_val=.true.)
302 CALL section_add_keyword(section, keyword)
303 CALL keyword_release(keyword)
304
305 CALL keyword_create(keyword, __location__, name="WFN_RESTART_FILE_NAME", &
306 variants=(/"RESTART_FILE_NAME"/), &
307 description="Root of the file names where to read the response functions from "// &
308 "which to restart the calculation of the linear response", &
309 usage="WFN_RESTART_FILE_NAME <FILENAME>", &
310 type_of_var=lchar_t)
311 CALL section_add_keyword(section, keyword)
312 CALL keyword_release(keyword)
313
314 IF (create_subsections) THEN
315 NULLIFY (subsection)
316
317 CALL create_localize_section(subsection)
318 CALL section_add_subsection(section, subsection)
319 CALL section_release(subsection)
320
321 CALL create_current_section(subsection)
322 CALL section_add_subsection(section, subsection)
323 CALL section_release(subsection)
324
325 CALL create_nmr_section(subsection)
326 CALL section_add_subsection(section, subsection)
327 CALL section_release(subsection)
328
329 CALL create_spin_spin_section(subsection)
330 CALL section_add_subsection(section, subsection)
331 CALL section_release(subsection)
332
333 CALL create_epr_section(subsection)
334 CALL section_add_subsection(section, subsection)
335 CALL section_release(subsection)
336
337 CALL create_polarizability_section(subsection)
338 CALL section_add_subsection(section, subsection)
339 CALL section_release(subsection)
340
341 CALL create_dcdr_section(subsection)
342 CALL section_add_subsection(section, subsection)
343 CALL section_release(subsection)
344
345 CALL create_vcd_section(subsection)
346 CALL section_add_subsection(section, subsection)
347 CALL section_release(subsection)
348
349 CALL section_create(subsection, __location__, name="PRINT", &
350 description="printing of information during the linear response calculation", &
351 repeats=.false.)
352
354 print_key, __location__, "program_run_info", &
355 description="Controls the printing of basic iteration information during the LINRES calculation", &
356 print_level=low_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
357 CALL section_add_subsection(subsection, print_key)
358 CALL section_release(print_key)
359
360 CALL cp_print_key_section_create(print_key, __location__, "RESTART", &
361 description="Controls the dumping of restart file of the response wavefunction. "// &
362 "For each set of response functions, i.e. for each perturbation, "// &
363 "one different restart file is dumped. These restart files should be "// &
364 "employed only to restart the same type of LINRES calculation, "// &
365 "i.e. with the same perturbation.", &
366 print_level=low_print_level, common_iter_levels=3, each_iter_names=s2a("ITER"), &
367 add_last=add_last_numeric, each_iter_values=(/3/), filename="")
368 CALL section_add_subsection(subsection, print_key)
369 CALL section_release(print_key)
370
371 CALL section_add_subsection(section, subsection)
372 CALL section_release(subsection)
373
374 END IF
375
376 END SUBROUTINE create_linres_section
377
378! **************************************************************************************************
379!> \brief creates the input structure used to activate
380!> calculation of position perturbation DFPT
381!> \param section ...
382!> \author Sandra Luber, Edward Ditler
383! **************************************************************************************************
384 SUBROUTINE create_dcdr_section(section)
385
386 TYPE(section_type), POINTER :: section
387
388 LOGICAL :: failure
389 TYPE(keyword_type), POINTER :: keyword
390 TYPE(section_type), POINTER :: print_key, subsection
391
392 failure = .false.
393 NULLIFY (keyword, print_key, subsection)
394
395 cpassert(.NOT. ASSOCIATED(section))
396
397 IF (.NOT. failure) THEN
398 CALL section_create(section, __location__, name="DCDR", &
399 description="Compute analytical gradients the dipole moments.", &
400 n_keywords=50, n_subsections=1, repeats=.false.)
401
402 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
403 description="controls the activation of the APT calculation", &
404 usage="&DCDR T", &
405 default_l_val=.false., &
406 lone_keyword_l_val=.true.)
407 CALL section_add_keyword(section, keyword)
408 CALL keyword_release(keyword)
409
410 CALL keyword_create(keyword, __location__, name="LIST_OF_ATOMS", &
411 description="Specifies a list of atoms.", &
412 usage="LIST_OF_ATOMS {integer} {integer} .. {integer}", repeats=.true., &
413 n_var=-1, type_of_var=integer_t)
414 CALL section_add_keyword(section, keyword)
415 CALL keyword_release(keyword)
416
417 CALL keyword_create(keyword, __location__, name="DISTRIBUTED_ORIGIN", &
418 variants=(/"DO_GAUGE"/), &
419 description="Use the distributed origin (DO) gauge?", &
420 usage="DISTRIBUTED_ORIGIN T", &
421 default_l_val=.false., lone_keyword_l_val=.true.)
422 CALL section_add_keyword(section, keyword)
423 CALL keyword_release(keyword)
424
425 CALL keyword_create(keyword, __location__, name="ORBITAL_CENTER", &
426 description="The orbital center.", &
427 usage="ORBITAL_CENTER WANNIER", &
428 default_i_val=current_orb_center_wannier, &
429 enum_c_vals=s2a("WANNIER", "COMMON", "ATOM", "BOX"), &
430 enum_desc=s2a("Use the Wannier centers.", &
431 "Use a common center (works only for an isolate molecule).", &
432 "Use the atoms as center.", &
433 "Boxing."), &
436 CALL section_add_keyword(section, keyword)
437 CALL keyword_release(keyword)
438
439 CALL keyword_create(keyword, __location__, name="REFERENCE", &
440 description="Gauge origin of the velocity gauge factor.", &
441 enum_c_vals=s2a("COM", "COAC", "USER_DEFINED", "ZERO"), &
442 enum_desc=s2a("Use Center of Mass", &
443 "Use Center of Atomic Charges", &
444 "Use User-defined Point", &
445 "Use Origin of Coordinate System"), &
446 enum_i_vals=(/use_mom_ref_com, &
450 default_i_val=use_mom_ref_zero)
451 CALL section_add_keyword(section, keyword)
452 CALL keyword_release(keyword)
453
454 CALL keyword_create(keyword, __location__, name="REFERENCE_POINT", &
455 description="User-defined reference point of the velocity gauge factor.", &
456 usage="REFERENCE_POINT x y z", &
457 repeats=.false., n_var=3, type_of_var=real_t, unit_str='bohr')
458 CALL section_add_keyword(section, keyword)
459 CALL keyword_release(keyword)
460
461 CALL keyword_create(keyword, __location__, name="Z_MATRIX_METHOD", &
462 description="Use Z_matrix method to solve the response equation", &
463 usage="Z_MATRIX_METHOD T", &
464 default_l_val=.false., lone_keyword_l_val=.true.)
465 CALL section_add_keyword(section, keyword)
466 CALL keyword_release(keyword)
467
468 CALL keyword_create(keyword, __location__, name="APT_FD", &
469 description="Use numerical differentiation to compute the APT, "// &
470 "switches off the calculation of dcdr analytical derivatives. "// &
471 "Requires RUN_TYPE = ENERGY_FORCE or MD.", &
472 usage="APT_FD T", &
473 default_l_val=.false., lone_keyword_l_val=.true.)
474 CALL section_add_keyword(section, keyword)
475 CALL keyword_release(keyword)
476
477 CALL keyword_create(keyword, __location__, name="APT_FD_DE", &
478 description="Electric field strength (atomic units) to use for finite differences", &
479 repeats=.false., &
480 n_var=1, &
481 type_of_var=real_t, &
482 default_r_val=0.0003_dp, &
483 usage="APT_FD_DE 1.0E-4")
484 CALL section_add_keyword(section, keyword)
485 CALL keyword_release(keyword)
486
487 CALL keyword_create(keyword, __location__, name="APT_FD_METHOD", &
488 description="Numerical differentiation method", &
489 usage="APT_FD_METHOD FD", &
490 default_i_val=1, &
491 !enum_c_vals=s2a("FD", "2PNT"), &
492 enum_c_vals=s2a("2PNT"), &
493 !enum_desc=s2a("Forward differences.", &
494 ! "Symmetric two-point differences."), &
495 enum_desc=s2a("Symmetric two-point differences."), &
496 !enum_i_vals=(/0, 1/))
497 enum_i_vals=(/1/))
498 CALL section_add_keyword(section, keyword)
499 CALL keyword_release(keyword)
500
501 NULLIFY (subsection)
502 CALL section_create(subsection, __location__, name="PRINT", &
503 description="print results of the magnetic dipole moment calculation", &
504 repeats=.false.)
505
506 CALL cp_print_key_section_create(print_key, __location__, "APT", &
507 description="Controls the printing of the electric dipole gradient", &
508 print_level=low_print_level, add_last=add_last_numeric, filename="")
509 CALL section_add_subsection(subsection, print_key)
510 CALL section_release(print_key)
511
512 CALL section_add_subsection(section, subsection)
513 CALL section_release(subsection)
514
515 NULLIFY (subsection)
516 CALL create_interp_section(subsection)
517 CALL section_add_subsection(section, subsection)
518 CALL section_release(subsection)
519
520 END IF
521
522 END SUBROUTINE create_dcdr_section
523
524! **************************************************************************************************
525!> \brief creates the input structure used to activate
526!> calculation of VCD spectra using DFPT
527!> \param section ...
528!> \author Sandra Luber, Tomas Zimmermann, Edward Ditler
529! **************************************************************************************************
530 SUBROUTINE create_vcd_section(section)
531
532 TYPE(section_type), POINTER :: section
533
534 TYPE(keyword_type), POINTER :: keyword
535 TYPE(section_type), POINTER :: print_key, subsection
536
537 NULLIFY (keyword, print_key, subsection)
538
539 cpassert(.NOT. ASSOCIATED(section))
540
541 CALL section_create(section, __location__, name="VCD", &
542 description="Carry out a VCD calculation.", &
543 n_keywords=50, n_subsections=1, repeats=.false.)
544
545 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
546 description="controls the activation of the APT/AAT calculation", &
547 usage="&VCD T", &
548 default_l_val=.false., &
549 lone_keyword_l_val=.true.)
550 CALL section_add_keyword(section, keyword)
551 CALL keyword_release(keyword)
552
553 CALL keyword_create(keyword, __location__, name="LIST_OF_ATOMS", &
554 description="Specifies a list of atoms.", &
555 usage="LIST_OF_ATOMS {integer} {integer} .. {integer}", repeats=.true., &
556 n_var=-1, type_of_var=integer_t)
557 CALL section_add_keyword(section, keyword)
558 CALL keyword_release(keyword)
559
560 CALL keyword_create(keyword, __location__, name="DISTRIBUTED_ORIGIN", &
561 variants=(/"DO_GAUGE"/), &
562 description="Use the distributed origin (DO) gauge?", &
563 usage="DISTRIBUTED_ORIGIN T", &
564 default_l_val=.false., lone_keyword_l_val=.true.)
565 CALL section_add_keyword(section, keyword)
566 CALL keyword_release(keyword)
567
568 CALL keyword_create(keyword, __location__, name="ORIGIN_DEPENDENT_MFP", &
569 description="Use the origin dependent MFP operator.", &
570 usage="ORIGIN_DEPENDENT_MFP T", &
571 default_l_val=.false., lone_keyword_l_val=.true.)
572 CALL section_add_keyword(section, keyword)
573 CALL keyword_release(keyword)
574
575 CALL keyword_create(keyword, __location__, name="ORBITAL_CENTER", &
576 description="The orbital center.", &
577 usage="ORBITAL_CENTER WANNIER", &
578 default_i_val=current_orb_center_wannier, &
579 enum_c_vals=s2a("WANNIER", "COMMON", "ATOM", "BOX"), &
580 enum_desc=s2a("Use the Wannier centers.", &
581 "Use a common center (works only for an isolate molecule).", &
582 "Use the atoms as center.", &
583 "Boxing."), &
586 CALL section_add_keyword(section, keyword)
587 CALL keyword_release(keyword)
588
589 ! The origin of the magnetic dipole operator (r - MAGNETIC_ORIGIN) x momentum
590 CALL keyword_create(keyword, __location__, name="MAGNETIC_ORIGIN", &
591 description="Gauge origin of the magnetic dipole operator.", &
592 enum_c_vals=s2a("COM", "COAC", "USER_DEFINED", "ZERO"), &
593 enum_desc=s2a("Use Center of Mass", &
594 "Use Center of Atomic Charges", &
595 "Use User-defined Point", &
596 "Use Origin of Coordinate System"), &
597 enum_i_vals=(/use_mom_ref_com, &
601 default_i_val=use_mom_ref_zero)
602 CALL section_add_keyword(section, keyword)
603 CALL keyword_release(keyword)
604
605 CALL keyword_create(keyword, __location__, name="MAGNETIC_ORIGIN_REFERENCE", &
606 description="User-defined reference point of the magnetic dipole operator.", &
607 usage="MAGNETIC_ORIGIN_REFERENCE x y z", &
608 repeats=.false., n_var=3, type_of_var=real_t, unit_str='bohr')
609 CALL section_add_keyword(section, keyword)
610 CALL keyword_release(keyword)
611
612 ! The origin of the coordinate system
613 CALL keyword_create(keyword, __location__, name="SPATIAL_ORIGIN", &
614 description="Gauge origin of the velocity gauge factor/spatial origin.", &
615 enum_c_vals=s2a("COM", "COAC", "USER_DEFINED", "ZERO"), &
616 enum_desc=s2a("Use Center of Mass", &
617 "Use Center of Atomic Charges", &
618 "Use User-defined Point", &
619 "Use Origin of Coordinate System"), &
620 enum_i_vals=(/use_mom_ref_com, &
624 default_i_val=use_mom_ref_zero)
625 CALL section_add_keyword(section, keyword)
626 CALL keyword_release(keyword)
627
628 CALL keyword_create(keyword, __location__, name="SPATIAL_ORIGIN_REFERENCE", &
629 description="User-defined reference point of the velocity gauge factor/spatial origin.", &
630 usage="SPATIAL_ORIGIN_REFERENCE x y z", &
631 repeats=.false., n_var=3, type_of_var=real_t, unit_str='bohr')
632 CALL section_add_keyword(section, keyword)
633 CALL keyword_release(keyword)
634
635 NULLIFY (subsection)
636 CALL section_create(subsection, __location__, name="PRINT", &
637 description="print results of the magnetic dipole moment calculation", &
638 repeats=.false.)
639
640 CALL cp_print_key_section_create(print_key, __location__, "VCD", &
641 description="Controls the printing of the APTs and AATs", &
642 print_level=low_print_level, add_last=add_last_numeric, filename="")
643 CALL section_add_subsection(subsection, print_key)
644 CALL section_release(print_key)
645
646 CALL section_add_subsection(section, subsection)
647 CALL section_release(subsection)
648
649 NULLIFY (subsection)
650 CALL create_interp_section(subsection)
651 CALL section_add_subsection(section, subsection)
652 CALL section_release(subsection)
653
654 END SUBROUTINE create_vcd_section
655
656! **************************************************************************************************
657!> \brief creates the input structure used to activate
658!> calculation of induced current DFPT
659!> Available properties : none
660!> \param section the section to create
661!> \author MI/VW
662! **************************************************************************************************
663 SUBROUTINE create_current_section(section)
664 TYPE(section_type), POINTER :: section
665
666 TYPE(keyword_type), POINTER :: keyword
667 TYPE(section_type), POINTER :: print_key, subsection
668
669 NULLIFY (keyword, print_key, subsection)
670
671 cpassert(.NOT. ASSOCIATED(section))
672 CALL section_create(section, __location__, name="current", &
673 description="The induced current density is calculated by DFPT.", &
674 n_keywords=4, n_subsections=1, repeats=.false., &
675 citations=(/sebastiani2001, weber2009/))
676
677 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
678 description="controls the activation of the induced current calculation", &
679 usage="&CURRENT T", &
680 default_l_val=.false., &
681 lone_keyword_l_val=.true.)
682 CALL section_add_keyword(section, keyword)
683 CALL keyword_release(keyword)
684
685 CALL keyword_create(keyword, __location__, name="GAUGE", &
686 description="The gauge used to compute the induced current within GAPW.", &
687 usage="GAUGE R", &
688 default_i_val=current_gauge_r_and_step_func, &
689 enum_c_vals=s2a("R", "R_AND_STEP_FUNCTION", "ATOM"), &
690 enum_desc=s2a("Position gauge (doesnt work well).", &
691 "Position and step function for the soft and the local parts, respectively.", &
692 "Atoms."), &
694 CALL section_add_keyword(section, keyword)
695 CALL keyword_release(keyword)
696
697 CALL keyword_create(keyword, __location__, name="GAUGE_ATOM_RADIUS", &
698 description="Build the gauge=atom using only the atoms within this radius.", &
699 usage="GAUGE_ATOM_RADIUS 10.0", &
700 type_of_var=real_t, &
701 default_r_val=cp_unit_to_cp2k(value=4.0_dp, unit_str="angstrom"), &
702 unit_str="angstrom")
703 CALL section_add_keyword(section, keyword)
704 CALL keyword_release(keyword)
705
706 CALL keyword_create(keyword, __location__, name="USE_OLD_GAUGE_ATOM", &
707 description="Use the old way to compute the gauge.", &
708 usage="USE_OLD_GAUGE_ATOM T", &
709 default_l_val=.true., lone_keyword_l_val=.true.)
710 CALL section_add_keyword(section, keyword)
711 CALL keyword_release(keyword)
712
713 CALL keyword_create(keyword, __location__, name="ORBITAL_CENTER", &
714 description="The orbital center.", &
715 usage="ORBITAL_CENTER WANNIER", &
716 default_i_val=current_orb_center_wannier, &
717 enum_c_vals=s2a("WANNIER", "COMMON", "ATOM", "BOX"), &
718 enum_desc=s2a("Use the Wannier centers.", &
719 "Use a common center (works only for an isolate molecule).", &
720 "Use the atoms as center.", &
721 "Boxing."), &
724 CALL section_add_keyword(section, keyword)
725 CALL keyword_release(keyword)
726
727 CALL keyword_create(keyword, __location__, name="COMMON_CENTER", &
728 description="The common center ", usage="COMMON_CENTER 0.0 1.0 0.0", &
729 n_var=3, default_r_vals=(/0.0_dp, 0.0_dp, 0.0_dp/), type_of_var=real_t, &
730 unit_str="angstrom")
731 CALL section_add_keyword(section, keyword)
732 CALL keyword_release(keyword)
733
734 CALL keyword_create(keyword, __location__, name="NBOX", &
735 description="How many boxes along each directions ", usage="NBOX 6 6 5", &
736 n_var=3, default_i_vals=(/4, 4, 4/), type_of_var=integer_t)
737 CALL section_add_keyword(section, keyword)
738 CALL keyword_release(keyword)
739
740 CALL keyword_create(keyword, __location__, name="CHI_PBC", &
741 description="Calculate the succeptibility correction to the shift with PBC", &
742 usage="CHI_PBC T", &
743 default_l_val=.false., lone_keyword_l_val=.true.)
744 CALL section_add_keyword(section, keyword)
745 CALL keyword_release(keyword)
746
747 CALL keyword_create(keyword, __location__, name="FORCE_NO_FULL", &
748 description="Avoid the calculation of the state dependent perturbation term, "// &
749 "even if the orbital centers are set at Wannier centers or at Atom centers", &
750 usage="FORCE_NO_FULL T", &
751 default_l_val=.false., 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="SELECTED_STATES_ON_ATOM_LIST", &
756 description="Indexes of the atoms for selecting"// &
757 " the states to be used for the response calculations.", &
758 usage="SELECTED_STATES_ON_ATOM_LIST 1 2 10", &
759 n_var=-1, type_of_var=integer_t, repeats=.true.)
760 CALL section_add_keyword(section, keyword)
761 CALL keyword_release(keyword)
762
763 CALL keyword_create(keyword, __location__, name="SELECTED_STATES_ATOM_RADIUS", &
764 description="Select all the states included in the given radius around each atoms "// &
765 "in SELECTED_STATES_ON_ATOM_LIST.", &
766 usage="SELECTED_STATES_ATOM_RADIUS 2.0", &
767 type_of_var=real_t, &
768 default_r_val=cp_unit_to_cp2k(value=4.0_dp, unit_str="angstrom"), &
769 unit_str="angstrom")
770 CALL section_add_keyword(section, keyword)
771 CALL keyword_release(keyword)
772
773 CALL keyword_create(keyword, __location__, name="RESTART_CURRENT", &
774 description="Restart the induced current density calculation"// &
775 " from a previous run (not working yet).", &
776 usage="RESTART_CURRENT", default_l_val=.false., &
777 lone_keyword_l_val=.true.)
778 CALL section_add_keyword(section, keyword)
779 CALL keyword_release(keyword)
780
781 NULLIFY (subsection)
782 CALL section_create(subsection, __location__, name="PRINT", &
783 description="print results of induced current density calculation", &
784 repeats=.false.)
785
786 CALL cp_print_key_section_create(print_key, __location__, "CURRENT_CUBES", &
787 description="Controls the printing of the induced current density (not working yet).", &
788 print_level=high_print_level, add_last=add_last_numeric, filename="")
789 CALL keyword_create(keyword, __location__, name="stride", &
790 description="The stride (X,Y,Z) used to write the cube file "// &
791 "(larger values result in smaller cube files). You can provide 3 numbers (for X,Y,Z) or"// &
792 " 1 number valid for all components (not working yet).", &
793 usage="STRIDE 2 2 2", n_var=-1, default_i_vals=(/2, 2, 2/), type_of_var=integer_t)
794 CALL section_add_keyword(print_key, keyword)
795 CALL keyword_release(keyword)
796 CALL keyword_create(keyword, __location__, name="APPEND", &
797 description="append the cube files when they already exist", &
798 default_l_val=.false., lone_keyword_l_val=.true.)
799 CALL section_add_keyword(print_key, keyword)
800 CALL keyword_release(keyword)
801
802 CALL section_add_subsection(subsection, print_key)
803 CALL section_release(print_key)
804
805 CALL cp_print_key_section_create(print_key, __location__, "RESPONSE_FUNCTION_CUBES", &
806 description="Controls the printing of the response functions (not working yet).", &
807 print_level=high_print_level, add_last=add_last_numeric, filename="")
808 CALL keyword_create(keyword, __location__, name="stride", &
809 description="The stride (X,Y,Z) used to write the cube file "// &
810 "(larger values result in smaller cube files). You can provide 3 numbers (for X,Y,Z) or"// &
811 " 1 number valid for all components (not working yet).", &
812 usage="STRIDE 2 2 2", n_var=-1, default_i_vals=(/2, 2, 2/), type_of_var=integer_t)
813 CALL section_add_keyword(print_key, keyword)
814 CALL keyword_release(keyword)
815
816 CALL keyword_create(keyword, __location__, name="CUBES_LU_BOUNDS", &
817 variants=(/"CUBES_LU"/), &
818 description="The lower and upper index of the states to be printed as cube (not working yet).", &
819 usage="CUBES_LU_BOUNDS integer integer", &
820 n_var=2, default_i_vals=(/0, -2/), type_of_var=integer_t)
821 CALL section_add_keyword(print_key, keyword)
822 CALL keyword_release(keyword)
823
824 CALL keyword_create(keyword, __location__, name="CUBES_LIST", &
825 description="Indexes of the states to be printed as cube files "// &
826 "This keyword can be repeated several times "// &
827 "(useful if you have to specify many indexes) (not working yet).", &
828 usage="CUBES_LIST 1 2", &
829 n_var=-1, type_of_var=integer_t, repeats=.true.)
830 CALL section_add_keyword(print_key, keyword)
831 CALL keyword_release(keyword)
832 CALL keyword_create(keyword, __location__, name="APPEND", &
833 description="append the cube files when they already exist", &
834 default_l_val=.false., lone_keyword_l_val=.true.)
835 CALL section_add_keyword(print_key, keyword)
836 CALL keyword_release(keyword)
837
838 CALL section_add_subsection(subsection, print_key)
839 CALL section_release(print_key)
840
841 CALL section_add_subsection(section, subsection)
842 CALL section_release(subsection)
843
844 NULLIFY (subsection)
845 CALL create_interp_section(subsection)
846 CALL section_add_subsection(section, subsection)
847 CALL section_release(subsection)
848
849 END SUBROUTINE create_current_section
850
851! **************************************************************************************************
852!> \brief creates the input structure used to activate
853!> calculation of NMR chemical shift using
854!> the induced current obtained from DFPT
855!> Available properties : none
856!> \param section the section to create
857!> \author MI/VW
858! **************************************************************************************************
859 SUBROUTINE create_nmr_section(section)
860 TYPE(section_type), POINTER :: section
861
862 TYPE(keyword_type), POINTER :: keyword
863 TYPE(section_type), POINTER :: print_key, subsection
864
865 NULLIFY (keyword, print_key, subsection)
866
867 cpassert(.NOT. ASSOCIATED(section))
868 CALL section_create(section, __location__, name="nmr", &
869 description="The chemical shift is calculated by DFPT.", &
870 n_keywords=5, n_subsections=1, repeats=.false., &
871 citations=(/weber2009/))
872
873 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
874 description="controls the activation of the nmr calculation", &
875 usage="&NMR T", &
876 default_l_val=.false., &
877 lone_keyword_l_val=.true.)
878 CALL section_add_keyword(section, keyword)
879 CALL keyword_release(keyword)
880
881 CALL keyword_create(keyword, __location__, name="INTERPOLATE_SHIFT", &
882 description="Calculate the soft part of the chemical shift by interpolation ", &
883 usage="INTERPOLATE_SHIFT T", &
884 default_l_val=.false., lone_keyword_l_val=.true.)
885 CALL section_add_keyword(section, keyword)
886 CALL keyword_release(keyword)
887
888 CALL keyword_create(keyword, __location__, name="NICS", &
889 description="Calculate the chemical shift in a set of points"// &
890 " given from an external file", usage="NICS", &
891 default_l_val=.false., lone_keyword_l_val=.true.)
892 CALL section_add_keyword(section, keyword)
893 CALL keyword_release(keyword)
894
895 CALL keyword_create(keyword, __location__, name="NICS_FILE_NAME", &
896 description="Name of the file with the NICS points coordinates", &
897 usage="NICS_FILE_NAME nics_file", &
898 default_lc_val="nics_file")
899 CALL section_add_keyword(section, keyword)
900 CALL keyword_release(keyword)
901
902 CALL keyword_create(keyword, __location__, name="RESTART_NMR", &
903 description="Restart the NMR calculation from a previous run (NOT WORKING YET)", &
904 usage="RESTART_NMR", default_l_val=.false., &
905 lone_keyword_l_val=.true.)
906 CALL section_add_keyword(section, keyword)
907 CALL keyword_release(keyword)
908
909 CALL keyword_create(keyword, __location__, name="SHIFT_GAPW_RADIUS", &
910 description="While computing the local part of the shift (GAPW), "// &
911 "the integration is restricted to nuclei that are within this radius.", &
912 usage="SHIFT_GAPW_RADIUS 20.0", &
913 type_of_var=real_t, &
914 default_r_val=cp_unit_to_cp2k(value=60.0_dp, unit_str="angstrom"), &
915 unit_str="angstrom")
916 CALL section_add_keyword(section, keyword)
917 CALL keyword_release(keyword)
918
919 NULLIFY (subsection)
920 CALL section_create(subsection, __location__, name="PRINT", &
921 description="print results of nmr calculation", &
922 repeats=.false.)
923
924 CALL cp_print_key_section_create(print_key, __location__, "RESPONSE_FUNCTION_CUBES", &
925 description="Controls the printing of the response functions ", &
926 print_level=high_print_level, add_last=add_last_numeric, filename="")
927 CALL keyword_create(keyword, __location__, name="stride", &
928 description="The stride (X,Y,Z) used to write the cube file "// &
929 "(larger values result in smaller cube files). You can provide 3 numbers (for X,Y,Z) or"// &
930 " 1 number valid for all components.", &
931 usage="STRIDE 2 2 2", n_var=-1, default_i_vals=(/2, 2, 2/), type_of_var=integer_t)
932 CALL section_add_keyword(print_key, keyword)
933 CALL keyword_release(keyword)
934
935 CALL keyword_create(keyword, __location__, name="CUBES_LU_BOUNDS", &
936 variants=(/"CUBES_LU"/), &
937 description="The lower and upper index of the states to be printed as cube", &
938 usage="CUBES_LU_BOUNDS integer integer", &
939 n_var=2, default_i_vals=(/0, -2/), type_of_var=integer_t)
940 CALL section_add_keyword(print_key, keyword)
941 CALL keyword_release(keyword)
942
943 CALL keyword_create(keyword, __location__, name="CUBES_LIST", &
944 description="Indexes of the states to be printed as cube files "// &
945 "This keyword can be repeated several times "// &
946 "(useful if you have to specify many indexes).", &
947 usage="CUBES_LIST 1 2", &
948 n_var=-1, type_of_var=integer_t, repeats=.true.)
949 CALL section_add_keyword(print_key, keyword)
950 CALL keyword_release(keyword)
951 CALL keyword_create(keyword, __location__, name="APPEND", &
952 description="append the cube files when they already exist", &
953 default_l_val=.false., lone_keyword_l_val=.true.)
954 CALL section_add_keyword(print_key, keyword)
955 CALL keyword_release(keyword)
956
957 CALL section_add_subsection(subsection, print_key)
958 CALL section_release(print_key)
959
960 CALL cp_print_key_section_create(print_key, __location__, "CHI_TENSOR", &
961 description="Controls the printing of susceptibility", &
962 print_level=high_print_level, add_last=add_last_numeric, filename="")
963 CALL section_add_subsection(subsection, print_key)
964 CALL section_release(print_key)
965
966 CALL cp_print_key_section_create(print_key, __location__, "SHIELDING_TENSOR", &
967 description="Controls the printing of the chemical shift", &
968 print_level=low_print_level, add_last=add_last_numeric, filename="")
969
970 CALL keyword_create(keyword, __location__, name="ATOMS_LU_BOUNDS", &
971 variants=(/"ATOMS_LU"/), &
972 description="The lower and upper atomic index for which the tensor is printed", &
973 usage="ATOMS_LU_BOUNDS integer integer", &
974 n_var=2, default_i_vals=(/0, -2/), type_of_var=integer_t)
975 CALL section_add_keyword(print_key, keyword)
976 CALL keyword_release(keyword)
977
978 CALL keyword_create(keyword, __location__, name="ATOMS_LIST", &
979 description="list of atoms for which the shift is printed into a file ", &
980 usage="ATOMS_LIST 1 2", n_var=-1, &
981 type_of_var=integer_t, repeats=.true.)
982 CALL section_add_keyword(print_key, keyword)
983 CALL keyword_release(keyword)
984
985 CALL section_add_subsection(subsection, print_key)
986 CALL section_release(print_key)
987
988 CALL section_add_subsection(section, subsection)
989 CALL section_release(subsection)
990
991 NULLIFY (subsection)
992 CALL create_interp_section(subsection)
993 CALL section_add_subsection(section, subsection)
994 CALL section_release(subsection)
995
996 END SUBROUTINE create_nmr_section
997
998! **************************************************************************************************
999!> \brief creates the input structure used to activate
1000!> calculation of NMR spin-spin coupling (implementation not operating)
1001!> Available properties : none
1002!> \param section the section to create
1003!> \author VW
1004! **************************************************************************************************
1005 SUBROUTINE create_spin_spin_section(section)
1006 TYPE(section_type), POINTER :: section
1007
1008 TYPE(keyword_type), POINTER :: keyword
1009 TYPE(section_type), POINTER :: print_key, subsection
1010
1011 NULLIFY (keyword, print_key, subsection)
1012
1013 cpassert(.NOT. ASSOCIATED(section))
1014 CALL section_create(section, __location__, name="spinspin", &
1015 description="Compute indirect spin-spin coupling constants.", &
1016 n_keywords=5, n_subsections=1, repeats=.false.)
1017
1018 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
1019 description="controls the activation of the nmr calculation", &
1020 usage="&SPINSPIN T", &
1021 default_l_val=.false., &
1022 lone_keyword_l_val=.true.)
1023 CALL section_add_keyword(section, keyword)
1024 CALL keyword_release(keyword)
1025
1026 CALL keyword_create(keyword, __location__, name="RESTART_SPINSPIN", &
1027 description="Restart the spin-spin calculation from a previous run (NOT WORKING YET)", &
1028 usage="RESTART_SPINSPIN", default_l_val=.false., &
1029 lone_keyword_l_val=.true.)
1030 CALL section_add_keyword(section, keyword)
1031 CALL keyword_release(keyword)
1032
1033 CALL keyword_create(keyword, __location__, name="ISSC_ON_ATOM_LIST", &
1034 description="Atoms for which the issc is computed.", &
1035 usage="ISSC_ON_ATOM_LIST 1 2 10", &
1036 n_var=-1, type_of_var=integer_t, repeats=.true.)
1037 CALL section_add_keyword(section, keyword)
1038 CALL keyword_release(keyword)
1039
1040 CALL keyword_create(keyword, __location__, name="DO_FC", &
1041 description="Compute the Fermi contact contribution", &
1042 usage="DO_FC F", &
1043 default_l_val=.true., lone_keyword_l_val=.true.)
1044 CALL section_add_keyword(section, keyword)
1045 CALL keyword_release(keyword)
1046
1047 CALL keyword_create(keyword, __location__, name="DO_SD", &
1048 description="Compute the spin-dipolar contribution", &
1049 usage="DO_SD F", &
1050 default_l_val=.true., lone_keyword_l_val=.true.)
1051 CALL section_add_keyword(section, keyword)
1052 CALL keyword_release(keyword)
1053
1054 CALL keyword_create(keyword, __location__, name="DO_PSO", &
1055 description="Compute the paramagnetic spin-orbit contribution", &
1056 usage="DO_PSO F", &
1057 default_l_val=.true., lone_keyword_l_val=.true.)
1058 CALL section_add_keyword(section, keyword)
1059 CALL keyword_release(keyword)
1060
1061 CALL keyword_create(keyword, __location__, name="DO_DSO", &
1062 description="Compute the diamagnetic spin-orbit contribution (NOT YET IMPLEMENTED)", &
1063 usage="DO_DSO F", &
1064 default_l_val=.true., lone_keyword_l_val=.true.)
1065 CALL section_add_keyword(section, keyword)
1066 CALL keyword_release(keyword)
1067
1068 NULLIFY (subsection)
1069 CALL section_create(subsection, __location__, name="PRINT", &
1070 description="print results of the indirect spin-spin calculation", &
1071 repeats=.false.)
1072
1073 CALL cp_print_key_section_create(print_key, __location__, "K_MATRIX", &
1074 description="Controls the printing of the indirect spin-spin matrix", &
1075 print_level=low_print_level, add_last=add_last_numeric, filename="")
1076
1077 CALL keyword_create(keyword, __location__, name="ATOMS_LIST", &
1078 description="list of atoms for which the indirect spin-spin is printed into a file ", &
1079 usage="ATOMS_LIST 1 2", n_var=-1, &
1080 type_of_var=integer_t, repeats=.true.)
1081 CALL section_add_keyword(print_key, keyword)
1082 CALL keyword_release(keyword)
1083
1084 CALL section_add_subsection(subsection, print_key)
1085 CALL section_release(print_key)
1086
1087 CALL section_add_subsection(section, subsection)
1088 CALL section_release(subsection)
1089
1090 NULLIFY (subsection)
1091 CALL create_interp_section(subsection)
1092 CALL section_add_subsection(section, subsection)
1093 CALL section_release(subsection)
1094
1095 END SUBROUTINE create_spin_spin_section
1096
1097! **************************************************************************************************
1098!> \brief creates the input structure used to activate
1099!> calculation of EPR using
1100!> the induced current obtained from DFPT
1101!> Available properties : none
1102!> \param section the section to create
1103!> \author VW
1104! **************************************************************************************************
1105 SUBROUTINE create_epr_section(section)
1106 TYPE(section_type), POINTER :: section
1107
1108 TYPE(keyword_type), POINTER :: keyword
1109 TYPE(section_type), POINTER :: print_key, subsection, subsubsection
1110
1111 NULLIFY (keyword, print_key, subsection, subsubsection)
1112
1113 cpassert(.NOT. ASSOCIATED(section))
1114 CALL section_create(section, __location__, name="EPR", &
1115 description="The g tensor is calculated by DFPT ", &
1116 n_keywords=5, n_subsections=1, repeats=.false., &
1117 citations=(/weber2009/))
1118
1119 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
1120 description="controls the activation of the epr calculation", &
1121 usage="&EPR T", &
1122 default_l_val=.false., &
1123 lone_keyword_l_val=.true.)
1124 CALL section_add_keyword(section, keyword)
1125 CALL keyword_release(keyword)
1126
1127 CALL keyword_create(keyword, __location__, name="RESTART_EPR", &
1128 description="Restart the EPR calculation from a previous run (NOT WORKING)", &
1129 usage="RESTART_EPR", default_l_val=.false., &
1130 lone_keyword_l_val=.true.)
1131 CALL section_add_keyword(section, keyword)
1132 CALL keyword_release(keyword)
1133
1134 NULLIFY (subsection)
1135 CALL section_create(subsection, __location__, name="PRINT", &
1136 description="print results of epr calculation", &
1137 repeats=.false.)
1138
1139 CALL cp_print_key_section_create(print_key, __location__, "NABLAVKS_CUBES", &
1140 description="Controls the printing of the components of nabla v_ks ", &
1141 print_level=high_print_level, add_last=add_last_numeric, filename="")
1142 CALL keyword_create(keyword, __location__, name="stride", &
1143 description="The stride (X,Y,Z) used to write the cube file "// &
1144 "(larger values result in smaller cube files). You can provide 3 numbers (for X,Y,Z) or"// &
1145 " 1 number valid for all components.", &
1146 usage="STRIDE 2 2 2", n_var=-1, default_i_vals=(/2, 2, 2/), type_of_var=integer_t)
1147 CALL section_add_keyword(print_key, keyword)
1148 CALL keyword_release(keyword)
1149 CALL keyword_create(keyword, __location__, name="APPEND", &
1150 description="append the cube files when they already exist", &
1151 default_l_val=.false., lone_keyword_l_val=.true.)
1152 CALL section_add_keyword(print_key, keyword)
1153 CALL keyword_release(keyword)
1154
1155 CALL section_add_subsection(subsection, print_key)
1156 CALL section_release(print_key)
1157
1158 CALL cp_print_key_section_create(print_key, __location__, "G_TENSOR", &
1159 description="Controls the printing of the g tensor", &
1160 print_level=high_print_level, add_last=add_last_numeric, filename="")
1161 CALL create_xc_section(subsubsection)
1162 CALL section_add_subsection(print_key, subsubsection)
1163 CALL section_release(subsubsection)
1164
1165 CALL keyword_create(keyword, __location__, name="GAPW_MAX_ALPHA", &
1166 description="Maximum alpha of GTH potentials allowed on the soft grids ", &
1167 usage="GAPW_MAX_ALPHA real", default_r_val=5.0_dp)
1168 CALL section_add_keyword(print_key, keyword)
1169 CALL keyword_release(keyword)
1170
1171 CALL keyword_create(keyword, __location__, name="SOO_RHO_HARD", &
1172 description="Whether or not to include the atomic parts of the density "// &
1173 "in the SOO part of the g tensor", usage="SOO_RHO_HARD", &
1174 default_l_val=.false., lone_keyword_l_val=.true.)
1175 CALL section_add_keyword(print_key, keyword)
1176 CALL keyword_release(keyword)
1177
1178 CALL section_add_subsection(subsection, print_key)
1179 CALL section_release(print_key)
1180
1181 CALL cp_print_key_section_create(print_key, __location__, "RESPONSE_FUNCTION_CUBES", &
1182 description="Controls the printing of the response functions ", &
1183 print_level=high_print_level, add_last=add_last_numeric, filename="")
1184 CALL keyword_create(keyword, __location__, name="stride", &
1185 description="The stride (X,Y,Z) used to write the cube file "// &
1186 "(larger values result in smaller cube files). You can provide 3 numbers (for X,Y,Z) or"// &
1187 " 1 number valid for all components.", &
1188 usage="STRIDE 2 2 2", n_var=-1, default_i_vals=(/2, 2, 2/), type_of_var=integer_t)
1189 CALL section_add_keyword(print_key, keyword)
1190 CALL keyword_release(keyword)
1191
1192 CALL keyword_create(keyword, __location__, name="CUBES_LU_BOUNDS", &
1193 variants=(/"CUBES_LU"/), &
1194 description="The lower and upper index of the states to be printed as cube", &
1195 usage="CUBES_LU_BOUNDS integer integer", &
1196 n_var=2, default_i_vals=(/0, -2/), type_of_var=integer_t)
1197 CALL section_add_keyword(print_key, keyword)
1198 CALL keyword_release(keyword)
1199
1200 CALL keyword_create(keyword, __location__, name="CUBES_LIST", &
1201 description="Indexes of the states to be printed as cube files "// &
1202 "This keyword can be repeated several times "// &
1203 "(useful if you have to specify many indexes).", &
1204 usage="CUBES_LIST 1 2", &
1205 n_var=-1, type_of_var=integer_t, repeats=.true.)
1206 CALL section_add_keyword(print_key, keyword)
1207 CALL keyword_release(keyword)
1208 CALL keyword_create(keyword, __location__, name="APPEND", &
1209 description="append the cube files when they already exist", &
1210 default_l_val=.false., lone_keyword_l_val=.true.)
1211 CALL section_add_keyword(print_key, keyword)
1212 CALL keyword_release(keyword)
1213
1214 CALL section_add_subsection(subsection, print_key)
1215 CALL section_release(print_key)
1216
1217 CALL section_add_subsection(section, subsection)
1218 CALL section_release(subsection)
1219
1220 NULLIFY (subsection)
1221 CALL create_interp_section(subsection)
1222 CALL section_add_subsection(section, subsection)
1223 CALL section_release(subsection)
1224
1225 END SUBROUTINE create_epr_section
1226
1227! **************************************************************************************************
1228!> \brief creates the input structure used to activate
1229!> calculation of polarizability tensor DFPT
1230!> Available properties : none
1231!> \param section the section to create
1232!> \author SL
1233! **************************************************************************************************
1234 SUBROUTINE create_polarizability_section(section)
1235
1236 TYPE(section_type), POINTER :: section
1237
1238 TYPE(keyword_type), POINTER :: keyword
1239 TYPE(section_type), POINTER :: print_key, subsection
1240
1241 NULLIFY (keyword, print_key, subsection)
1242
1243 cpassert(.NOT. ASSOCIATED(section))
1244 CALL section_create(section, __location__, name="POLAR", &
1245 description="Compute polarizabilities.", &
1246 n_keywords=5, n_subsections=1, repeats=.false., &
1247 citations=(/putrino2002/))
1248
1249 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
1250 description="controls the activation of the polarizability calculation", &
1251 usage="&POLAR T", &
1252 default_l_val=.false., &
1253 lone_keyword_l_val=.true.)
1254 CALL section_add_keyword(section, keyword)
1255 CALL keyword_release(keyword)
1256
1257 CALL keyword_create(keyword, __location__, name="DO_RAMAN", &
1258 description="Compute the electric-dipole--electric-dipole polarizability", &
1259 usage="DO_RAMAN F", &
1260 citations=(/luber2014/), &
1261 default_l_val=.true., lone_keyword_l_val=.true.)
1262 CALL section_add_keyword(section, keyword)
1263 CALL keyword_release(keyword)
1264
1265 CALL keyword_create(keyword, __location__, name="PERIODIC_DIPOLE_OPERATOR", &
1266 description="Type of dipole operator: Berry phase(T) or Local(F)", &
1267 usage="PERIODIC_DIPOLE_OPERATOR T", &
1268 default_l_val=.true., lone_keyword_l_val=.true.)
1269 CALL section_add_keyword(section, keyword)
1270 CALL keyword_release(keyword)
1271
1272 NULLIFY (subsection)
1273 CALL section_create(subsection, __location__, name="PRINT", &
1274 description="print results of the polarizability calculation", &
1275 repeats=.false.)
1276
1277 CALL cp_print_key_section_create(print_key, __location__, "POLAR_MATRIX", &
1278 description="Controls the printing of the polarizabilities", &
1279 print_level=low_print_level, add_last=add_last_numeric, filename="")
1280
1281 CALL section_add_subsection(subsection, print_key)
1282 CALL section_release(print_key)
1283 CALL section_add_subsection(section, subsection)
1284 CALL section_release(subsection)
1285
1286 NULLIFY (subsection)
1287 CALL create_interp_section(subsection)
1288 CALL section_add_subsection(section, subsection)
1289 CALL section_release(subsection)
1290
1291 END SUBROUTINE create_polarizability_section
1292
1293! **************************************************************************************************
1294!> \brief creates the section for electron transfer coupling
1295!> \param section ...
1296!> \author fschiff
1297! **************************************************************************************************
1298 SUBROUTINE create_et_coupling_section(section)
1299 TYPE(section_type), POINTER :: section
1300
1301 TYPE(keyword_type), POINTER :: keyword
1302 TYPE(section_type), POINTER :: print_key, subsection
1303
1304 NULLIFY (keyword)
1305 cpassert(.NOT. ASSOCIATED(section))
1306 CALL section_create(section, __location__, name="ET_COUPLING", &
1307 description="specifies the two constraints/restraints for extracting ET coupling elements", &
1308 n_keywords=1, n_subsections=4, repeats=.false., citations=(/kondov2007, futera2017/))
1309
1310 NULLIFY (subsection)
1311 CALL create_ddapc_restraint_section(subsection, "DDAPC_RESTRAINT_A")
1312 CALL section_add_subsection(section, subsection)
1313 CALL section_release(subsection)
1314
1315 NULLIFY (subsection)
1316 CALL create_ddapc_restraint_section(subsection, "DDAPC_RESTRAINT_B")
1317 CALL section_add_subsection(section, subsection)
1318 CALL section_release(subsection)
1319
1320 NULLIFY (subsection)
1321 CALL create_projection(subsection, "PROJECTION")
1322 CALL section_add_subsection(section, subsection)
1323 CALL section_release(subsection)
1324
1325 CALL keyword_create(keyword, __location__, name="TYPE_OF_CONSTRAINT", &
1326 description="Specifies the type of constraint", &
1327 usage="TYPE_OF_CONSTRAINT DDAPC", &
1328 enum_c_vals=s2a("NONE", "DDAPC"), &
1329 enum_i_vals=(/do_no_et, do_et_ddapc/), &
1330 enum_desc=s2a("NONE", "DDAPC Constraint"), &
1331 default_i_val=do_no_et)
1332 CALL section_add_keyword(section, keyword)
1333 CALL keyword_release(keyword)
1334
1335 NULLIFY (print_key)
1336 CALL cp_print_key_section_create(print_key, __location__, "PROGRAM_RUN_INFO", &
1337 description="Controls the printing basic info about the method", &
1338 print_level=low_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
1339 CALL section_add_subsection(section, print_key)
1340 CALL section_release(print_key)
1341
1342 END SUBROUTINE create_et_coupling_section
1343
1344! **************************************************************************************************
1345!> \brief defines input sections for specification of Hilbert space partitioning
1346!> in projection-operator approach of electronic coupling calulation
1347!> \param section pointer to the section data structure
1348!> \param section_name name of the projection section
1349!> \author Z. Futera (02.2017)
1350! **************************************************************************************************
1351 SUBROUTINE create_projection(section, section_name)
1352
1353 ! Routine arguments
1354 TYPE(section_type), POINTER :: section
1355 CHARACTER(len=*), INTENT(in) :: section_name
1356
1357 TYPE(keyword_type), POINTER :: keyword
1358 TYPE(section_type), POINTER :: print_key, section_block, section_print
1359
1360! Routine name for dubug purposes
1361
1362 ! Sanity check
1363 cpassert(.NOT. ASSOCIATED(section))
1364
1365 ! Initialization
1366 NULLIFY (keyword)
1367 NULLIFY (print_key)
1368 NULLIFY (section_block)
1369 NULLIFY (section_print)
1370
1371 ! Input-file section definition
1372 CALL section_create(section, __location__, name=trim(adjustl(section_name)), &
1373 description="Projection-operator approach fo ET coupling calculation", &
1374 n_keywords=0, n_subsections=2, repeats=.false.)
1375
1376 ! Subsection #0: Log printing
1377 CALL cp_print_key_section_create(print_key, __location__, 'PROGRAM_RUN_INFO', &
1378 description="Controls printing of data and informations to log file", &
1379 print_level=low_print_level, filename="__STD_OUT__")
1380 CALL section_add_subsection(section, print_key)
1381 CALL section_release(print_key)
1382
1383 ! Subsection #1: Atomic blocks
1384 CALL section_create(section_block, __location__, name='BLOCK', &
1385 description="Part of the system (donor, acceptor, bridge,...)", &
1386 n_keywords=2, n_subsections=1, repeats=.true.)
1387 CALL section_add_subsection(section, section_block)
1388
1389 ! S#1 - Keyword #1: Atom IDs defining a Hilbert space block
1390 CALL keyword_create(keyword, __location__, name='ATOMS', &
1391 description="Array of atom IDs in the system part", &
1392 usage="ATOMS {integer} {integer} .. {integer}", &
1393 n_var=-1, type_of_var=integer_t, repeats=.false.)
1394 CALL section_add_keyword(section_block, keyword)
1395 CALL keyword_release(keyword)
1396
1397 ! S#1 - Keyword #1: Atom IDs defining a Hilbert space block
1398 CALL keyword_create(keyword, __location__, name='NELECTRON', &
1399 description="Number of electrons expected in the system part", &
1400 usage="NELECTRON {integer}", default_i_val=0)
1401 CALL section_add_keyword(section_block, keyword)
1402 CALL keyword_release(keyword)
1403
1404 ! S#1 - Subsection #1: Printing setting
1405 CALL section_create(section_print, __location__, name='PRINT', &
1406 description="Possible printing options in ET system part", &
1407 n_keywords=0, n_subsections=0, repeats=.false.)
1408 CALL section_add_subsection(section_block, section_print)
1409
1410 ! S#1 - S#1 - Keyword #1: MO coefficient on specific atom
1411 CALL keyword_create(keyword, __location__, name='MO_COEFF_ATOM', &
1412 description="Print out MO coeffiecients on given atom", &
1413 usage="MO_COEFF_ATOM {integer} {integer} .. {integer}", &
1414 type_of_var=integer_t, n_var=-1, repeats=.true.)
1415 CALL section_add_keyword(section_print, keyword)
1416 CALL keyword_release(keyword)
1417
1418 ! S#1 - S#1 - Keyword #1: MO coefficient of specific state
1419 CALL keyword_create(keyword, __location__, name='MO_COEFF_ATOM_STATE', &
1420 description="Print out MO coeffiecients of specific state", &
1421 usage="MO_COEFF_ATOM_STATE {integer} {integer} .. {integer}", &
1422 type_of_var=integer_t, n_var=-1, repeats=.true.)
1423 CALL section_add_keyword(section_print, keyword)
1424 CALL keyword_release(keyword)
1425
1426 ! S#1 - S#1 - Subsection #1: Saving MOs to CUBE files
1427 CALL cp_print_key_section_create(print_key, __location__, 'MO_CUBES', &
1428 description="Controls saving of MO cube files", &
1429 print_level=high_print_level, filename="")
1430
1431 ! S#1 - S#1 - S#1 - Keyword #1: Stride
1432 CALL keyword_create(keyword, __location__, name='STRIDE', &
1433 description="The stride (X,Y,Z) used to write the cube file", &
1434 usage="STRIDE {integer} {integer} {integer}", n_var=-1, &
1435 default_i_vals=(/2, 2, 2/), type_of_var=integer_t)
1436 CALL section_add_keyword(print_key, keyword)
1437 CALL keyword_release(keyword)
1438
1439 ! S#1 - S#1 - S#1 - Keyword #2: List of MO IDs
1440 CALL keyword_create(keyword, __location__, name='MO_LIST', &
1441 description="Indices of molecular orbitals to save", &
1442 usage="MO_LIST {integer} {integer} .. {integer}", &
1443 type_of_var=integer_t, n_var=-1, repeats=.true.)
1444 CALL section_add_keyword(print_key, keyword)
1445 CALL keyword_release(keyword)
1446
1447 ! S#1 - S#1 - S#1 - Keyword #2: Number of unoccupied states
1448 CALL keyword_create(keyword, __location__, name='NLUMO', &
1449 description="Number of unoccupied molecular orbitals to save", &
1450 usage="NLUMO {integer}", default_i_val=1)
1451 CALL section_add_keyword(print_key, keyword)
1452 CALL keyword_release(keyword)
1453
1454 ! S#1 - S#1 - S#1 - Keyword #3: Number of occupied states
1455 CALL keyword_create(keyword, __location__, name='NHOMO', &
1456 description="Number of occupied molecular orbitals to save", &
1457 usage="NHOMO {integer}", default_i_val=1)
1458 CALL section_add_keyword(print_key, keyword)
1459 CALL keyword_release(keyword)
1460
1461 CALL section_add_subsection(section_print, print_key)
1462 CALL section_release(print_key)
1463
1464 ! S#1 - S#1 - Clean
1465 CALL section_release(section_print)
1466
1467 ! S#1 - Clean
1468 CALL section_release(section_block)
1469
1470 ! S#1 - Subsection #1: Printing setting
1471 CALL section_create(section_print, __location__, name='PRINT', &
1472 description="Possible printing options in ET", &
1473 n_keywords=0, n_subsections=0, repeats=.false.)
1474 CALL section_add_subsection(section, section_print)
1475
1476 ! Print couplings
1477 CALL cp_print_key_section_create(print_key, __location__, 'COUPLINGS', &
1478 description="Controls printing couplings onto file", &
1479 print_level=low_print_level, filename="")
1480
1481 CALL keyword_create(keyword, __location__, name="APPEND", &
1482 description="append the files when they already exist", &
1483 default_l_val=.false., lone_keyword_l_val=.true.)
1484 CALL section_add_keyword(print_key, keyword)
1485 CALL keyword_release(keyword)
1486
1487 CALL section_add_subsection(section_print, print_key)
1488 CALL section_release(print_key)
1489
1490 CALL section_release(section_print)
1491
1492 END SUBROUTINE create_projection
1493
1494! **************************************************************************************************
1495!> \brief creates an input section for tddfpt calculation
1496!> \param section section to create
1497!> \par History
1498!> * 05.2016 forked from create_tddfpt_section [Sergey Chulkov]
1499!> * 08.2016 moved from module input_cp2k_dft [Sergey Chulkov]
1500! **************************************************************************************************
1501 SUBROUTINE create_tddfpt2_section(section)
1502 TYPE(section_type), POINTER :: section
1503
1504 TYPE(keyword_type), POINTER :: keyword
1505 TYPE(section_type), POINTER :: print_key, subsection
1506
1507 cpassert(.NOT. ASSOCIATED(section))
1508 CALL section_create(section, __location__, name="TDDFPT", &
1509 description="Parameters needed to set up the Time-Dependent "// &
1510 "Density Functional Perturbation Theory. "// &
1511 "Current implementation works for hybrid functionals. "// &
1512 "Can be used with Gaussian and Plane Waves (GPW) method only.", &
1513 n_keywords=14, n_subsections=4, repeats=.false., &
1514 citations=(/iannuzzi2005/))
1515
1516 NULLIFY (keyword, print_key, subsection)
1517
1518 CALL keyword_create(keyword, __location__, &
1519 name="_SECTION_PARAMETERS_", &
1520 description="Controls the activation of the TDDFPT procedure", &
1521 default_l_val=.false., &
1522 lone_keyword_l_val=.true.)
1523 CALL section_add_keyword(section, keyword)
1524 CALL keyword_release(keyword)
1525
1526 ! Integer
1527 CALL keyword_create(keyword, __location__, name="NSTATES", &
1528 description="Number of excited states to converge.", &
1529 n_var=1, type_of_var=integer_t, &
1530 default_i_val=1)
1531 CALL section_add_keyword(section, keyword)
1532 CALL keyword_release(keyword)
1533
1534 CALL keyword_create(keyword, __location__, name="MAX_ITER", &
1535 description="Maximal number of iterations to be performed.", &
1536 n_var=1, type_of_var=integer_t, &
1537 default_i_val=50)
1538 CALL section_add_keyword(section, keyword)
1539 CALL keyword_release(keyword)
1540
1541 CALL keyword_create(keyword, __location__, name="MAX_KV", &
1542 description="Maximal number of Krylov space vectors. "// &
1543 "Davidson iterations will be restarted upon reaching this limit.", &
1544 n_var=1, type_of_var=integer_t, &
1545 default_i_val=5000)
1546 CALL section_add_keyword(section, keyword)
1547 CALL keyword_release(keyword)
1548
1549 CALL keyword_create(keyword, __location__, name="NLUMO", &
1550 description="Number of unoccupied orbitals to consider. "// &
1551 "Default is to use all unoccupied orbitals (-1).", &
1552 n_var=1, type_of_var=integer_t, &
1553 default_i_val=-1)
1554 CALL section_add_keyword(section, keyword)
1555 CALL keyword_release(keyword)
1556
1557 CALL keyword_create(keyword, __location__, name="NPROC_STATE", &
1558 description="Number of MPI processes to be used per excited state. "// &
1559 "Default is to use all processors (0).", &
1560 n_var=1, type_of_var=integer_t, &
1561 default_i_val=0)
1562 CALL section_add_keyword(section, keyword)
1563 CALL keyword_release(keyword)
1564
1565 ! kernel type
1566 CALL keyword_create(keyword, __location__, name="KERNEL", &
1567 description="Options to compute the kernel", &
1568 usage="KERNEL FULL", &
1569 enum_c_vals=s2a("FULL", "sTDA", "NONE"), &
1571 default_i_val=tddfpt_kernel_full)
1572 CALL section_add_keyword(section, keyword)
1573 CALL keyword_release(keyword)
1574
1575 ! spin-flip TDDFPT options
1576 CALL keyword_create(keyword, __location__, name="SPINFLIP", &
1577 description="Selects the type of spin-flip TDDFPT kernel", &
1578 usage="SPINFLIP NONCOLLINEAR", &
1579 enum_c_vals=s2a("NONE", "COLLINEAR", "NONCOLLINEAR"), &
1580 enum_i_vals=(/no_sf_tddfpt, tddfpt_sf_col, tddfpt_sf_noncol/), &
1581 default_i_val=no_sf_tddfpt)
1582 CALL section_add_keyword(section, keyword)
1583 CALL keyword_release(keyword)
1584
1585 CALL keyword_create(keyword, __location__, name="OE_CORR", &
1586 description="Orbital energy correction potential.", &
1587 enum_c_vals=s2a("NONE", "LB94", "GLLB", "SAOP", "SHIFT"), &
1588 enum_i_vals=(/oe_none, oe_lb, oe_gllb, oe_saop, oe_shift/), &
1589 enum_desc=s2a("No orbital correction scheme is used", &
1590 "van Leeuwen and Baerends. PRA, 49:2421, 1994", &
1591 "Gritsenko, van Leeuwen, van Lenthe, Baerends. PRA, 51:1944, 1995", &
1592 "Gritsenko, Schipper, Baerends. Chem. Phys. Lett., 302:199, 1999", &
1593 "Constant shift of virtual and/or open-shell orbitals"), &
1594 default_i_val=oe_none)
1595 CALL section_add_keyword(section, keyword)
1596 CALL keyword_release(keyword)
1597
1598 ! SHIFTS
1599 CALL keyword_create(keyword, __location__, name="EV_SHIFT", &
1600 variants=s2a("VIRTUAL_SHIFT"), &
1601 description="Constant shift of virtual state eigenvalues.", &
1602 usage="EV_SHIFT 0.500", &
1603 n_var=1, type_of_var=real_t, &
1604 unit_str="eV", &
1605 default_r_val=0.0_dp)
1606 CALL section_add_keyword(section, keyword)
1607 CALL keyword_release(keyword)
1608 !
1609 CALL keyword_create(keyword, __location__, name="EOS_SHIFT", &
1610 variants=s2a("OPEN_SHELL_SHIFT"), &
1611 description="Constant shift of open shell eigenvalues.", &
1612 usage="EOS_SHIFT 0.200", &
1613 n_var=1, type_of_var=real_t, &
1614 unit_str="eV", &
1615 default_r_val=0.0_dp)
1616 CALL section_add_keyword(section, keyword)
1617 CALL keyword_release(keyword)
1618
1619 ! Real
1620 CALL keyword_create(keyword, __location__, name="CONVERGENCE", &
1621 description="Target accuracy for excited state energies.", &
1622 n_var=1, type_of_var=real_t, unit_str="hartree", &
1623 default_r_val=1.0e-5_dp)
1624 CALL section_add_keyword(section, keyword)
1625 CALL keyword_release(keyword)
1626
1627 CALL keyword_create(keyword, __location__, name="MIN_AMPLITUDE", &
1628 description="The smallest excitation amplitude to print.", &
1629 n_var=1, type_of_var=real_t, &
1630 default_r_val=5.0e-2_dp)
1631 CALL section_add_keyword(section, keyword)
1632 CALL keyword_release(keyword)
1633
1634 CALL keyword_create(keyword, __location__, name="ORTHOGONAL_EPS", &
1635 description="The largest possible overlap between the ground state and "// &
1636 "orthogonalised excited state wave-functions. Davidson iterations "// &
1637 "will be restarted when the overlap goes beyond this threshold in "// &
1638 "order to prevent numerical instability.", &
1639 n_var=1, type_of_var=real_t, &
1640 default_r_val=1.0e-4_dp)
1641 CALL section_add_keyword(section, keyword)
1642 CALL keyword_release(keyword)
1643
1644 ! Logical
1645 CALL keyword_create(keyword, __location__, name="RESTART", &
1646 description="Restart the TDDFPT calculation if a restart file exists", &
1647 n_var=1, type_of_var=logical_t, &
1648 default_l_val=.false., lone_keyword_l_val=.true.)
1649 CALL section_add_keyword(section, keyword)
1650 CALL keyword_release(keyword)
1651
1652 CALL keyword_create(keyword, __location__, name="RKS_TRIPLETS", &
1653 description="Compute triplet excited states using spin-unpolarised molecular orbitals.", &
1654 n_var=1, type_of_var=logical_t, &
1655 default_l_val=.false.)
1656 CALL section_add_keyword(section, keyword)
1657 CALL keyword_release(keyword)
1658
1659 CALL keyword_create(keyword, __location__, name="ADMM_KERNEL_XC_CORRECTION", &
1660 description="Use/Ignore ADMM correction xc functional for TD kernel. "// &
1661 "XC correction functional is defined in ground state XC section.", &
1662 n_var=1, type_of_var=logical_t, &
1663 default_l_val=.true., lone_keyword_l_val=.true.)
1664 CALL section_add_keyword(section, keyword)
1665 CALL keyword_release(keyword)
1666
1667 CALL keyword_create(keyword, __location__, name="ADMM_KERNEL_CORRECTION_SYMMETRIC", &
1668 description="ADMM correction functional in kernel is applied symmetrically. "// &
1669 "Original implementation is using a non-symmetric formula.", &
1670 n_var=1, type_of_var=logical_t, &
1671 default_l_val=.true., lone_keyword_l_val=.true.)
1672 CALL section_add_keyword(section, keyword)
1673 CALL keyword_release(keyword)
1674
1675 CALL keyword_create(keyword, __location__, name="DO_LRIGPW", &
1676 description="Local resolution of identity for Coulomb contribution.", &
1677 n_var=1, type_of_var=logical_t, &
1678 default_l_val=.false.)
1679 CALL section_add_keyword(section, keyword)
1680 CALL keyword_release(keyword)
1681
1682 CALL keyword_create(keyword, __location__, name="AUTO_BASIS", &
1683 description="Specify size of automatically generated auxiliary basis sets: "// &
1684 "Options={small,medium,large,huge}", &
1685 usage="AUTO_BASIS {basis_type} {basis_size}", &
1686 type_of_var=char_t, repeats=.true., n_var=-1, default_c_vals=(/"X", "X"/))
1687 CALL section_add_keyword(section, keyword)
1688 CALL keyword_release(keyword)
1689
1690 CALL keyword_create(keyword, __location__, name="DO_SMEARING", &
1691 description="Implying smeared occupation. ", &
1692 n_var=1, type_of_var=logical_t, &
1693 default_l_val=.false., lone_keyword_l_val=.true.)
1694 CALL section_add_keyword(section, keyword)
1695 CALL keyword_release(keyword)
1696
1697 CALL keyword_create(keyword, __location__, name="EXCITON_DESCRIPTORS", &
1698 description="Compute exciton descriptors. "// &
1699 "Details given in Manual section about Bethe Salpeter equation.", &
1700 n_var=1, type_of_var=logical_t, &
1701 default_l_val=.false.)
1702 CALL section_add_keyword(section, keyword)
1703 CALL keyword_release(keyword)
1704
1705 CALL keyword_create(keyword, __location__, name="DIRECTIONAL_EXCITON_DESCRIPTORS", &
1706 description="Print cartesian components of exciton descriptors.", &
1707 n_var=1, type_of_var=logical_t, &
1708 default_l_val=.false.)
1709 CALL section_add_keyword(section, keyword)
1710 CALL keyword_release(keyword)
1711
1712 ! Strings
1713 CALL keyword_create(keyword, __location__, name="WFN_RESTART_FILE_NAME", &
1714 variants=(/"RESTART_FILE_NAME"/), &
1715 description="Name of the wave function restart file, may include a path."// &
1716 " If no file is specified, the default is to open the file as generated by"// &
1717 " the wave function restart print key.", &
1718 usage="WFN_RESTART_FILE_NAME <FILENAME>", &
1719 type_of_var=lchar_t)
1720 CALL section_add_keyword(section, keyword)
1721 CALL keyword_release(keyword)
1722
1723 ! DIPOLE subsection
1724 CALL section_create(subsection, __location__, name="DIPOLE_MOMENTS", &
1725 description="Parameters to compute oscillator strengths in the dipole approximation.", &
1726 n_keywords=3, n_subsections=0, repeats=.false.)
1727
1728 CALL keyword_create(keyword, __location__, name="DIPOLE_FORM", &
1729 description="Form of dipole transition integrals.", &
1730 enum_c_vals=s2a("BERRY", "LENGTH", "VELOCITY"), &
1731 enum_desc=s2a("Based on Berry phase formula (valid for fully periodic molecular systems only)", &
1732 "Length form &lang; i | r | j &rang; (valid for non-periodic molecular systems only)", &
1733 "Velocity form &lang; i | d/dr | j &rang;"), &
1735 default_i_val=tddfpt_dipole_velocity)
1736 CALL section_add_keyword(subsection, keyword)
1737 CALL keyword_release(keyword)
1738
1739 CALL keyword_create(keyword, __location__, name="REFERENCE", &
1740 description="Reference point to calculate electric "// &
1741 "dipole moments using the dipole integrals in the length form.", &
1742 enum_c_vals=s2a("COM", "COAC", "USER_DEFINED", "ZERO"), &
1743 enum_desc=s2a("Use Center of Mass", &
1744 "Use Center of Atomic Charges", &
1745 "Use User-defined Point", &
1746 "Use Origin of Coordinate System"), &
1747 enum_i_vals=(/use_mom_ref_com, &
1750 use_mom_ref_zero/), &
1751 default_i_val=use_mom_ref_com)
1752 CALL section_add_keyword(subsection, keyword)
1753 CALL keyword_release(keyword)
1754
1755 CALL keyword_create(keyword, __location__, name="REFERENCE_POINT", &
1756 description="User-defined reference point.", &
1757 usage="REFERENCE_POINT x y z", &
1758 repeats=.false., n_var=3, type_of_var=real_t, unit_str='bohr')
1759 CALL section_add_keyword(subsection, keyword)
1760 CALL keyword_release(keyword)
1761
1762 CALL section_add_subsection(section, subsection)
1763 CALL section_release(subsection)
1764
1765 ! SOC functional
1766
1767 CALL section_create(subsection, __location__, name="SOC", &
1768 description="Is jet to be implemented", &
1769 n_keywords=2, n_subsections=0, repeats=.false.)
1770
1771 CALL keyword_create(keyword, __location__, name="EPS_FILTER", &
1772 variants=s2a("EPS_FILTER_MATRIX"), &
1773 description="The threshold used for sparse matrix operations", &
1774 usage="EPS_FILTER {real}", &
1775 type_of_var=real_t, &
1776 default_r_val=1.0e-10_dp)
1777 CALL section_add_keyword(subsection, keyword)
1778 CALL keyword_release(keyword)
1779
1780 CALL keyword_create(keyword, __location__, name="GRID", &
1781 variants=(/"ATOMIC_GRID"/), &
1782 description="Specification of the atomic angular and radial grids for "// &
1783 "a atomic kind. This keyword must be repeated for all kinds! "// &
1784 "Usage: GRID < LEBEDEV_GRID > < RADIAL_GRID >", &
1785 usage="GRID {string} {integer} {integer}", &
1786 n_var=3, type_of_var=char_t, repeats=.true.)
1787 CALL section_add_keyword(subsection, keyword)
1788 CALL keyword_release(keyword)
1789
1790 CALL section_add_subsection(section, subsection)
1791 CALL section_release(subsection)
1792
1793 ! kernel XC functional
1794 CALL create_xc_section(subsection)
1795 CALL section_add_subsection(section, subsection)
1796 CALL section_release(subsection)
1797
1798 ! MGRID subsection
1799 CALL create_mgrid_section(subsection, create_subsections=.false.)
1800 CALL section_add_subsection(section, subsection)
1801 CALL section_release(subsection)
1802
1803 ! sTDA subsection
1804 CALL create_stda_section(subsection)
1805 CALL section_add_subsection(section, subsection)
1806 CALL section_release(subsection)
1807
1808 CALL keyword_create(keyword, __location__, name="DO_BSE", &
1809 description="Choosing BSE kernel.", &
1810 usage="DO_BSE", default_l_val=.false., lone_keyword_l_val=.true.)
1811 CALL section_add_keyword(section, keyword)
1812 CALL keyword_release(keyword)
1813
1814 ! LRI subsection
1815 CALL create_lrigpw_section(subsection)
1816 CALL section_add_subsection(section, subsection)
1817 CALL section_release(subsection)
1818
1819 ! LINRES section
1820 CALL create_linres_section(subsection, create_subsections=.false., default_set_tdlr=.true.)
1821 CALL section_add_subsection(section, subsection)
1822 CALL section_release(subsection)
1823
1824 ! PRINT subsection
1825 CALL section_create(subsection, __location__, name="PRINT", &
1826 description="Printing of information during the TDDFT run.", repeats=.false.)
1827
1828 CALL cp_print_key_section_create(print_key, __location__, name="PROGRAM_BANNER", &
1829 description="Controls the printing of the banner for TDDFPT program", &
1830 print_level=silent_print_level, filename="__STD_OUT__")
1831 CALL section_add_subsection(subsection, print_key)
1832 CALL section_release(print_key)
1833
1834 CALL cp_print_key_section_create(print_key, __location__, name="GUESS_VECTORS", &
1835 description="Controls the printing of initial guess vectors.", &
1836 print_level=low_print_level, filename="__STD_OUT__")
1837 CALL section_add_subsection(subsection, print_key)
1838 CALL section_release(print_key)
1839
1840 CALL cp_print_key_section_create(print_key, __location__, name="ITERATION_INFO", &
1841 description="Controls the printing of basic iteration information "// &
1842 "during the TDDFT run.", &
1843 print_level=low_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
1844 CALL section_add_subsection(subsection, print_key)
1845 CALL section_release(print_key)
1846
1847 CALL cp_print_key_section_create(print_key, __location__, name="DETAILED_ENERGY", &
1848 description="Controls the printing of detailed energy information "// &
1849 "during the TDDFT run.", &
1850 print_level=medium_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
1851 CALL section_add_subsection(subsection, print_key)
1852 CALL section_release(print_key)
1853
1854 CALL cp_print_key_section_create(print_key, __location__, name="BASIS_SET_FILE", &
1855 description="Controls the printing of a file with all basis sets used.", &
1856 print_level=debug_print_level, filename="BASIS_SETS")
1857 CALL section_add_subsection(subsection, print_key)
1858 CALL section_release(print_key)
1859
1860 CALL cp_print_key_section_create(print_key, __location__, name="RESTART", &
1861 description="Controls the dumping of the MO restart file during TDDFPT. "// &
1862 "By default keeps a short history of three restarts.", &
1863 print_level=low_print_level, common_iter_levels=3, &
1864 each_iter_names=s2a("TDDFT_SCF"), each_iter_values=(/10/), &
1865 add_last=add_last_numeric, filename="RESTART")
1866 CALL keyword_create(keyword, __location__, name="BACKUP_COPIES", &
1867 description="Specifies the maximum number of backup copies.", &
1868 usage="BACKUP_COPIES {int}", &
1869 default_i_val=1)
1870 CALL section_add_keyword(print_key, keyword)
1871 CALL keyword_release(keyword)
1872 CALL section_add_subsection(subsection, print_key)
1873 CALL section_release(print_key)
1874
1875 CALL cp_print_key_section_create(print_key, __location__, name="NTO_ANALYSIS", &
1876 description="Perform a natural transition orbital analysis.", &
1877 print_level=medium_print_level)
1878 CALL keyword_create(keyword, __location__, name="THRESHOLD", &
1879 description="Threshold for sum of NTO eigenvalues considered", &
1880 usage="Threshold 0.95", &
1881 n_var=1, &
1882 type_of_var=real_t, &
1883 default_r_val=0.975_dp)
1884 CALL section_add_keyword(print_key, keyword)
1885 CALL keyword_release(keyword)
1886 CALL keyword_create(keyword, __location__, name="INTENSITY_THRESHOLD", &
1887 description="Threshold for oscillator strength to screen states.", &
1888 usage="Intensity_threshold 0.01", &
1889 n_var=1, &
1890 type_of_var=real_t, &
1891 default_r_val=0.0_dp)
1892 CALL section_add_keyword(print_key, keyword)
1893 CALL keyword_release(keyword)
1894 CALL keyword_create(keyword, __location__, name="STATE_LIST", &
1895 description="Specifies a list of states for the NTO calculations.", &
1896 usage="STATE_LIST {integer} {integer} .. {integer}", &
1897 n_var=-1, type_of_var=integer_t)
1898 CALL section_add_keyword(print_key, keyword)
1899 CALL keyword_release(keyword)
1900 CALL keyword_create(keyword, __location__, name="CUBE_FILES", &
1901 description="Print NTOs on Cube Files", &
1902 usage="CUBE_FILES {logical}", repeats=.false., n_var=1, &
1903 default_l_val=.false., lone_keyword_l_val=.true., type_of_var=logical_t)
1904 CALL section_add_keyword(print_key, keyword)
1905 CALL keyword_release(keyword)
1906 CALL keyword_create(keyword, __location__, name="STRIDE", &
1907 description="The stride (X,Y,Z) used to write the cube file "// &
1908 "(larger values result in smaller cube files). Provide 3 numbers (for X,Y,Z) or"// &
1909 " 1 number valid for all components.", &
1910 usage="STRIDE 2 2 2", n_var=-1, default_i_vals=(/2, 2, 2/), type_of_var=integer_t)
1911 CALL section_add_keyword(print_key, keyword)
1912 CALL keyword_release(keyword)
1913 CALL keyword_create(keyword, __location__, name="APPEND", &
1914 description="append the cube files when they already exist", &
1915 default_l_val=.false., lone_keyword_l_val=.true.)
1916 CALL section_add_keyword(print_key, keyword)
1917 CALL keyword_release(keyword)
1918 CALL section_add_subsection(subsection, print_key)
1919 CALL section_release(print_key)
1920
1921 CALL cp_print_key_section_create(print_key, __location__, "MOS_MOLDEN", &
1922 description="Write the NTO in Molden file format, for visualisation.", &
1923 print_level=debug_print_level + 1, add_last=add_last_numeric, filename="MOS")
1924 CALL keyword_create(keyword, __location__, name="NDIGITS", &
1925 description="Specifies the number of significant digits retained. 3 is OK for visualization.", &
1926 usage="NDIGITS {int}", &
1927 default_i_val=3)
1928 CALL section_add_keyword(print_key, keyword)
1929 CALL keyword_release(keyword)
1930 CALL keyword_create(keyword, __location__, name="GTO_KIND", &
1931 description="Representation of Gaussian-type orbitals", &
1932 default_i_val=gto_spherical, &
1933 enum_c_vals=s2a("CARTESIAN", "SPHERICAL"), &
1934 enum_desc=s2a( &
1935 "Cartesian Gaussian orbitals. Use with caution", &
1936 "Spherical Gaussian orbitals. Incompatible with VMD"), &
1937 enum_i_vals=(/gto_cartesian, gto_spherical/))
1938 CALL section_add_keyword(print_key, keyword)
1939 CALL keyword_release(keyword)
1940 CALL section_add_subsection(subsection, print_key)
1941 CALL section_release(print_key)
1942
1943 CALL cp_print_key_section_create(print_key, __location__, name="NAMD_PRINT", &
1944 description="Controls the printout required for NAMD with NEWTONX.", &
1945 print_level=debug_print_level + 1, filename="CP2K_NEWTONX")
1946 CALL keyword_create(keyword, __location__, name="PRINT_VIRTUALS", &
1947 description="Print occupied AND virtual molecular orbital coefficients", &
1948 default_l_val=.false., lone_keyword_l_val=.true.)
1949 CALL section_add_keyword(print_key, keyword)
1950 CALL keyword_release(keyword)
1951 CALL keyword_create(keyword, __location__, name="PRINT_PHASES", &
1952 description="Print phases of occupied and virtuals MOs.", &
1953 default_l_val=.false., lone_keyword_l_val=.true.)
1954 CALL section_add_keyword(print_key, keyword)
1955 CALL keyword_release(keyword)
1956 CALL keyword_create(keyword, __location__, name="SCALE_WITH_PHASES", &
1957 description="Scale ES eigenvectors with phases of occupied and virtuals MOs.", &
1958 default_l_val=.false., lone_keyword_l_val=.true.)
1959 CALL section_add_keyword(print_key, keyword)
1960 CALL keyword_release(keyword)
1961 CALL section_add_subsection(subsection, print_key)
1962 CALL section_release(print_key)
1963
1964 !! SOC PRINT SECTION
1965 CALL cp_print_key_section_create(print_key, __location__, name="SOC_PRINT", &
1966 description="Controls the printout of the tddfpt2_soc modul", &
1967 print_level=debug_print_level + 1, filename="SOC")
1968 CALL keyword_create(keyword, __location__, name="UNIT_eV", &
1969 description="Will detrement if output in eVolt will be printef.", &
1970 default_l_val=.true., lone_keyword_l_val=.true.)
1971 CALL section_add_keyword(print_key, keyword)
1972 CALL keyword_release(keyword)
1973 CALL keyword_create(keyword, __location__, name="UNIT_wn", &
1974 description="Will detrement if output in wavenumbers will be printed.", &
1975 default_l_val=.false., lone_keyword_l_val=.true.)
1976 CALL section_add_keyword(print_key, keyword)
1977 CALL keyword_release(keyword)
1978 CALL keyword_create(keyword, __location__, name="SPLITTING", &
1979 description="Will add the SOC-Splitting as additional output", &
1980 default_l_val=.false., lone_keyword_l_val=.true.)
1981 CALL section_add_keyword(print_key, keyword)
1982 CALL keyword_release(keyword)
1983 CALL keyword_create(keyword, __location__, name="SOME", &
1984 description="Will add the SOC-Matrix as additional output in a different file", &
1985 default_l_val=.false., lone_keyword_l_val=.true.)
1986 CALL section_add_keyword(print_key, keyword)
1987 CALL keyword_release(keyword)
1988 CALL section_add_subsection(subsection, print_key)
1989 CALL section_release(print_key)
1990
1991 CALL cp_print_key_section_create(print_key, __location__, name="FORCES", &
1992 description="Controls the calculation and printing of excited state forces. "// &
1993 "This needs a RUN_TYPE that includes force evaluation, e.g. ENERGY_FORCE", &
1994 print_level=debug_print_level, filename="TDFORCE")
1995 CALL keyword_create(keyword, __location__, name="LIST", &
1996 description="Specifies a list of states for the force calculations.", &
1997 usage="LIST {integer} {integer} .. {integer}", repeats=.true., &
1998 n_var=-1, type_of_var=integer_t)
1999 CALL section_add_keyword(print_key, keyword)
2000 CALL keyword_release(keyword)
2001 CALL keyword_create(keyword, __location__, name="THRESHOLD", &
2002 description="Threshold for oszillator strength to screen states.", &
2003 usage="Threshold 0.01", &
2004 n_var=1, &
2005 type_of_var=real_t, &
2006 default_r_val=0.0_dp)
2007 CALL section_add_keyword(print_key, keyword)
2008 CALL keyword_release(keyword)
2009 CALL section_add_subsection(subsection, print_key)
2010 CALL section_release(print_key)
2011
2012 CALL section_add_subsection(section, subsection)
2013 CALL section_release(subsection)
2014
2015 END SUBROUTINE create_tddfpt2_section
2016
2017! **************************************************************************************************
2018!> \brief creates the stda input section (simplified Tamm Dancoff Approximation)
2019!> \param section the section to create
2020! **************************************************************************************************
2021 SUBROUTINE create_stda_section(section)
2022 TYPE(section_type), POINTER :: section
2023
2024 TYPE(keyword_type), POINTER :: keyword
2025
2026 cpassert(.NOT. ASSOCIATED(section))
2027 CALL section_create(section, __location__, name="sTDA", &
2028 description="parameters needed and setup for sTDA calculations", &
2029 n_keywords=3, n_subsections=0, repeats=.false.)
2030 NULLIFY (keyword)
2031
2032 CALL keyword_create(keyword, __location__, name="FRACTION", &
2033 variants=(/"HFX_FRACTION"/), &
2034 description="The fraction of TB Hartree-Fock exchange to use in the Kernel. "// &
2035 "0.0 implies no HFX part is used in the kernel. ", &
2036 usage="FRACTION 0.0", default_r_val=0.0_dp)
2037 CALL section_add_keyword(section, keyword)
2038 CALL keyword_release(keyword)
2039
2040 ! even if scaling parameter for exchange FRACTION (see above) is zero, the semi-empirical electron repulsion
2041 ! operator for exchange is not, so that a keyword is required to switch off sTDA exchange (if wanted)
2042 CALL keyword_create(keyword, __location__, name="DO_EXCHANGE", &
2043 description="Explicitly including or switching off sTDA exchange", &
2044 usage="DO_EXCHANGE", default_l_val=.true., lone_keyword_l_val=.true.)
2045 CALL section_add_keyword(section, keyword)
2046 CALL keyword_release(keyword)
2047
2048 CALL keyword_create(keyword, __location__, name="DO_EWALD", &
2049 description="Use Ewald type method for Coulomb interaction", &
2050 usage="DO_EWALD", default_l_val=.false., lone_keyword_l_val=.true.)
2051 CALL section_add_keyword(section, keyword)
2052 CALL keyword_release(keyword)
2053
2054 CALL keyword_create(keyword, __location__, name="EPS_TD_FILTER", &
2055 description="Threshold for filtering the transition density matrix", &
2056 usage="EPS_TD_FILTER epsf", default_r_val=1.e-10_dp)
2057 CALL section_add_keyword(section, keyword)
2058 CALL keyword_release(keyword)
2059
2060 CALL keyword_create(keyword, __location__, name="MATAGA_NISHIMOTO_CEXP", &
2061 description="Exponent used in Mataga-Nishimoto formula for Coulomb (alpha). "// &
2062 "Default value is method dependent!", &
2063 usage="MATAGA_NISHIMOTO_CEXP cexp", default_r_val=-99.0_dp)
2064 CALL section_add_keyword(section, keyword)
2065 CALL keyword_release(keyword)
2066
2067 CALL keyword_create(keyword, __location__, name="MATAGA_NISHIMOTO_XEXP", &
2068 description="Exponent used in Mataga-Nishimoto formula for Exchange (beta). "// &
2069 "Default value is method dependent!", &
2070 usage="MATAGA_NISHIMOTO_XEXP xexp", default_r_val=-99.0_dp)
2071 CALL section_add_keyword(section, keyword)
2072 CALL keyword_release(keyword)
2073
2074 CALL keyword_create(keyword, __location__, name="COULOMB_SR_CUT", &
2075 description="Maximum range of short range part of Coulomb interaction.", &
2076 usage="COULOMB_SR_CUT rcut", default_r_val=20.0_dp)
2077 CALL section_add_keyword(section, keyword)
2078 CALL keyword_release(keyword)
2079
2080 CALL keyword_create(keyword, __location__, name="COULOMB_SR_EPS", &
2081 description="Threshold for short range part of Coulomb interaction.", &
2082 usage="COULOMB_SR_EPS sreps", default_r_val=1.e-03_dp)
2083 CALL section_add_keyword(section, keyword)
2084 CALL keyword_release(keyword)
2085
2086 END SUBROUTINE create_stda_section
2087
2088! **************************************************************************************************
2089!> \brief creates an input section for electronic band structure calculations
2090!> \param section section to create
2091!> \par History
2092!> * 07.2023 created [Jan Wilhelm]
2093! **************************************************************************************************
2094 SUBROUTINE create_bandstructure_section(section)
2095 TYPE(section_type), POINTER :: section
2096
2097 TYPE(keyword_type), POINTER :: keyword
2098 TYPE(section_type), POINTER :: subsection
2099
2100 cpassert(.NOT. ASSOCIATED(section))
2101 CALL section_create(section, __location__, name="BANDSTRUCTURE", &
2102 description="Parameters needed to set up a calculation for "// &
2103 "electronic level energies of molecules and the electronic band "// &
2104 "structure of materials from post-SCF schemes (GW, perturbative "// &
2105 "spin-orbit coupling). Also, the density of states (DOS), "// &
2106 "projected density of states (PDOS), local density of states (LDOS), "// &
2107 "local valence band maximum (LVBM), local conduction band minimum "// &
2108 "(LCBM) and local band gap can be calculated. Please note that "// &
2109 "all methods in this section start from a Gamma-only DFT SCF. "// &
2110 "You need to make sure that the cell chosen in the DFT SCF is "// &
2111 "converged in the cell size. Band structures are computed "// &
2112 "for the primitive cell (i.e. the smallest possible unit cell of "// &
2113 "the input structure which is detected automatically). Moreover, "// &
2114 "spin-orbit coupling (SOC) on eigenvalues and band structures is "// &
2115 "available using Hartwigsen-Goedecker-Hutter "// &
2116 "pseudopotentials.", &
2117 n_keywords=1, n_subsections=1, repeats=.false.)
2118
2119 NULLIFY (keyword, subsection)
2120 CALL keyword_create(keyword, __location__, &
2121 name="_SECTION_PARAMETERS_", &
2122 description="Controls the activation of the band structure calculation.", &
2123 default_l_val=.false., &
2124 lone_keyword_l_val=.true.)
2125 CALL section_add_keyword(section, keyword)
2126 CALL keyword_release(keyword)
2127
2128 ! here we generate a subsection for getting a k-point path for the bandstructure
2129 CALL create_kpoint_set_section(subsection, "BANDSTRUCTURE_PATH")
2130 CALL section_add_subsection(section, subsection)
2131 CALL section_release(subsection)
2132
2133 CALL create_gw_section(subsection)
2134 CALL section_add_subsection(section, subsection)
2135 CALL section_release(subsection)
2136
2137 CALL create_soc_section(subsection)
2138 CALL section_add_subsection(section, subsection)
2139 CALL section_release(subsection)
2140
2141 CALL create_dos_section(subsection)
2142 CALL section_add_subsection(section, subsection)
2143 CALL section_release(subsection)
2144
2145 END SUBROUTINE create_bandstructure_section
2146
2147! **************************************************************************************************
2148!> \brief creates an input section for a GW calculation for the electronic band structure
2149!> \param section section to create
2150!> \par History
2151!> * 07.2023 created [Jan Wilhelm]
2152! **************************************************************************************************
2153 SUBROUTINE create_gw_section(section)
2154 TYPE(section_type), POINTER :: section
2155
2156 TYPE(keyword_type), POINTER :: keyword
2157 TYPE(section_type), POINTER :: print_key, subsection
2158
2159 cpassert(.NOT. ASSOCIATED(section))
2160 CALL section_create(section, __location__, name="GW", &
2161 description="Parameters needed to set up a GW calculation for "// &
2162 "electronic level energies $\varepsilon_{n\mathbf{k}}^{G_0W_0}$ "// &
2163 "of molecules and the band structure of materials: "// &
2164 "$\varepsilon_{n\mathbf{k}}^{G_0W_0}= "// &
2165 "\varepsilon_{n\mathbf{k}}^\text{DFT}+\Sigma_{n\mathbf{k}} "// &
2166 "-v^\text{xc}_{n\mathbf{k}}$. "// &
2167 "For the GW algorithm for molecules, see "// &
2168 "<https://doi.org/10.1021/acs.jctc.0c01282>. "// &
2169 "For 2D materials, see <https://doi.org/10.1021/acs.jctc.3c01230>.", &
2170 n_keywords=1, n_subsections=1, repeats=.false.)
2171
2172 NULLIFY (keyword)
2173 CALL keyword_create(keyword, __location__, &
2174 name="_SECTION_PARAMETERS_", &
2175 description="Controls the activation of the GW calculation.", &
2176 default_l_val=.false., &
2177 lone_keyword_l_val=.true.)
2178 CALL section_add_keyword(section, keyword)
2179 CALL keyword_release(keyword)
2180
2181 CALL keyword_create(keyword, __location__, name="NUM_TIME_FREQ_POINTS", &
2182 description="Number of discrete points for the imaginary-time "// &
2183 "grid and the imaginary-frequency grid. The more points, the more "// &
2184 "precise is the calculation. Typically, 10 points are good "// &
2185 "for 0.1 eV precision of band structures and molecular energy "// &
2186 "levels, 20 points for 0.03 eV precision, "// &
2187 "and 30 points for 0.01 eV precision, see Table I in "// &
2188 "<https://doi.org/10.1021/acs.jctc.0c01282>. GW computation time "// &
2189 "increases linearly with `NUM_TIME_FREQ_POINTS`.", &
2190 usage="NUM_TIME_FREQ_POINTS 30", &
2191 default_i_val=30)
2192 CALL section_add_keyword(section, keyword)
2193 CALL keyword_release(keyword)
2194
2195 CALL keyword_create(keyword, __location__, name="EPS_FILTER", &
2196 description="Determines a threshold for the DBCSR based sparse "// &
2197 "multiplications. Normally, `EPS_FILTER` determines accuracy "// &
2198 "and timing of low-scaling GW calculations. (Lower filter means "// &
2199 "higher numerical precision, but higher computational cost.)", &
2200 usage="EPS_FILTER 1.0E-6", &
2201 default_r_val=1.0e-8_dp)
2202 CALL section_add_keyword(section, keyword)
2203 CALL keyword_release(keyword)
2204
2205 CALL keyword_create(keyword, __location__, name="REGULARIZATION_MINIMAX", &
2206 description="Parameter to regularize the Fourier transformation with minimax grids. "// &
2207 "In case the parameter 0.0 is chosen, no regularization is performed.", &
2208 usage="REGULARIZATION_MINIMAX 1.0E-4", &
2209 default_r_val=-1.0_dp)
2210 CALL section_add_keyword(section, keyword)
2211 CALL keyword_release(keyword)
2212
2213 CALL keyword_create(keyword, __location__, name="REGULARIZATION_RI", &
2214 description="Parameter for RI regularization, setting a negative "// &
2215 "value triggers the default value. Affects RI basis set convergence "// &
2216 "but in any case large RI basis will give RI basis set convergence.", &
2217 usage="REGULARIZATION_RI 1.0E-4", &
2218 default_r_val=-1.0_dp)
2219 CALL section_add_keyword(section, keyword)
2220 CALL keyword_release(keyword)
2221
2222 CALL keyword_create(keyword, __location__, name="CUTOFF_RADIUS_RI", &
2223 description="The cutoff radius (in Angstrom) for the truncated "// &
2224 "Coulomb operator. The larger the cutoff radius, the faster "// &
2225 "converges the resolution of the identity (RI) with respect to the "// &
2226 "RI basis set size. Larger cutoff radius means higher computational "// &
2227 "cost.", &
2228 usage="CUTOFF_RADIUS_RI 3.0", &
2229 default_r_val=cp_unit_to_cp2k(value=3.0_dp, unit_str="angstrom"), &
2230 type_of_var=real_t, unit_str="angstrom")
2231 CALL section_add_keyword(section, keyword)
2232 CALL keyword_release(keyword)
2233
2234 CALL keyword_create(keyword, __location__, name="MEMORY_PER_PROC", &
2235 description="Specify the available memory per MPI process. Set "// &
2236 "`MEMORY_PER_PROC` as accurately as possible for good performance. If "// &
2237 "`MEMORY_PER_PROC` is set lower as the actually available "// &
2238 "memory per MPI process, the performance will be "// &
2239 "bad; if `MEMORY_PER_PROC` is set higher as the actually "// &
2240 "available memory per MPI process, the program might run out of "// &
2241 "memory. You can calculate `MEMORY_PER_PROC` as follows: "// &
2242 "Get the memory per node on your machine, mem_per_node "// &
2243 "(for example, from a supercomputer website, typically between "// &
2244 "100 GB and 2 TB), get the number of "// &
2245 "MPI processes per node, n_MPI_proc_per_node"// &
2246 " (for example from your run-script; if you "// &
2247 "use slurm, the number behind '--ntasks-per-node' is the number "// &
2248 "of MPI processes per node). Then calculate "// &
2249 "`MEMORY_PER_PROC` = mem_per_node / n_MPI_proc_per_node "// &
2250 "(typically between 2 GB and 50 GB). Unit of keyword: Gigabyte (GB).", &
2251 usage="MEMORY_PER_PROC 16", &
2252 default_r_val=2.0_dp)
2253 CALL section_add_keyword(section, keyword)
2254 CALL keyword_release(keyword)
2255
2256 CALL keyword_create(keyword, __location__, name="APPROX_KP_EXTRAPOL", &
2257 description="If true, use only a 4x4 kpoint mesh for frequency "// &
2258 "points $\omega_j, j \ge 2$ (instead of a 4x4 and 6x6 k-point mesh). "// &
2259 "The k-point extrapolation of $W_{PQ}(i\omega_j,\mathbf{q})$ "// &
2260 "is done approximately from $W_{PQ}(i\omega_1,\mathbf{q})$.", &
2261 usage="APPROX_KP_EXTRAPOL", &
2262 default_l_val=.false., lone_keyword_l_val=.true.)
2263 CALL section_add_keyword(section, keyword)
2264 CALL keyword_release(keyword)
2265
2266 CALL keyword_create(keyword, __location__, name="SIZE_LATTICE_SUM", &
2267 description="Parameter determines how many neighbor cells $\mathbf{R}$ "// &
2268 "are used for computing "// &
2269 "$V_{PQ}(\mathbf{k}) = "// &
2270 "\sum_{\mathbf{R}} e^{i\mathbf{k}\cdot\mathbf{R}}\,\langle P, "// &
2271 "\text{cell}{=}\mathbf{0}|1/r|Q,\text{cell}{=}\mathbf{R}\rangle$. "// &
2272 "Normally, parameter does not need to be touched.", &
2273 usage="SIZE_LATTICE_SUM 4", &
2274 default_i_val=3)
2275 CALL section_add_keyword(section, keyword)
2276 CALL keyword_release(keyword)
2277
2278 CALL keyword_create( &
2279 keyword, __location__, name="KPOINTS_W", &
2280 description="Monkhorst-Pack k-point mesh of size N_x, N_y, N_z for calculating "// &
2281 "$W_{PQ}^\mathbf{R}=\int_\text{BZ}\frac{d\mathbf{k}}{\Omega_\text{BZ}}\, "// &
2282 "e^{-i\mathbf{k}\cdot\mathbf{R}}\,W_{PQ}(\mathbf{k})$. "// &
2283 αα"For non-periodic directions , choose N_ = 1. "// &
2284 "Automatic choice of the k-point mesh for negative "// &
2285 "values, i.e. KPOINTS_W -1 -1 -1. "// &
2286 "K-point extrapolation of W is automatically switched on.", &
2287 usage="KPOINTS_W N_x N_y N_z", &
2288 n_var=3, type_of_var=integer_t, default_i_vals=(/-1, -1, -1/))
2289 CALL section_add_keyword(section, keyword)
2290 CALL keyword_release(keyword)
2291
2292 CALL keyword_create(keyword, __location__, name="HEDIN_SHIFT", &
2293 description="If true, use Hedin's shift in G0W0, evGW and evGW0. "// &
2294 "Details see in Li et al. JCTC 18, 7570 "// &
2295 "(2022), Figure 1. G0W0 with Hedin's shift should give "// &
2296 "similar GW eigenvalues as evGW0; at a lower "// &
2297 "computational cost.", &
2298 usage="HEDIN_SHIFT", &
2299 default_l_val=.false., &
2300 lone_keyword_l_val=.true.)
2301 CALL section_add_keyword(section, keyword)
2302 CALL keyword_release(keyword)
2303
2304 CALL keyword_create(keyword, __location__, name="FREQ_MAX_FIT", &
2305 description=Σω"For analytic continuation, a fit on (i) is performed. "// &
2306 Σω"This fit is then evaluated at a real frequency, (), which is used "// &
2307 "in the quasiparticle equation "// &
2308 "$\varepsilon_{n\mathbf{k}}^{G_0W_0}= "// &
2309 "\varepsilon_{n\mathbf{k}}^\text{DFT}+\Sigma_{n\mathbf{k}} "// &
2310 "-v^\text{xc}_{n\mathbf{k}}$. The keyword FREQ_MAX_FIT "// &
2311 Σω"determines fitting range for the self-energy (i) on "// &
2312 ωω"imaginary axis: i*[0, _max] for empty orbitals/bands, i*[-_max,0] "// &
2313 ω"for occ orbitals. A smaller _max might lead to better numerical "// &
2314 "stability (i.e., if you observe clearly wrong GW eigenvalues/bands "// &
2315 ω"around HOMO/LUMO, decreasing _max might fix this issue). "// &
2316 ω"A small benchmark of _max is contained in Fig. 5 of "// &
2317 "J. Wilhelm et al., JCTC 12, 3623-3635 (2016). "// &
2318 ω"Note that we used _max = 1 Ha = 27.211 eV in the benchmark "// &
2319 "M. Azizi et al., PRB 109, 245101 (2024).", &
2320 unit_str="eV", &
2321 usage="FREQ_MAX_FIT 20.0", &
2322 default_r_val=cp_unit_to_cp2k(value=10.0_dp, unit_str="eV"))
2323 CALL section_add_keyword(section, keyword)
2324 CALL keyword_release(keyword)
2325
2326 NULLIFY (subsection, print_key)
2327 CALL section_create(subsection, __location__, name="PRINT", &
2328 description="Printing of GW restarts.", &
2329 n_keywords=0, n_subsections=1, repeats=.false.)
2330 CALL cp_print_key_section_create(print_key, __location__, "RESTART", &
2331 description="Controls the printing of restart files "// &
2332 χΣ"for , W, .", &
2333 filename="", print_level=low_print_level, &
2334 common_iter_levels=3)
2335 CALL section_add_subsection(subsection, print_key)
2336 CALL section_release(print_key)
2337
2338 CALL section_add_subsection(section, subsection)
2339 CALL section_release(subsection)
2340
2341 END SUBROUTINE create_gw_section
2342
2343! **************************************************************************************************
2344!> \brief creates an input section for calculation SOC for the electronic band structure
2345!> \param section section to create
2346!> \par History
2347!> * 09.2023 created [Jan Wilhelm]
2348! **************************************************************************************************
2349 SUBROUTINE create_soc_section(section)
2350 TYPE(section_type), POINTER :: section
2351
2352 TYPE(keyword_type), POINTER :: keyword
2353
2354 cpassert(.NOT. ASSOCIATED(section))
2355 CALL section_create(section, __location__, name="SOC", &
2356 description="Switch on or off spin-orbit coupling. Use SOC "// &
2357 "parameters from non-local pseudopotentials as given in "// &
2358 "Hartwigsen, Goedecker, Hutter, Eq.(18), (19), "// &
2359 "<https://doi.org/10.1103/PhysRevB.58.3641>, "// &
2360 "$V_{\mu\nu}^{\mathrm{SOC}, (\alpha)} = "// &
2361 "(\hbar/2) \langle \phi_\mu | \sum_l \Delta "// &
2362 "V_l^\mathrm{SO}(\mathbf{r},\mathbf{r}') "// &
2363 "L^{(\alpha)} | \phi_\nu \rangle, "// &
2364 "\alpha = x, y, z$.", &
2365 n_keywords=1, n_subsections=1, repeats=.false.)
2366
2367 NULLIFY (keyword)
2368 CALL keyword_create(keyword, __location__, &
2369 name="_SECTION_PARAMETERS_", &
2370 description="Controls the activation of the SOC calculation.", &
2371 default_l_val=.false., &
2372 lone_keyword_l_val=.true.)
2373 CALL section_add_keyword(section, keyword)
2374 CALL keyword_release(keyword)
2375
2376 CALL keyword_create(keyword, __location__, name="ENERGY_WINDOW", &
2377 description="Apply SOC only for states with eigenvalues in the "// &
2378 "interval $[\varepsilon_\mathrm{VBM}-E_\mathrm{window}/2, "// &
2379 "\varepsilon_\mathrm{CBM}+E_\mathrm{window}/2]$. Might be necessary "// &
2380 "to use for large systems to prevent numerical instabilities.", &
2381 usage="ENERGY_WINDOW 5.0", &
2382 default_r_val=cp_unit_to_cp2k(value=40.0_dp, unit_str="eV"), &
2383 unit_str="eV")
2384 CALL section_add_keyword(section, keyword)
2385 CALL keyword_release(keyword)
2386
2387 END SUBROUTINE create_soc_section
2388
2389! **************************************************************************************************
2390!> \brief input section for computing the density of states and the projected density of states
2391!> \param section section to create
2392!> \par History
2393!> * 09.2023 created [Jan Wilhelm]
2394! **************************************************************************************************
2395 SUBROUTINE create_dos_section(section)
2396 TYPE(section_type), POINTER :: section
2397
2398 TYPE(keyword_type), POINTER :: keyword
2399 TYPE(section_type), POINTER :: subsection
2400
2401 cpassert(.NOT. ASSOCIATED(section))
2402 CALL section_create(section, __location__, name="DOS", &
2403 description="Parameters needed to calculate the density of states "// &
2404 "(DOS) and the projected density of states (PDOS).", &
2405 n_keywords=1, n_subsections=1, repeats=.false.)
2406
2407 NULLIFY (keyword)
2408 CALL keyword_create(keyword, __location__, &
2409 name="_SECTION_PARAMETERS_", &
2410 description="Controls the activation of the DOS calculation.", &
2411 default_l_val=.false., &
2412 lone_keyword_l_val=.true.)
2413 CALL section_add_keyword(section, keyword)
2414 CALL keyword_release(keyword)
2415
2416 CALL keyword_create(keyword, __location__, name="ENERGY_WINDOW", &
2417 description="Print DOS and PDOS in the energy window "// &
2418 "$[\varepsilon_\mathrm{VBM}-E_\mathrm{window}/2, "// &
2419 "\varepsilon_\mathrm{CBM}+E_\mathrm{window}/2]$,"// &
2420 " where VBM is the valence "// &
2421 "band maximum (or highest occupied molecular orbital, HOMO, for "// &
2422 "molecules) and CBM the conduction band minimum (or lowest "// &
2423 "unoccupied molecular orbital, LUMO, for molecules).", &
2424 usage="ENERGY_WINDOW 5.0", &
2425 default_r_val=cp_unit_to_cp2k(value=10.0_dp, unit_str="eV"), &
2426 unit_str="eV")
2427 CALL section_add_keyword(section, keyword)
2428 CALL keyword_release(keyword)
2429
2430 CALL keyword_create(keyword, __location__, name="ENERGY_STEP", &
2431 description="Resolution of the energy E when computing the $\rho(E)$.", &
2432 usage="ENERGY_STEP 0.01", &
2433 default_r_val=cp_unit_to_cp2k(value=0.01_dp, unit_str="eV"), &
2434 unit_str="eV")
2435 CALL section_add_keyword(section, keyword)
2436 CALL keyword_release(keyword)
2437
2438 CALL keyword_create(keyword, __location__, name="BROADENING", &
2439 description=α"Broadening in Gaussians used in the DOS; "// &
2440 "$\rho(E) = \sum_n \exp(((E-\varepsilon_n)/\alpha)^2)/("// &
2441 " \sqrt{2\pi} \alpha)$.", &
2442 usage="BROADENING 0.01", &
2443 default_r_val=cp_unit_to_cp2k(value=0.01_dp, unit_str="eV"), &
2444 unit_str="eV")
2445 CALL section_add_keyword(section, keyword)
2446 CALL keyword_release(keyword)
2447
2448 CALL keyword_create( &
2449 keyword, __location__, name="KPOINTS", &
2450 description="Monkhorst-Pack k-point mesh of size N_x, N_y, N_z for calculating "// &
2451 "the density of states (DOS). In GW, the KPOINT_DOS mesh is thus used as k-point "// &
2452 αα"mesh for the self-energy. For non-periodic directions , choose N_ = 1. "// &
2453 "Automatic choice of the k-point mesh for negative "// &
2454 α"values, i.e. KPOINTS_DOS -1 -1 -1 (automatic choice: N_ = 1 in non-periodic "// &
2455 "direction, 8 k-points in periodic direction). If you like to compute a "// &
2456 "band structure along a k-path, you can specify the k-path in "// &
2457 "&KPOINT_SET.", &
2458 usage="KPOINTS N_x N_y N_z", &
2459 n_var=3, type_of_var=integer_t, default_i_vals=(/-1, -1, -1/))
2460 CALL section_add_keyword(section, keyword)
2461 CALL keyword_release(keyword)
2462
2463 NULLIFY (subsection)
2464 CALL create_ldos_section(subsection)
2465 CALL section_add_subsection(section, subsection)
2466 CALL section_release(subsection)
2467
2468 END SUBROUTINE create_dos_section
2469
2470! **************************************************************************************************
2471!> \brief ...
2472!> \param section ...
2473! **************************************************************************************************
2474 SUBROUTINE create_ldos_section(section)
2475 TYPE(section_type), POINTER :: section
2476
2477 TYPE(keyword_type), POINTER :: keyword
2478
2479 cpassert(.NOT. ASSOCIATED(section))
2480 CALL section_create(section, __location__, name="LDOS", &
2481 description="Parameters needed to calculate the local density "// &
2482 "of states (LDOS). "// &
2483 "The LDOS is computed as $\rho(\mathbf{r},E) = "// &
2484 "\sum\limits_{n,\mathbf{k}}"// &
2485 " |\psi_{n\mathbf{k}}(r)|^2\, w_\mathbf{k}\, g(E-\varepsilon_{n\mathbf{k}})$ "// &
2486 "using the Gaussian weight function "// &
2487 "$g(x) = \exp(x^2/\alpha^2)/(\sqrt{2\pi}\alpha)$, $\alpha$ is the broadening "// &
2488 "from the &DOS section, and the k-point weight "// &
2489 "$w_\mathbf{k}$. The k-mesh is taken from the &DOS section.", &
2490 n_keywords=2, repeats=.false.)
2491
2492 NULLIFY (keyword)
2493 CALL keyword_create(keyword, __location__, &
2494 name="_SECTION_PARAMETERS_", &
2495 description="Activates the local VBM CBM gap calculation.", &
2496 default_l_val=.false., &
2497 lone_keyword_l_val=.true.)
2498 CALL section_add_keyword(section, keyword)
2499 CALL keyword_release(keyword)
2500
2501 CALL keyword_create(keyword, __location__, name="INTEGRATION", &
2502 description="Defines whether the LDOS is integrated along a "// &
2503 "coordinate. As an example, for INTEGRATION Z, the LDOS "// &
2504 "$\rho(x,y,E) = \int dz\, \rho(x,y,z,E)$ is computed.", &
2505 usage="INTEGRATION Z", &
2506 enum_c_vals=s2a("X", "Y", "Z", "NONE"), &
2507 enum_i_vals=(/int_ldos_x, int_ldos_y, int_ldos_z, int_ldos_none/), &
2508 enum_desc=s2a("Integrate over x coordinate (not yet implemented).", &
2509 "Integrate over y coordinate (not yet implemented).", &
2510 "Integrate over z coordinate.", &
2511 "No integration, print cube file as function "// &
2512 "of x,y,z (not yet implemented)."), &
2513 default_i_val=int_ldos_z)
2514 CALL section_add_keyword(section, keyword)
2515 CALL keyword_release(keyword)
2516
2517 CALL keyword_create( &
2518 keyword, __location__, name="BIN_MESH", &
2519 description="Mesh of size n x m for binning the space coordinates x and y of "// &
2520 "the LDOS $\rho(x,y,E)$. If -1, no binning is performed and the "// &
2521 "fine x, y resolution of the electron density from SCF is used.", &
2522 usage="BIN_MESH n m", &
2523 n_var=2, type_of_var=integer_t, default_i_vals=(/10, 10/))
2524 CALL section_add_keyword(section, keyword)
2525 CALL keyword_release(keyword)
2526
2527 END SUBROUTINE create_ldos_section
2528
2529! **************************************************************************************************
2530!> \brief creates an input section for a tip scan calculation
2531!> \param section section to create
2532!> \par History
2533!> * 04.2021 created [JGH]
2534! **************************************************************************************************
2535 SUBROUTINE create_tipscan_section(section)
2536 TYPE(section_type), POINTER :: section
2537
2538 TYPE(keyword_type), POINTER :: keyword
2539
2540 cpassert(.NOT. ASSOCIATED(section))
2541 CALL section_create(section, __location__, name="TIP_SCAN", &
2542 description="Parameters needed to set up a Tip Scan. "// &
2543 "Needs external definition of tip induced field.", &
2544 n_keywords=1, n_subsections=1, repeats=.false.)
2545
2546 NULLIFY (keyword)
2547
2548 CALL keyword_create(keyword, __location__, &
2549 name="_SECTION_PARAMETERS_", &
2550 description="Controls the activation of the Tip Scan procedure", &
2551 default_l_val=.false., &
2552 lone_keyword_l_val=.true.)
2553 CALL section_add_keyword(section, keyword)
2554 CALL keyword_release(keyword)
2555
2556 CALL keyword_create(keyword, __location__, name="SCAN_DIRECTION", &
2557 description="Defines scan direction and scan type(line, plane).", &
2558 usage="SCAN_DIRECTION XY", &
2559 enum_c_vals=s2a("X", "Y", "Z", "XY", "XZ", "YZ", "XYZ"), &
2560 enum_i_vals=(/scan_x, scan_y, scan_z, scan_xy, scan_xz, scan_yz, scan_xyz/), &
2561 default_i_val=scan_xy)
2562 CALL section_add_keyword(section, keyword)
2563 CALL keyword_release(keyword)
2564
2565 CALL keyword_create(keyword, __location__, name="REFERENCE_POINT", &
2566 description="The reference point to define the absolute position of the scan. ", &
2567 usage="REFERENCE_POINT 0.0 0.0 1.0", &
2568 n_var=3, default_r_vals=(/0.0_dp, 0.0_dp, 0.0_dp/), type_of_var=real_t, &
2569 unit_str="angstrom")
2570 CALL section_add_keyword(section, keyword)
2571 CALL keyword_release(keyword)
2572
2573 CALL keyword_create(keyword, __location__, name="SCAN_POINTS", &
2574 description="Number of points calculated for each scan direction.", &
2575 usage="SCAN_POINTS 20 20", &
2576 n_var=-1, type_of_var=integer_t)
2577 CALL section_add_keyword(section, keyword)
2578 CALL keyword_release(keyword)
2579
2580 CALL keyword_create(keyword, __location__, name="SCAN_STEP", &
2581 description="Step size for each scan direction.", &
2582 usage="SCAN_STEP 0.01 0.01", &
2583 n_var=-1, type_of_var=real_t, unit_str="angstrom")
2584 CALL section_add_keyword(section, keyword)
2585 CALL keyword_release(keyword)
2586
2587 CALL keyword_create(keyword, __location__, name="TIP_FILENAME", &
2588 description="Filename of tip potential defined in cube file format.", &
2589 usage="TIP_FILENAME <filename>", &
2590 type_of_var=lchar_t)
2591 CALL section_add_keyword(section, keyword)
2592 CALL keyword_release(keyword)
2593
2594 END SUBROUTINE create_tipscan_section
2595
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public putrino2000
integer, save, public weber2009
integer, save, public kondov2007
integer, save, public luber2014
integer, save, public iannuzzi2005
integer, save, public sebastiani2001
integer, save, public putrino2002
integer, save, public vazdacruz2021
integer, save, public futera2017
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 low_print_level
integer, parameter, public medium_print_level
integer, parameter, public high_print_level
integer, parameter, public add_last_numeric
integer, parameter, public silent_print_level
subroutine, public cp_print_key_section_create(print_key_section, location, name, description, print_level, each_iter_names, each_iter_values, add_last, filename, common_iter_levels, citations, unit_str)
creates a print_key section
unit conversion facility
Definition cp_units.F:30
real(kind=dp) function, public cp_unit_to_cp2k(value, unit_str, defaults, power)
converts to the internal cp2k units to the given unit
Definition cp_units.F:1149
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public int_ldos_y
integer, parameter, public use_mom_ref_coac
integer, parameter, public tddfpt_dipole_berry
integer, parameter, public oe_saop
integer, parameter, public do_no_et
integer, parameter, public current_orb_center_wannier
integer, parameter, public scan_x
integer, parameter, public scan_xyz
integer, parameter, public use_mom_ref_user
integer, parameter, public tddfpt_sf_col
integer, parameter, public tddfpt_kernel_none
integer, parameter, public use_mom_ref_com
integer, parameter, public gto_cartesian
integer, parameter, public gto_spherical
integer, parameter, public current_gauge_atom
integer, parameter, public current_gauge_r
integer, parameter, public oe_none
integer, parameter, public ot_precond_full_kinetic
integer, parameter, public current_gauge_r_and_step_func
integer, parameter, public oe_shift
integer, parameter, public int_ldos_x
integer, parameter, public current_orb_center_box
integer, parameter, public tddfpt_dipole_velocity
integer, parameter, public current_orb_center_common
integer, parameter, public ot_precond_full_single
integer, parameter, public scan_xy
integer, parameter, public scan_xz
integer, parameter, public tddfpt_kernel_full
integer, parameter, public ot_precond_none
integer, parameter, public scan_y
integer, parameter, public int_ldos_z
integer, parameter, public ot_precond_full_single_inverse
integer, parameter, public current_orb_center_atom
integer, parameter, public scan_z
integer, parameter, public tddfpt_sf_noncol
integer, parameter, public do_spin_density
integer, parameter, public tddfpt_dipole_length
integer, parameter, public use_mom_ref_zero
integer, parameter, public oe_lb
integer, parameter, public tddfpt_kernel_stda
integer, parameter, public int_ldos_none
integer, parameter, public no_sf_tddfpt
integer, parameter, public do_et_ddapc
integer, parameter, public ot_precond_s_inverse
integer, parameter, public do_full_density
integer, parameter, public scan_yz
integer, parameter, public oe_gllb
integer, parameter, public ot_precond_full_all
input section for atomic properties
subroutine, public create_atprop_section(section)
Creates the ATOMIC section.
function that build the dft section of the input
subroutine, public create_mgrid_section(section, create_subsections)
creates the multigrid
subroutine, public create_interp_section(section)
creates the interpolation section
function that build the kpoints section of the input
subroutine, public create_kpoint_set_section(section, section_name)
...
subroutine, public create_localize_section(section)
parameters fo the localization of wavefunctions
function that build the dft section of the input
subroutine, public create_properties_section(section)
Create the PROPERTIES section.
function that build the QS section of the input
subroutine, public create_lrigpw_section(section)
input section for optional parameters for LRIGPW LRI: local resolution of identity
subroutine, public create_ddapc_restraint_section(section, section_name)
...
function that builds the resp section of the input
subroutine, public create_resp_section(section)
Creates the RESP section.
function that build the XAS section of the input
subroutine, public create_xas_tdp_section(section)
makes the input section for core-level spectroscopy simulations using linear response TDDFT
function that build the xc section of the input
subroutine, public create_xc_section(section)
creates the input section for the xc part
represents keywords in an input
subroutine, public keyword_release(keyword)
releases the given keyword (see doc/ReferenceCounting.html)
subroutine, public keyword_create(keyword, location, name, description, usage, type_of_var, n_var, repeats, variants, default_val, default_l_val, default_r_val, default_lc_val, default_c_val, default_i_val, default_l_vals, default_r_vals, default_c_vals, default_i_vals, lone_keyword_val, lone_keyword_l_val, lone_keyword_r_val, lone_keyword_c_val, lone_keyword_i_val, lone_keyword_l_vals, lone_keyword_r_vals, lone_keyword_c_vals, lone_keyword_i_vals, enum_c_vals, enum_i_vals, enum, enum_strict, enum_desc, unit_str, citations, deprecation_notice, removed)
creates a keyword object
objects that represent the structure of input sections and the data contained in an input section
subroutine, public section_create(section, location, name, description, n_keywords, n_subsections, repeats, citations, deprecation_notice)
creates a list of keywords
subroutine, public section_add_keyword(section, keyword)
adds a keyword to the given section
subroutine, public section_add_subsection(section, subsection)
adds a subsection to the given section
recursive subroutine, public section_release(section)
releases the given keyword list (see doc/ReferenceCounting.html)
a wrapper for basic fortran types.
integer, parameter, public real_t
integer, parameter, public lchar_t
integer, parameter, public logical_t
integer, parameter, public char_t
integer, parameter, public integer_t
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
Utilities for string manipulations.
represent a keyword in the input
represent a section of the input file