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