(git:d43dca4)
Loading...
Searching...
No Matches
qs_dispersion_utils.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2026 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Set disperson types for DFT calculations
10!> \author JGH (04.2014)
11! **************************************************************************************************
13
20 USE input_constants, ONLY: &
26 USE kinds, ONLY: default_string_length,&
27 dp
28 USE physcon, ONLY: bohr,&
29 kjmol
34 USE qs_kind_types, ONLY: get_qs_kind,&
37#include "./base/base_uses.f90"
38
39 IMPLICIT NONE
40
41 PRIVATE
42
43 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_dispersion_utils'
44
46 PUBLIC :: cellhash
47
48! **************************************************************************************************
49CONTAINS
50! **************************************************************************************************
51!> \brief ...
52!> \param dispersion_env ...
53!> \param xc_section ...
54! **************************************************************************************************
55 SUBROUTINE qs_dispersion_env_set(dispersion_env, xc_section)
56 TYPE(qs_dispersion_type), POINTER :: dispersion_env
57 TYPE(section_vals_type), POINTER :: xc_section
58
59 LOGICAL :: exfun, explicit
60 REAL(dp), POINTER :: params(:), scal(:)
61 TYPE(section_vals_type), POINTER :: nl_section, pp_section, vdw_section, &
62 xc_fun_section
63
64 cpassert(ASSOCIATED(dispersion_env))
65
66 ! set general defaults
67 dispersion_env%doabc = .false.
68 dispersion_env%c9cnst = .false.
69 dispersion_env%lrc = .false.
70 dispersion_env%srb = .false.
71 dispersion_env%verbose = .false.
72 dispersion_env%nd3_exclude_pair = 0
73 NULLIFY (dispersion_env%c6ab, dispersion_env%maxci, dispersion_env%r0ab, dispersion_env%rcov, &
74 dispersion_env%r2r4, dispersion_env%cn, dispersion_env%cnkind, dispersion_env%cnlist, &
75 dispersion_env%d3_exclude_pair)
76 NULLIFY (dispersion_env%q_mesh, dispersion_env%kernel, dispersion_env%d2phi_dk2, &
77 dispersion_env%d2y_dx2)
78 NULLIFY (dispersion_env%sab_vdw, dispersion_env%sab_cn)
79 NULLIFY (dispersion_env%dftd_section)
80 NULLIFY (vdw_section, xc_fun_section)
81 vdw_section => section_vals_get_subs_vals(xc_section, "vdw_potential")
82 xc_fun_section => section_vals_get_subs_vals(xc_section, "XC_FUNCTIONAL")
83 CALL section_vals_val_get(vdw_section, "POTENTIAL_TYPE", i_val=dispersion_env%type)
84 IF (dispersion_env%type == xc_vdw_fun_pairpot) THEN
85 NULLIFY (pp_section)
86 pp_section => section_vals_get_subs_vals(vdw_section, "PAIR_POTENTIAL")
87 CALL section_vals_val_get(pp_section, "VERBOSE_OUTPUT", l_val=dispersion_env%verbose)
88 CALL section_vals_val_get(pp_section, "TYPE", i_val=dispersion_env%pp_type)
89 IF (dispersion_env%pp_type == vdw_pairpot_dftd2) THEN
90 ! functional parameters for Grimme D2 type
91 CALL section_vals_val_get(pp_section, "EXP_PRE", r_val=dispersion_env%exp_pre)
92 CALL section_vals_val_get(pp_section, "SCALING", explicit=explicit)
93 IF (.NOT. explicit) THEN
94 CALL section_vals_val_get(pp_section, "REFERENCE_FUNCTIONAL", explicit=exfun)
95 IF (.NOT. exfun) THEN
96 CALL cp_abort(__location__, "D2 vdW REFERENCE_FUNCTIONAL expected a second parameter "// &
97 "(example REFERENCE_FUNCTIONAL PBE). "// &
98 "Go to https://www.chemie.uni-bonn.de/grimme/de/software/dft-d3/ "// &
99 "for a full list of supported functionals")
100 END IF
101 CALL qs_scaling_dftd2(dispersion_env%scaling, vdw_section)
102 ELSE
103 CALL section_vals_val_get(pp_section, "SCALING", r_val=dispersion_env%scaling)
104 END IF
105 ELSE
106 dispersion_env%exp_pre = 0._dp
107 dispersion_env%scaling = 0._dp
108 END IF
109 IF (dispersion_env%pp_type == vdw_pairpot_dftd3 .OR. &
110 dispersion_env%pp_type == vdw_pairpot_dftd3bj) THEN
111 ! functional parameters for Grimme DFT-D3 type
112 CALL section_vals_val_get(pp_section, "EPS_CN", r_val=dispersion_env%eps_cn)
113 CALL section_vals_val_get(pp_section, "CALCULATE_C9_TERM", l_val=dispersion_env%doabc)
114 CALL section_vals_val_get(pp_section, "REFERENCE_C9_TERM", l_val=dispersion_env%c9cnst)
115 CALL section_vals_val_get(pp_section, "LONG_RANGE_CORRECTION", l_val=dispersion_env%lrc)
116 CALL section_vals_val_get(pp_section, "SHORT_RANGE_CORRECTION", l_val=dispersion_env%srb)
117 CALL section_vals_val_get(pp_section, "SHORT_RANGE_CORRECTION_PARAMETERS", r_vals=params)
118 dispersion_env%srb_params(1:4) = params(1:4)
119 ! KG corrections
120 CALL section_vals_val_get(pp_section, "MOLECULE_CORRECTION", l_val=dispersion_env%domol)
121 CALL section_vals_val_get(pp_section, "MOLECULE_CORRECTION_C8", r_val=dispersion_env%kgc8)
122 IF (dispersion_env%pp_type == vdw_pairpot_dftd3) THEN
123 CALL section_vals_val_get(pp_section, "D3_SCALING", explicit=explicit)
124 ELSE IF (dispersion_env%pp_type == vdw_pairpot_dftd3bj) THEN
125 CALL section_vals_val_get(pp_section, "D3BJ_SCALING", explicit=explicit)
126 END IF
127 IF (.NOT. explicit) THEN
128 CALL section_vals_val_get(pp_section, "REFERENCE_FUNCTIONAL", explicit=exfun)
129 IF (.NOT. exfun) THEN
130 CALL cp_abort(__location__, "D3 vdW REFERENCE_FUNCTIONAL expected a second parameter "// &
131 "(example REFERENCE_FUNCTIONAL PBE). "// &
132 "Go to https://www.chemie.uni-bonn.de/grimme/de/software/dft-d3/ for a full "// &
133 "list of supported functionals")
134 ELSE
135 CALL section_vals_val_get(vdw_section, &
136 "PAIR_POTENTIAL%REFERENCE_FUNCTIONAL", &
137 c_val=dispersion_env%ref_functional)
138 END IF
139 IF (dispersion_env%pp_type == vdw_pairpot_dftd3) THEN
140 CALL qs_scaling_dftd3(dispersion_env%s6, dispersion_env%sr6, &
141 dispersion_env%s8, vdw_section)
142 ELSE IF (dispersion_env%pp_type == vdw_pairpot_dftd3bj) THEN
143 CALL qs_scaling_dftd3bj(dispersion_env%s6, dispersion_env%a1, dispersion_env%s8, &
144 dispersion_env%a2, vdw_section)
145 END IF
146 ELSE
147 IF (dispersion_env%pp_type == vdw_pairpot_dftd3) THEN
148 ! zero damping
149 CALL section_vals_val_get(pp_section, "D3_SCALING", r_vals=scal)
150 dispersion_env%s6 = scal(1)
151 dispersion_env%sr6 = scal(2)
152 dispersion_env%s8 = scal(3)
153 dispersion_env%a1 = 0.0_dp
154 dispersion_env%a2 = 0.0_dp
155 ELSE IF (dispersion_env%pp_type == vdw_pairpot_dftd3bj) THEN
156 ! BJ damping
157 CALL section_vals_val_get(pp_section, "D3BJ_SCALING", r_vals=scal)
158 dispersion_env%s6 = scal(1)
159 dispersion_env%a1 = scal(2)
160 dispersion_env%s8 = scal(3)
161 dispersion_env%a2 = scal(4)
162 dispersion_env%sr6 = 0.0_dp
163 END IF
164 END IF
165 ELSE
166 dispersion_env%s6 = 0._dp
167 dispersion_env%sr6 = 0._dp
168 dispersion_env%s8 = 0._dp
169 dispersion_env%s9 = 0._dp
170 dispersion_env%a1 = 0._dp
171 dispersion_env%a2 = 0._dp
172 dispersion_env%eps_cn = 0._dp
173 END IF
174 IF (dispersion_env%pp_type == vdw_pairpot_dftd4) THEN
175 CALL section_vals_val_get(pp_section, "D4_SCALING", explicit=explicit)
176 IF (.NOT. explicit) THEN
177 CALL section_vals_val_get(pp_section, "REFERENCE_FUNCTIONAL", explicit=exfun)
178 IF (.NOT. exfun) THEN
179 cpabort("Missing REFERENCE_FUNCTIONAL or D4_SCALING for D4")
180 ELSE
181 CALL section_vals_val_get(vdw_section, &
182 "PAIR_POTENTIAL%REFERENCE_FUNCTIONAL", &
183 c_val=dispersion_env%ref_functional)
184 END IF
185 ELSE
186 CALL section_vals_val_get(pp_section, "D4_SCALING", r_vals=scal)
187 dispersion_env%s6 = scal(1)
188 dispersion_env%a1 = scal(2)
189 dispersion_env%s8 = scal(3)
190 dispersion_env%a2 = scal(4)
191 dispersion_env%sr6 = 0.0_dp
192 dispersion_env%ref_functional = "none"
193 END IF
194 IF (trim(adjustl(dispersion_env%ref_functional)) == "none") THEN
195 cpabort("Missing REFERENCE_FUNCTIONAL or D4_SCALING for D4")
196 END IF
197 CALL section_vals_val_get(pp_section, "EPS_CN", r_val=dispersion_env%eps_cn)
198 CALL section_vals_val_get(pp_section, "D4_REFERENCE_CODE", &
199 l_val=dispersion_env%d4_reference_code)
200 CALL section_vals_val_get(pp_section, "D4_DEBUG", l_val=dispersion_env%d4_debug)
201 CALL section_vals_val_get(pp_section, "D4_CUTOFF", r_val=dispersion_env%rc_d4)
202 CALL section_vals_val_get(pp_section, "D4_CN_CUTOFF", r_val=dispersion_env%rc_cn)
203 CALL section_vals_val_get(pp_section, "FACTOR_S9_TERM", r_val=dispersion_env%s9)
204 !C9 term default=T for D4
205 CALL section_vals_val_get(pp_section, "CALCULATE_C9_TERM", explicit=exfun)
206 IF (exfun) THEN
207 CALL section_vals_val_get(pp_section, "CALCULATE_C9_TERM", l_val=dispersion_env%doabc)
208 ELSE
209 dispersion_env%doabc = .true.
210 END IF
211 END IF
212 CALL section_vals_val_get(pp_section, "R_CUTOFF", r_val=dispersion_env%rc_disp)
213 CALL section_vals_val_get(pp_section, "PARAMETER_FILE_NAME", &
214 c_val=dispersion_env%parameter_file_name)
215 ! set DFTD section for output handling
216 dispersion_env%dftd_section => pp_section
217 ELSE IF (dispersion_env%type == xc_vdw_fun_nonloc) THEN
218 NULLIFY (nl_section)
219 nl_section => section_vals_get_subs_vals(vdw_section, "NON_LOCAL")
220 CALL section_vals_val_get(nl_section, "VERBOSE_OUTPUT", l_val=dispersion_env%verbose)
221 CALL section_vals_val_get(nl_section, "KERNEL_FILE_NAME", &
222 c_val=dispersion_env%kernel_file_name)
223 CALL section_vals_val_get(nl_section, "TYPE", i_val=dispersion_env%nl_type)
224 CALL section_vals_val_get(nl_section, "CUTOFF", r_val=dispersion_env%pw_cutoff)
225 CALL section_vals_val_get(nl_section, "PARAMETERS", r_vals=params)
226 CALL section_vals_val_get(nl_section, "SCALE", r_val=dispersion_env%scale_rvv10)
227 dispersion_env%b_value = params(1)
228 dispersion_env%c_value = params(2)
229 END IF
230 END SUBROUTINE qs_dispersion_env_set
231
232! **************************************************************************************************
233!> \brief ...
234!> \param qs_env ...
235!> \param dispersion_env ...
236!> \param ounit ...
237! **************************************************************************************************
238 SUBROUTINE qs_write_dispersion(qs_env, dispersion_env, ounit)
239 TYPE(qs_environment_type), POINTER :: qs_env
240 TYPE(qs_dispersion_type), POINTER :: dispersion_env
241 INTEGER, INTENT(in), OPTIONAL :: ounit
242
243 CHARACTER(LEN=2) :: symbol
244 INTEGER :: i, ikind, nkind, output_unit
245 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
246 TYPE(cp_logger_type), POINTER :: logger
247 TYPE(qs_atom_dispersion_type), POINTER :: disp
248 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
249 TYPE(section_vals_type), POINTER :: dft_section
250
251 IF (PRESENT(ounit)) THEN
252 output_unit = ounit
253 ELSE
254 NULLIFY (logger)
255 logger => cp_get_default_logger()
256
257 dft_section => section_vals_get_subs_vals(qs_env%input, "DFT")
258 output_unit = cp_print_key_unit_nr(logger, dft_section, &
259 "PRINT%DFT_CONTROL_PARAMETERS", extension=".Log")
260 END IF
261
262 IF (output_unit > 0) THEN
263 ! vdW type specific output
264 IF (dispersion_env%type == xc_vdw_fun_pairpot) THEN
265 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T67,'Pair Potential')")
266 ! Pair potentials
267 IF (dispersion_env%pp_type == vdw_pairpot_dftd2) THEN
268 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T35,'DFT-D2')")
269 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T35,'Potential Form: S. Grimme, JCC 27: 1787 (2006)')")
270 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T35,'Cutoff Radius [Bohr]:',T73,F8.2)") dispersion_env%rc_disp
271 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T35,'Scaling Factor:',T73,F8.4)") dispersion_env%scaling
272 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T35,'Exp Prefactor for Damping:',T73,F8.1)") dispersion_env%exp_pre
273 CALL get_qs_env(qs_env, atomic_kind_set=atomic_kind_set, qs_kind_set=qs_kind_set)
274 nkind = SIZE(atomic_kind_set)
275 DO ikind = 1, nkind
276 CALL get_atomic_kind(atomic_kind_set(ikind), element_symbol=symbol)
277 CALL get_qs_kind(qs_kind_set(ikind), dispersion=disp)
278 IF (disp%defined) THEN
279 WRITE (output_unit, fmt="(' vdW PARAMETER| ',T18,'Atom=',A2, "// &
280 "T28,'C6[J*nm^6*mol^-1]=',F8.4,T63,'r(vdW)[A]=',F8.4)") &
281 symbol, disp%c6/(1000._dp*bohr**6/kjmol), disp%vdw_radii/bohr
282 ELSE
283 WRITE (output_unit, fmt="(' vdW PARAMETER| ',T20,'Atom=',A2,T70,'not defined')")
284 END IF
285 END DO
286 ELSE IF (dispersion_env%pp_type == vdw_pairpot_dftd3) THEN
287 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'DFT-D3 (Version 3.1)')")
288 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'Potential Form: S. Grimme et al, JCP 132: 154104 (2010)')")
289 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'Zero Damping')")
290 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'Cutoff Radius [Bohr]:',T73,F8.2)") dispersion_env%rc_disp
291 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'s6 Scaling Factor:',T73,F8.4)") dispersion_env%s6
292 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'sr6 Scaling Factor:',T73,F8.4)") dispersion_env%sr6
293 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'s8 Scaling Factor:',T73,F8.4)") dispersion_env%s8
294 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'Cutoff for CN calculation:',T69,E12.4)") dispersion_env%eps_cn
295 IF (dispersion_env%nd3_exclude_pair > 0) THEN
296 DO i = 1, dispersion_env%nd3_exclude_pair
297 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'Excluded Pairs: ',T76,I2,' ',I2)") &
298 dispersion_env%d3_exclude_pair(i, :)
299 END DO
300 END IF
301 ELSE IF (dispersion_env%pp_type == vdw_pairpot_dftd3bj) THEN
302 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'DFT-D3 (Version 3.1)')")
303 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'Potential Form: S. Grimme et al, JCP 132: 154104 (2010)')")
304 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'BJ Damping: S. Grimme et al, JCC 32: 1456 (2011)')")
305 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'Cutoff Radius [Bohr]:',T73,F8.2)") dispersion_env%rc_disp
306 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'s6 Scaling Factor:',T73,F8.4)") dispersion_env%s6
307 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'a1 Damping Factor:',T73,F8.4)") dispersion_env%a1
308 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'s8 Scaling Factor:',T73,F8.4)") dispersion_env%s8
309 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'a2 Damping Factor:',T73,F8.4)") dispersion_env%a2
310 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'Cutoff for CN calculation:',T69,E12.4)") dispersion_env%eps_cn
311 IF (dispersion_env%nd3_exclude_pair > 0) THEN
312 DO i = 1, dispersion_env%nd3_exclude_pair
313 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'Excluded Kind Pairs: ',T76,I2,' ',I2)") &
314 dispersion_env%d3_exclude_pair(i, :)
315 END DO
316 END IF
317 ELSE IF (dispersion_env%pp_type == vdw_pairpot_dftd4) THEN
318#if defined(__DFTD4_V3)
319 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'DFT-D4(Version 3.7)')")
320#else
321 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'DFT-D4(Version 4.0)')")
322#endif
323 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'see https://github.com/dftd4/dftd4')")
324 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'E. Caldeweyher et al, PCCP 22: 8499 (2020)')")
325 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'E. Caldeweyher et al, JCP 150: 154122 (2019)')")
326 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T26,'E. Caldeweyher et al, JCP 147: 034112 (2017)')")
327 END IF
328 ELSE IF (dispersion_env%type == xc_vdw_fun_nonloc) THEN
329 WRITE (output_unit, fmt="(' vdW POTENTIAL| ',T61,'Non-local Functional')")
330 WRITE (output_unit, &
331 fmt="(' vdW POTENTIAL| ','Implementation: G. Roman-Perez, J. Soler, PRL 103: 096102 (2009)')")
332 WRITE (output_unit, &
333 fmt="(' vdW POTENTIAL| ',T38,' T. Thonhauser et al, PRB 76: 125112 (2007)')")
334 WRITE (output_unit, &
335 fmt="(' vdW POTENTIAL| ',T22,' R. Sabatini et al, J.Phys:Condens Matter 24: 424209 (2012)')")
336 WRITE (output_unit, &
337 fmt="(' vdW POTENTIAL| ',T16,' Based on QE implementation by Brian Kolb, Timo Thonhauser (2009)')")
338 SELECT CASE (dispersion_env%nl_type)
339 CASE DEFAULT
340 ! unknown functional
341 cpabort("")
342 CASE (vdw_nl_drsll)
343 WRITE (output_unit, &
344 fmt="(' vdW POTENTIAL| ','DRSLL Functional: M. Dion et al, PRL 92: 246401 (2004)')")
345 CASE (vdw_nl_lmkll)
346 WRITE (output_unit, &
347 fmt="(' vdW POTENTIAL| ','LMKLL Functional: K. Lee et al, PRB 82: 081101 (2010)')")
348 CASE (vdw_nl_rvv10)
349 WRITE (output_unit, &
350 fmt="(' vdW POTENTIAL| ','RVV10 Functional: R. Sabatini et al, PRB 87: 041108(R) (2013)')")
351 END SELECT
352 IF (dispersion_env%verbose) THEN
353 WRITE (output_unit, &
354 fmt="(' vdW POTENTIAL| ',' Carrying out vdW-DF run using the following parameters:')")
355 WRITE (output_unit, fmt="(' vdW POTENTIAL| ','Nqs =',I8,' Nr_points =',I8,' r_max =',F10.3)") &
356 dispersion_env%nqs, dispersion_env%nr_points, dispersion_env%r_max
357 WRITE (output_unit, fmt="(' vdW POTENTIAL| ','q_mesh =')")
358 WRITE (output_unit, fmt="(8X,4F18.8)") (dispersion_env%q_mesh(i), i=1, dispersion_env%nqs)
359 WRITE (output_unit, &
360 fmt="(' vdW POTENTIAL| ','Density cutoff for convolution [a.u.]:',T71,F10.1)") &
361 dispersion_env%pw_cutoff
362 END IF
363 END IF
364 END IF
365 IF (.NOT. PRESENT(ounit)) THEN
366 CALL cp_print_key_finished_output(output_unit, logger, dft_section, &
367 "PRINT%DFT_CONTROL_PARAMETERS")
368 END IF
369
370 END SUBROUTINE qs_write_dispersion
371
372! **************************************************************************************************
373!> \brief ...
374!> \param scaling ...
375!> \param vdw_section ...
376! **************************************************************************************************
377 SUBROUTINE qs_scaling_dftd2(scaling, vdw_section)
378 REAL(kind=dp), INTENT(inout) :: scaling
379 TYPE(section_vals_type), POINTER :: vdw_section
380
381 CHARACTER(LEN=default_string_length) :: functional
382
383 CALL section_vals_val_get(vdw_section, "PAIR_POTENTIAL%REFERENCE_FUNCTIONAL", c_val=functional)
384
385 SELECT CASE (trim(functional))
386 CASE DEFAULT
387 ! unknown functional
388 cpabort("No DFT-D2 s6 value available for this functional:"//trim(functional))
389 CASE ("BLYP")
390 scaling = 1.20_dp
391 CASE ("B3LYP")
392 scaling = 1.05_dp
393 CASE ("TPSS")
394 scaling = 1.00_dp
395 CASE ("PBE")
396 scaling = 0.75_dp
397 CASE ("PBE0")
398 scaling = 0.6_dp
399 CASE ("B2PLYP")
400 scaling = 0.55_dp
401 CASE ("BP86")
402 scaling = 1.05_dp
403 CASE ("B97")
404 scaling = 1.25_dp
405 END SELECT
406
407 END SUBROUTINE qs_scaling_dftd2
408
409! **************************************************************************************************
410!> \brief ...
411!> \param s6 ...
412!> \param sr6 ...
413!> \param s8 ...
414!> \param vdw_section ...
415! **************************************************************************************************
416 SUBROUTINE qs_scaling_dftd3(s6, sr6, s8, vdw_section)
417
418 REAL(kind=dp), INTENT(inout) :: s6, sr6, s8
419 TYPE(section_vals_type), POINTER :: vdw_section
420
421 CHARACTER(LEN=default_string_length) :: functional
422
423 CALL section_vals_val_get(vdw_section, "PAIR_POTENTIAL%REFERENCE_FUNCTIONAL", c_val=functional)
424 CALL uppercase(functional)
425 ! values for different functionals from:
426 ! https://www.chemie.uni-bonn.de/grimme/de/software/dft-d3
427 ! L. Goerigk et al. PCCP 2017, 32147-32744, SI
428 SELECT CASE (trim(functional))
429 CASE DEFAULT
430 ! unknown functional
431 cpabort("No DFT-D3 values available for this functional:"//trim(functional))
432 CASE ("B1B95")
433 s6 = 1.000_dp
434 sr6 = 1.613_dp
435 s8 = 1.868_dp
436 CASE ("B2GPPLYP")
437 ! L. Goerigk and S. Grimme
438 ! J. Chem. Theory Comput. 2011, 7, 291-309; doi:10.1021/ct100466k
439 s6 = 0.56_dp
440 sr6 = 1.586_dp
441 s8 = 0.760_dp
442 CASE ("B2PLYP")
443 ! L. Goerigk and S. Grimme
444 ! J. Chem. Theory Comput. 2011, 7, 291-309; doi:10.1021/ct100466k
445 s6 = 0.64_dp
446 sr6 = 1.427_dp
447 s8 = 1.022_dp
448 CASE ("DSD-BLYP")
449 ! L. Goerigk and S. Grimme
450 ! J. Chem. Theory Comput. 2011, 7, 291-309; doi:10.1021/ct100466k
451 s6 = 0.50_dp
452 sr6 = 1.569_dp
453 s8 = 0.705_dp
454 CASE ("B3LYP")
455 ! S. Grimme, J. Antony, S. Ehrlich, and H. Krieg
456 ! J. Chem. Phys. 132,154104 (2010); doi:10.1063/1.3382344
457 s6 = 1.000_dp
458 sr6 = 1.261_dp
459 s8 = 1.703_dp
460 CASE ("B97-D")
461 ! S. Grimme, J. Antony, S. Ehrlich, and H. Krieg
462 ! J. Chem. Phys. 132,154104 (2010); doi:10.1063/1.3382344
463 s6 = 1.000_dp
464 sr6 = 0.892_dp
465 s8 = 0.909_dp
466 CASE ("BHLYP")
467 s6 = 1.000_dp
468 sr6 = 1.370_dp
469 s8 = 1.442_dp
470 CASE ("BLYP")
471 ! S. Grimme, J. Antony, S. Ehrlich, and H. Krieg
472 ! J. Chem. Phys. 132,154104 (2010); doi:10.1063/1.3382344
473 s6 = 1.000_dp
474 sr6 = 1.094_dp
475 s8 = 1.682_dp
476 CASE ("BP86")
477 ! S. Grimme, J. Antony, S. Ehrlich, and H. Krieg
478 ! J. Chem. Phys. 132,154104 (2010); doi:10.1063/1.3382344
479 s6 = 1.000_dp
480 sr6 = 1.139_dp
481 s8 = 1.683_dp
482 CASE ("BPBE")
483 s6 = 1.000_dp
484 sr6 = 1.087_dp
485 s8 = 2.033_dp
486 CASE ("MPWLYP")
487 s6 = 1.000_dp
488 sr6 = 1.239_dp
489 s8 = 1.098_dp
490 CASE ("PBE")
491 ! S. Grimme, J. Antony, S. Ehrlich, and H. Krieg
492 ! J. Chem. Phys. 132,154104 (2010); doi:10.1063/1.3382344
493 s6 = 1.000_dp
494 sr6 = 1.217_dp
495 s8 = 0.722_dp
496 CASE ("PBEHPBE")
497 s6 = 1.000_dp
498 sr6 = 1.5703_dp
499 s8 = 1.4010_dp
500 CASE ("PBE0")
501 ! S. Grimme, J. Antony, S. Ehrlich, and H. Krieg
502 ! J. Chem. Phys. 132,154104 (2010); doi:10.1063/1.3382344
503 s6 = 1.000_dp
504 sr6 = 1.287_dp
505 s8 = 0.928_dp
506 CASE ("PW6B95")
507 ! S. Grimme, J. Antony, S. Ehrlich, and H. Krieg
508 ! J. Chem. Phys. 132,154104 (2010); doi:10.1063/1.3382344
509 s6 = 1.000_dp
510 sr6 = 1.532_dp
511 s8 = 0.862_dp
512 CASE ("PWB6K")
513 s6 = 1.000_dp
514 sr6 = 1.660_dp
515 s8 = 0.550_dp
516 CASE ("REVPBE")
517 ! S. Grimme, J. Antony, S. Ehrlich, and H. Krieg
518 ! J. Chem. Phys. 132,154104 (2010); doi:10.1063/1.3382344
519 s6 = 1.000_dp
520 sr6 = 0.923_dp
521 s8 = 1.010_dp
522 CASE ("RPBE")
523 s6 = 1.000_dp
524 sr6 = 0.872_dp
525 s8 = 0.514_dp
526 CASE ("TPSS")
527 ! S. Grimme, J. Antony, S. Ehrlich, and H. Krieg
528 ! J. Chem. Phys. 132,154104 (2010); doi:10.1063/1.3382344
529 s6 = 1.000_dp
530 sr6 = 1.166_dp
531 s8 = 1.105_dp
532 CASE ("TPSS0")
533 ! S. Grimme, J. Antony, S. Ehrlich, and H. Krieg
534 ! J. Chem. Phys. 132,154104 (2010); doi:10.1063/1.3382344
535 s6 = 1.000_dp
536 sr6 = 1.252_dp
537 s8 = 1.242_dp
538 CASE ("TPSSH")
539 s6 = 1.000_dp
540 sr6 = 1.223_dp
541 s8 = 1.219_dp
542 CASE ("B1LYP")
543 s6 = 1.000_dp
544 sr6 = 1.3725_dp
545 s8 = 1.9467_dp
546 CASE ("B1P86")
547 s6 = 1.000_dp
548 sr6 = 1.1815_dp
549 s8 = 1.1209_dp
550 CASE ("B3P86")
551 s6 = 1.000_dp
552 sr6 = 1.1897_dp
553 s8 = 1.1961_dp
554 CASE ("B3PW91")
555 s6 = 1.000_dp
556 sr6 = 1.176_dp
557 s8 = 1.775_dp
558 CASE ("BMK")
559 s6 = 1.000_dp
560 sr6 = 1.931_dp
561 s8 = 2.168_dp
562 CASE ("CAMB3LYP")
563 s6 = 1.000_dp
564 sr6 = 1.378_dp
565 s8 = 1.217_dp
566 CASE ("LCWPBE")
567 s6 = 1.000_dp
568 sr6 = 1.355_dp
569 s8 = 1.279_dp
570 CASE ("M052X")
571 s6 = 1.000_dp
572 sr6 = 1.417_dp
573 s8 = 0.000_dp
574 CASE ("M05")
575 s6 = 1.000_dp
576 sr6 = 1.373_dp
577 s8 = 0.595_dp
578 CASE ("M062X")
579 s6 = 1.000_dp
580 sr6 = 1.619_dp
581 s8 = 0.000_dp
582 CASE ("M06HF")
583 s6 = 1.000_dp
584 sr6 = 1.446_dp
585 s8 = 0.000_dp
586 CASE ("M06L")
587 s6 = 1.000_dp
588 sr6 = 1.581_dp
589 s8 = 0.000_dp
590 CASE ("M06N")
591 s6 = 1.000_dp
592 sr6 = 1.325_dp
593 s8 = 0.000_dp
594 CASE ("HCTH120")
595 s6 = 1.000_dp
596 sr6 = 1.221_dp
597 s8 = 1.206_dp
598 CASE ("HCTH407")
599 s6 = 1.000_dp
600 sr6 = 4.0426_dp
601 s8 = 2.7694_dp
602 CASE ("MPW2PLYP")
603 s6 = 1.000_dp
604 sr6 = 1.5527_dp
605 s8 = 0.7529_dp
606 CASE ("PKZB")
607 s6 = 1.000_dp
608 sr6 = 0.6327_dp
609 s8 = 0.000_dp
610 CASE ("PTPSS")
611 s6 = 0.750_dp
612 sr6 = 1.541_dp
613 s8 = 0.879_dp
614 CASE ("PWPB95")
615 s6 = 0.820_dp
616 sr6 = 1.557_dp
617 s8 = 0.705_dp
618 CASE ("OLYP")
619 s6 = 1.000_dp
620 sr6 = 0.806_dp
621 s8 = 1.764_dp
622 CASE ("OPBE")
623 s6 = 1.000_dp
624 sr6 = 0.837_dp
625 s8 = 2.055_dp
626 CASE ("OTPSS")
627 s6 = 1.000_dp
628 sr6 = 1.128_dp
629 s8 = 1.494_dp
630 CASE ("PBE1KCIS")
631 s6 = 1.000_dp
632 sr6 = 3.6355_dp
633 s8 = 1.7934_dp
634 CASE ("PBE38")
635 s6 = 1.000_dp
636 sr6 = 1.333_dp
637 s8 = 0.998_dp
638 CASE ("PBEH1PBE")
639 s6 = 1.000_dp
640 sr6 = 1.3719_dp
641 s8 = 1.0430_dp
642 CASE ("PBESOL")
643 s6 = 1.000_dp
644 sr6 = 1.345_dp
645 s8 = 0.612_dp
646 CASE ("REVSSB")
647 s6 = 1.000_dp
648 sr6 = 1.221_dp
649 s8 = 0.560_dp
650 CASE ("REVTPSS")
651 s6 = 1.000_dp
652 sr6 = 1.3491_dp
653 s8 = 1.3666_dp
654 CASE ("SSB")
655 s6 = 1.000_dp
656 sr6 = 1.215_dp
657 s8 = 0.663_dp
658 CASE ("B97-1")
659 s6 = 1.000_dp
660 sr6 = 3.7924_dp
661 s8 = 1.6418_dp
662 CASE ("B97-2")
663 s6 = 1.000_dp
664 sr6 = 1.7066_dp
665 s8 = 2.4661_dp
666 CASE ("B98")
667 s6 = 1.000_dp
668 sr6 = 2.6895_dp
669 s8 = 1.9078_dp
670 CASE ("BOP")
671 s6 = 1.000_dp
672 sr6 = 0.929_dp
673 s8 = 1.975_dp
674 CASE ("HISS")
675 s6 = 1.000_dp
676 sr6 = 1.3338_dp
677 s8 = 0.7615_dp
678 CASE ("HSE03")
679 s6 = 1.000_dp
680 sr6 = 1.3944_dp
681 s8 = 1.0156_dp
682 CASE ("HSE06")
683 s6 = 1.000_dp
684 sr6 = 1.129_dp
685 s8 = 0.109_dp
686 CASE ("M08HX")
687 s6 = 1.000_dp
688 sr6 = 1.6247_dp
689 s8 = 0.000_dp
690 CASE ("MN15L")
691 s6 = 1.000_dp
692 sr6 = 3.3388_dp
693 s8 = 0.000_dp
694 CASE ("MPWPW91")
695 s6 = 1.0000_dp
696 sr6 = 1.3725_dp
697 s8 = 1.9467_dp
698 CASE ("MPW1B95")
699 s6 = 1.000_dp
700 sr6 = 1.605_dp
701 s8 = 1.118_dp
702 CASE ("MPW1KCIS")
703 s6 = 1.000_dp
704 sr6 = 1.7231_dp
705 s8 = 2.2917_dp
706 CASE ("MPW1LYP")
707 s6 = 1.000_dp
708 sr6 = 2.0512_dp
709 s8 = 1.9529_dp
710 CASE ("MPW1PW91")
711 s6 = 1.000_dp
712 sr6 = 1.2892_dp
713 s8 = 1.4758_dp
714 CASE ("MPWB1K")
715 s6 = 1.000_dp
716 sr6 = 1.671_dp
717 s8 = 1.061_dp
718 CASE ("MPWKCIS1K")
719 s6 = 1.000_dp
720 sr6 = 1.4853_dp
721 s8 = 1.7553_dp
722 CASE ("O3LYP")
723 s6 = 1.000_dp
724 sr6 = 1.4060_dp
725 s8 = 1.8058_dp
726 CASE ("PW1PW")
727 s6 = 1.000_dp
728 sr6 = 1.4968_dp
729 s8 = 1.1786_dp
730 CASE ("PW91P86")
731 s6 = 1.0000_dp
732 sr6 = 2.1040_dp
733 s8 = 0.8747_dp
734 CASE ("REVPBE0")
735 s6 = 1.000_dp
736 sr6 = 0.949_dp
737 s8 = 0.792_dp
738 CASE ("REVPBE38")
739 s6 = 1.000_dp
740 sr6 = 1.021_dp
741 s8 = 0.862_dp
742 CASE ("REVTPSSh")
743 s6 = 1.000_dp
744 sr6 = 1.3224_dp
745 s8 = 1.2504_dp
746 CASE ("REVTPSS0")
747 s6 = 1.000_dp
748 sr6 = 1.2881_dp
749 s8 = 1.0649_dp
750 CASE ("TPSS1KCIS")
751 s6 = 1.000_dp
752 sr6 = 1.7729_dp
753 s8 = 2.0902_dp
754 CASE ("THCTHHYB")
755 s6 = 1.000_dp
756 sr6 = 1.5001_dp
757 s8 = 1.6302_dp
758 CASE ("RPW86PBE")
759 s6 = 1.000_dp
760 sr6 = 1.224_dp
761 s8 = 0.901_dp
762 CASE ("SCAN")
763 s6 = 1.000_dp
764 sr6 = 1.324_dp
765 s8 = 0.000_dp
766 CASE ("THCTH")
767 s6 = 1.000_dp
768 sr6 = 0.932_dp
769 s8 = 0.5662_dp
770 CASE ("XLYP")
771 s6 = 1.0000_dp
772 sr6 = 0.9384_dp
773 s8 = 0.7447_dp
774 CASE ("X3LYP")
775 s6 = 1.000_dp
776 sr6 = 1.0000_dp
777 s8 = 0.2990_dp
778 END SELECT
779
780 END SUBROUTINE qs_scaling_dftd3
781
782! **************************************************************************************************
783!> \brief ...
784!> \param s6 ...
785!> \param a1 ...
786!> \param s8 ...
787!> \param a2 ...
788!> \param vdw_section ...
789! **************************************************************************************************
790 SUBROUTINE qs_scaling_dftd3bj(s6, a1, s8, a2, vdw_section)
791 REAL(kind=dp), INTENT(inout) :: s6, a1, s8, a2
792 TYPE(section_vals_type), POINTER :: vdw_section
793
794 CHARACTER(LEN=default_string_length) :: functional
795
796 CALL section_vals_val_get(vdw_section, "PAIR_POTENTIAL%REFERENCE_FUNCTIONAL", c_val=functional)
797
798 ! values for different functionals from:
799 ! http://www.thch.uni-bonn.de/tc/downloads/DFT-D3/functionalsbj.html
800 ! L. Goerigk et al. PCCP 2017, 32147-32744, SI
801 SELECT CASE (trim(functional))
802 CASE DEFAULT
803 ! unknown functional
804 cpabort("No DFT-D3(BJ) values available for this functional:"//trim(functional))
805 CASE ("B1B95")
806 s6 = 1.0000_dp
807 a1 = 0.2092_dp
808 s8 = 1.4507_dp
809 a2 = 5.5545_dp
810 CASE ("B2GPPLYP")
811 s6 = 0.5600_dp
812 a1 = 0.0000_dp
813 s8 = 0.2597_dp
814 a2 = 6.3332_dp
815 CASE ("B3PW91")
816 s6 = 1.0000_dp
817 a1 = 0.4312_dp
818 s8 = 2.8524_dp
819 a2 = 4.4693_dp
820 CASE ("BHLYP")
821 s6 = 1.0000_dp
822 a1 = 0.2793_dp
823 s8 = 1.0354_dp
824 a2 = 4.9615_dp
825 CASE ("BMK")
826 s6 = 1.0000_dp
827 a1 = 0.1940_dp
828 s8 = 2.0860_dp
829 a2 = 5.9197_dp
830 CASE ("BOP")
831 s6 = 1.0000_dp
832 a1 = 0.4870_dp
833 s8 = 3.2950_dp
834 a2 = 3.5043_dp
835 CASE ("BPBE")
836 s6 = 1.0000_dp
837 a1 = 0.4567_dp
838 s8 = 4.0728_dp
839 a2 = 4.3908_dp
840 CASE ("B97-3c")
841 s6 = 1.0000_dp
842 a1 = 0.3700_dp
843 s8 = 1.5000_dp
844 a2 = 4.1000_dp
845 CASE ("CAMB3LYP")
846 s6 = 1.0000_dp
847 a1 = 0.3708_dp
848 s8 = 2.0674_dp
849 a2 = 5.4743_dp
850 CASE ("DSDBLYP")
851 s6 = 0.5000_dp
852 a1 = 0.0000_dp
853 s8 = 0.2130_dp
854 a2 = 6.0519_dp
855 CASE ("DSDPBEP86")
856 s6 = 0.4180_dp
857 a1 = 0.0000_dp
858 s8 = 0.0000_dp
859 a2 = 5.6500_dp
860 CASE ("DSDPBEB95")
861 s6 = 0.6100_dp
862 a1 = 0.0000_dp
863 s8 = 0.0000_dp
864 a2 = 6.2000_dp
865 CASE ("LCWPBE")
866 s6 = 1.0000_dp
867 a1 = 0.3919_dp
868 s8 = 1.8541_dp
869 a2 = 5.0897_dp
870 CASE ("LCWhPBE")
871 s6 = 1.0000_dp
872 a1 = 0.2746_dp
873 s8 = 1.1908_dp
874 a2 = 5.3157_dp
875 CASE ("MPW1B95")
876 s6 = 1.0000_dp
877 a1 = 0.1955_dp
878 s8 = 1.0508_dp
879 a2 = 6.4177_dp
880 CASE ("MPW2PLYP")
881 s6 = 0.6600_dp
882 a1 = 0.4105_dp
883 s8 = 0.6223_dp
884 a2 = 5.0136_dp
885 CASE ("MPWB1K")
886 s6 = 1.0000_dp
887 a1 = 0.1474_dp
888 s8 = 0.9499_dp
889 a2 = 6.6223_dp
890 CASE ("MPWLYP")
891 s6 = 1.0000_dp
892 a1 = 0.4831_dp
893 s8 = 2.0077_dp
894 a2 = 4.5323_dp
895 CASE ("OLYP")
896 s6 = 1.0000_dp
897 a1 = 0.5299_dp
898 s8 = 2.6205_dp
899 a2 = 2.8065_dp
900 CASE ("OPBE")
901 s6 = 1.0000_dp
902 a1 = 0.5512_dp
903 s8 = 3.3816_dp
904 a2 = 2.9444_dp
905 CASE ("OTPSS")
906 s6 = 1.0000_dp
907 a1 = 0.4634_dp
908 s8 = 2.7495_dp
909 a2 = 4.3153_dp
910 CASE ("PBE38")
911 s6 = 1.0000_dp
912 a1 = 0.3995_dp
913 s8 = 1.4623_dp
914 a2 = 5.1405_dp
915 CASE ("PBEsol")
916 s6 = 1.0000_dp
917 a1 = 0.4466_dp
918 s8 = 2.9491_dp
919 a2 = 6.1742_dp
920 CASE ("PTPSS")
921 s6 = 0.7500_dp
922 a1 = 0.0000_dp
923 s8 = 0.2804_dp
924 a2 = 6.5745_dp
925 CASE ("PWB6K")
926 s6 = 1.0000_dp
927 a1 = 0.1805_dp
928 s8 = 0.9383_dp
929 a2 = 7.7627_dp
930 CASE ("revSSB")
931 s6 = 1.0000_dp
932 a1 = 0.4720_dp
933 s8 = 0.4389_dp
934 a2 = 4.0986_dp
935 CASE ("SSB")
936 s6 = 1.0000_dp
937 a1 = -0.0952_dp
938 s8 = -0.1744_dp
939 a2 = 5.2170_dp
940 CASE ("TPSSh")
941 s6 = 1.0000_dp
942 a1 = 0.4529_dp
943 s8 = 2.2382_dp
944 a2 = 4.6550_dp
945 CASE ("HCTH120")
946 s6 = 1.0000_dp
947 a1 = 0.3563_dp
948 s8 = 1.0821_dp
949 a2 = 4.3359_dp
950 CASE ("B2PLYP")
951 s6 = 0.6400_dp
952 a1 = 0.3065_dp
953 s8 = 0.9147_dp
954 a2 = 5.0570_dp
955 CASE ("B1LYP")
956 s6 = 1.0000_dp
957 a1 = 0.1986_dp
958 s8 = 2.1167_dp
959 a2 = 5.3875_dp
960 CASE ("B1P86")
961 s6 = 1.0000_dp
962 a1 = 0.4724_dp
963 s8 = 3.5681_dp
964 a2 = 4.9858_dp
965 CASE ("B3LYP")
966 s6 = 1.0000_dp
967 a1 = 0.3981_dp
968 s8 = 1.9889_dp
969 a2 = 4.4211_dp
970 CASE ("B3P86")
971 s6 = 1.0000_dp
972 a1 = 0.4601_dp
973 s8 = 3.3211_dp
974 a2 = 4.9294_dp
975 CASE ("B97-1")
976 s6 = 1.0000_dp
977 a1 = 0.0000_dp
978 s8 = 0.4814_dp
979 a2 = 6.2279_dp
980 CASE ("B97-2")
981 s6 = 1.0000_dp
982 a1 = 0.0000_dp
983 s8 = 0.9448_dp
984 a2 = 5.4603_dp
985 CASE ("B97-D")
986 s6 = 1.0000_dp
987 a1 = 0.5545_dp
988 s8 = 2.2609_dp
989 a2 = 3.2297_dp
990 CASE ("B98")
991 s6 = 1.0000_dp
992 a1 = 0.0000_dp
993 s8 = 0.7086_dp
994 a2 = 6.0672_dp
995 CASE ("BLYP")
996 s6 = 1.0000_dp
997 a1 = 0.4298_dp
998 s8 = 2.6996_dp
999 a2 = 4.2359_dp
1000 CASE ("BP86")
1001 s6 = 1.0000_dp
1002 a1 = 0.3946_dp
1003 s8 = 3.2822_dp
1004 a2 = 4.8516_dp
1005 CASE ("DSD-BLYP")
1006 s6 = 0.5000_dp
1007 a1 = 0.0000_dp
1008 s8 = 0.2130_dp
1009 a2 = 6.0519_dp
1010 CASE ("HCTH407")
1011 s6 = 1.0000_dp
1012 a1 = 0.0000_dp
1013 s8 = 0.6490_dp
1014 a2 = 4.8162_dp
1015 CASE ("HISS")
1016 s6 = 1.0000_dp
1017 a1 = 0.0000_dp
1018 s8 = 1.6112_dp
1019 a2 = 7.3539_dp
1020 CASE ("HSE03")
1021 s6 = 1.0000_dp
1022 a1 = 0.0000_dp
1023 s8 = 1.1243_dp
1024 a2 = 6.8889_dp
1025 CASE ("HSE06")
1026 s6 = 1.0000_dp
1027 a1 = 0.3830_dp
1028 s8 = 2.3100_dp
1029 a2 = 5.6850_dp
1030 CASE ("M11")
1031 s6 = 1.0000_dp
1032 a1 = 0.0000_dp
1033 s8 = 2.8112_dp
1034 a2 = 10.1389_dp
1035 CASE ("MN12SX")
1036 s6 = 1.0000_dp
1037 a1 = 0.0983_dp
1038 s8 = 1.1674_dp
1039 a2 = 8.0259_dp
1040 CASE ("MN15")
1041 s6 = 1.0000_dp
1042 a1 = 2.0971_dp
1043 s8 = 0.7862_dp
1044 a2 = 7.5923_dp
1045 CASE ("mPWPW91")
1046 s6 = 1.0000_dp
1047 a1 = 0.3168_dp
1048 s8 = 1.7974_dp
1049 a2 = 4.7732_dp
1050 CASE ("MPW1PW91")
1051 s6 = 1.0000_dp
1052 a1 = 0.3342_dp
1053 s8 = 1.8744_dp
1054 a2 = 4.9819_dp
1055 CASE ("MPW1KCIS")
1056 s6 = 1.0000_dp
1057 a1 = 0.0576_dp
1058 s8 = 1.0893_dp
1059 a2 = 5.5314_dp
1060 CASE ("MPWKCIS1K")
1061 s6 = 1.0000_dp
1062 a1 = 0.0855_dp
1063 s8 = 1.2875_dp
1064 a2 = 5.8961_dp
1065 CASE ("N12SX")
1066 s6 = 1.0000_dp
1067 a1 = 0.3283_dp
1068 s8 = 2.4900_dp
1069 a2 = 5.7898_dp
1070 CASE ("O3LYP")
1071 s6 = 1.0000_dp
1072 a1 = 0.0963_dp
1073 s8 = 1.8171_dp
1074 a2 = 5.9940_dp
1075 CASE ("PBE0")
1076 s6 = 1.0000_dp
1077 a1 = 0.4145_dp
1078 s8 = 1.2177_dp
1079 a2 = 4.8593_dp
1080 CASE ("PBE")
1081 s6 = 1.0000_dp
1082 a1 = 0.4289_dp
1083 s8 = 0.7875_dp
1084 a2 = 4.4407_dp
1085 CASE ("PBEhPBE")
1086 s6 = 1.0000_dp
1087 a1 = 0.0000_dp
1088 s8 = 1.1152_dp
1089 a2 = 6.7184_dp
1090 CASE ("PBEh1PBE")
1091 s6 = 1.0000_dp
1092 a1 = 0.0000_dp
1093 s8 = 1.4877_dp
1094 a2 = 7.0385_dp
1095 CASE ("PBE1KCIS")
1096 s6 = 1.0000_dp
1097 a1 = 0.0000_dp
1098 s8 = 0.7688_dp
1099 a2 = 6.2794_dp
1100 CASE ("PW6B95")
1101 s6 = 1.0000_dp
1102 a1 = 0.2076_dp
1103 s8 = 0.7257_dp
1104 a2 = 6.3750_dp
1105 CASE ("PWPB95")
1106 s6 = 0.8200_dp
1107 a1 = 0.0000_dp
1108 s8 = 0.2904_dp
1109 a2 = 7.3141_dp
1110 CASE ("revPBE0")
1111 s6 = 1.0000_dp
1112 a1 = 0.4679_dp
1113 s8 = 1.7588_dp
1114 a2 = 3.7619_dp
1115 CASE ("revPBE38")
1116 s6 = 1.0000_dp
1117 a1 = 0.4309_dp
1118 s8 = 1.4760_dp
1119 a2 = 3.9446_dp
1120 CASE ("revPBE")
1121 s6 = 1.0000_dp
1122 a1 = 0.5238_dp
1123 s8 = 2.3550_dp
1124 a2 = 3.5016_dp
1125 CASE ("revTPSS")
1126 s6 = 1.0000_dp
1127 a1 = 0.4426_dp
1128 s8 = 1.4023_dp
1129 a2 = 4.4723_dp
1130 CASE ("revTPSS0")
1131 s6 = 1.0000_dp
1132 a1 = 0.2218_dp
1133 s8 = 1.6151_dp
1134 a2 = 5.7985_dp
1135 CASE ("revTPSSh")
1136 s6 = 1.0000_dp
1137 a1 = 0.2660_dp
1138 s8 = 1.4076_dp
1139 a2 = 5.3761_dp
1140 CASE ("RPBE")
1141 s6 = 1.0000_dp
1142 a1 = 0.1820_dp
1143 s8 = 0.8318_dp
1144 a2 = 4.0094_dp
1145 CASE ("RPW86PBE")
1146 s6 = 1.0000_dp
1147 a1 = 0.4613_dp
1148 s8 = 1.3845_dp
1149 a2 = 4.5062_dp
1150 CASE ("SCAN")
1151 s6 = 1.0000_dp
1152 a1 = 0.538_dp
1153 s8 = 0.0000_dp
1154 a2 = 5.420_dp
1155 CASE ("SOGGA11X")
1156 s6 = 1.0000_dp
1157 a1 = 0.1330_dp
1158 s8 = 1.1426_dp
1159 a2 = 5.7381_dp
1160 CASE ("TPSS0")
1161 s6 = 1.0000_dp
1162 a1 = 0.3768_dp
1163 s8 = 1.2576_dp
1164 a2 = 4.5865_dp
1165 CASE ("TPSS1KCIS")
1166 s6 = 1.0000_dp
1167 a1 = 0.0000_dp
1168 s8 = 1.0542_dp
1169 a2 = 6.0201_dp
1170 CASE ("TPSS")
1171 s6 = 1.0000_dp
1172 a1 = 0.4535_dp
1173 s8 = 1.9435_dp
1174 a2 = 4.4752_dp
1175 CASE ("tHCTH")
1176 s6 = 1.0000_dp
1177 a1 = 0.0000_dp
1178 s8 = 1.2626_dp
1179 a2 = 5.6162_dp
1180 CASE ("tHCTHhyb")
1181 s6 = 1.0000_dp
1182 a1 = 0.0000_dp
1183 s8 = 0.9585_dp
1184 a2 = 6.2303_dp
1185 CASE ("XLYP")
1186 s6 = 1.0000_dp
1187 a1 = 0.0809_dp
1188 s8 = 1.5669_dp
1189 a2 = 5.3166_dp
1190 CASE ("X3LYP")
1191 s6 = 1.0000_dp
1192 a1 = 0.2022_dp
1193 s8 = 1.5744_dp
1194 a2 = 5.4184_dp
1195 END SELECT
1196
1197 END SUBROUTINE qs_scaling_dftd3bj
1198
1199! **************************************************************************************************
1200!> \brief ...
1201!> \param cell ...
1202!> \param ncell ...
1203!> \return ...
1204! **************************************************************************************************
1205 FUNCTION cellhash(cell, ncell) RESULT(hash)
1206 INTEGER, DIMENSION(3), INTENT(IN) :: cell, ncell
1207 INTEGER :: hash
1208
1209 INTEGER :: ix, iy, iz, nx, ny, nz
1210
1211 cpassert(all(abs(cell) <= ncell))
1212
1213 ix = cell(1)
1214 IF (ix /= 0) THEN
1215 ix = 2*abs(ix) - (1 + sign(1, ix))/2
1216 END IF
1217 iy = cell(2)
1218 IF (iy /= 0) THEN
1219 iy = 2*abs(iy) - (1 + sign(1, iy))/2
1220 END IF
1221 iz = cell(3)
1222 IF (iz /= 0) THEN
1223 iz = 2*abs(iz) - (1 + sign(1, iz))/2
1224 END IF
1225
1226 nx = 2*ncell(1) + 1
1227 ny = 2*ncell(2) + 1
1228 nz = 2*ncell(3) + 1
1229
1230 hash = ix*ny*nz + iy*nz + iz + 1
1231
1232 END FUNCTION cellhash
1233! **************************************************************************************************
1234
1235END MODULE qs_dispersion_utils
1236
static unsigned int hash(const dbm_task_t task)
Private hash function based on Szudzik's elegant pairing. Using unsigned int to return a positive num...
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind(atomic_kind, fist_potential, element_symbol, name, mass, kind_number, natom, atom_list, rcov, rvdw, z, qeff, apol, cpol, mm_radius, shell, shell_active, damping)
Get attributes of an atomic kind.
various routines to log and control the output. The idea is that decisions about where to log should ...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer function, public cp_print_key_unit_nr(logger, basis_section, print_key_path, extension, middle_name, local, log_filename, ignore_should_output, file_form, file_position, file_action, file_status, do_backup, on_file, is_new_file, mpi_io, fout)
...
subroutine, public cp_print_key_finished_output(unit_nr, logger, basis_section, print_key_path, local, ignore_should_output, on_file, mpi_io)
should be called after you finish working with a unit obtained with cp_print_key_unit_nr,...
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public vdw_nl_rvv10
integer, parameter, public vdw_pairpot_dftd3
integer, parameter, public xc_vdw_fun_nonloc
integer, parameter, public vdw_pairpot_dftd4
integer, parameter, public vdw_nl_drsll
integer, parameter, public vdw_nl_lmkll
integer, parameter, public vdw_pairpot_dftd2
integer, parameter, public xc_vdw_fun_pairpot
integer, parameter, public vdw_pairpot_dftd3bj
objects that represent the structure of input sections and the data contained in an input section
recursive type(section_vals_type) function, pointer, public section_vals_get_subs_vals(section_vals, subsection_name, i_rep_section, can_return_null)
returns the values of the requested subsection
subroutine, public section_vals_val_get(section_vals, keyword_name, i_rep_section, i_rep_val, n_rep_val, val, l_val, i_val, r_val, c_val, l_vals, i_vals, r_vals, c_vals, explicit)
returns the requested value
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
integer, parameter, public default_string_length
Definition kinds.F:57
Definition of physical constants:
Definition physcon.F:68
real(kind=dp), parameter, public kjmol
Definition physcon.F:168
real(kind=dp), parameter, public bohr
Definition physcon.F:147
Definition of disperson types for DFT calculations.
Set disperson types for DFT calculations.
subroutine, public qs_dispersion_env_set(dispersion_env, xc_section)
...
integer function, public cellhash(cell, ncell)
...
subroutine, public qs_write_dispersion(qs_env, dispersion_env, ounit)
...
subroutine, public get_qs_env(qs_env, atomic_kind_set, qs_kind_set, cell, super_cell, cell_ref, use_ref_cell, kpoints, dft_control, mos, sab_orb, sab_all, qmmm, qmmm_periodic, mimic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, sab_cneo, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, xcint_weights, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, rhoz_cneo_set, ecoul_1c, rho0_s_rs, rho0_s_gs, rhoz_cneo_s_rs, rhoz_cneo_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, harris_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs, do_rixs, tb_tblite)
Get the QUICKSTEP environment.
Define the quickstep kind type and their sub types.
subroutine, public get_qs_kind(qs_kind, basis_set, basis_type, ncgf, nsgf, all_potential, tnadd_potential, gth_potential, sgp_potential, upf_potential, cneo_potential, se_parameter, dftb_parameter, xtb_parameter, dftb3_param, zatom, zeff, elec_conf, mao, lmax_dftb, alpha_core_charge, ccore_charge, core_charge, core_charge_radius, paw_proj_set, paw_atom, hard_radius, hard0_radius, max_rad_local, covalent_radius, vdw_radius, gpw_type_forced, harmonics, max_iso_not0, max_s_harm, grid_atom, ngrid_ang, ngrid_rad, lmax_rho0, dft_plus_u_atom, l_of_dft_plus_u, n_of_dft_plus_u, u_minus_j, u_of_dft_plus_u, j_of_dft_plus_u, alpha_of_dft_plus_u, beta_of_dft_plus_u, j0_of_dft_plus_u, occupation_of_dft_plus_u, dispersion, bs_occupation, magnetization, no_optimize, addel, laddel, naddel, orbitals, max_scf, eps_scf, smear, u_ramping, u_minus_j_target, eps_u_ramping, init_u_ramping_each_scf, reltmat, ghost, monovalent, floating, name, element_symbol, pao_basis_size, pao_model_file, pao_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
Utilities for string manipulations.
elemental subroutine, public uppercase(string)
Convert all lower case characters in a string to upper case.
Provides all information about an atomic kind.
type of a logger, at the moment it contains just a print level starting at which level it should be l...
Provides all information about a quickstep kind.