44 #include "./base/base_uses.f90"
50 INTEGER,
PARAMETER,
PRIVATE :: nelem = 106
61 INTEGER,
DIMENSION(0:nelem), &
62 PARAMETER,
PRIVATE :: zval = (/-1, &
64 1, 2, 3, 4, 5, 6, 7, 8, &
65 1, 2, 3, 4, 5, 6, 7, 8, &
66 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 2, 3, 4, 5, 6, 7, 8, &
67 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 2, 3, 4, 5, 6, 7, 8, &
68 1, 2, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, &
69 4, 5, 6, 7, 8, 9, 10, 11, 2, 3, 4, 5, 6, 7, 8, &
70 -1, -1, -1, 4, -1, 6, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1/)
75 REAL(KIND=
dp),
DIMENSION(0:nelem), &
76 PARAMETER,
PRIVATE :: eneg = (/0.00_dp, &
78 0.98_dp, 1.57_dp, 2.04_dp, 2.55_dp, 3.04_dp, 3.44_dp, 3.98_dp, 4.50_dp, &
79 0.93_dp, 1.31_dp, 1.61_dp, 1.90_dp, 2.19_dp, 2.58_dp, 3.16_dp, 3.50_dp, &
80 0.82_dp, 1.00_dp, 1.36_dp, 1.54_dp, 1.63_dp, 1.66_dp, 1.55_dp, 1.83_dp, &
81 1.88_dp, 1.91_dp, 1.90_dp, 1.65_dp, 1.81_dp, 2.01_dp, 2.18_dp, 2.55_dp, 2.96_dp, 3.00_dp, &
82 0.82_dp, 0.95_dp, 1.22_dp, 1.33_dp, 1.60_dp, 2.16_dp, 1.90_dp, 2.20_dp, &
83 2.28_dp, 2.20_dp, 1.93_dp, 1.69_dp, 1.78_dp, 1.96_dp, 2.05_dp, 2.10_dp, 2.66_dp, 2.60_dp, &
84 0.79_dp, 0.89_dp, 1.10_dp, &
85 1.12_dp, 1.13_dp, 1.14_dp, 1.15_dp, 1.17_dp, 1.18_dp, 1.20_dp, 1.21_dp, &
86 1.22_dp, 1.23_dp, 1.24_dp, 1.25_dp, 1.26_dp, 1.27_dp, &
87 1.30_dp, 1.50_dp, 2.36_dp, 1.90_dp, 2.20_dp, 2.20_dp, 2.28_dp, 2.54_dp, &
88 2.00_dp, 2.04_dp, 2.33_dp, 2.02_dp, 2.00_dp, 2.20_dp, 2.20_dp, &
89 0.70_dp, 0.89_dp, 1.10_dp, &
90 1.30_dp, 1.50_dp, 1.38_dp, 1.36_dp, 1.28_dp, 1.30_dp, 1.30_dp, 1.30_dp, &
91 1.30_dp, 1.30_dp, 1.30_dp, 1.30_dp, 1.30_dp, 1.50_dp, &
92 1.50_dp, 1.50_dp, 1.50_dp/)
97 INTEGER,
DIMENSION(1:5, 0:nelem) :: occupation = reshape((/0,0,0,0,0, &
98 1,0,0,0,0, 2,0,0,0,0, &
99 1,0,0,0,0, 2,0,0,0,0, 2,1,0,0,0, 2,2,0,0,0, 2,3,0,0,0, 2,4,0,0,0, 2,5,0,0,0, 2,6,0,0,0, &
100 1,0,0,0,0, 2,0,0,0,0, 2,1,0,0,0, 2,2,0,0,0, 2,3,0,0,0, 2,4,0,0,0, 2,5,0,0,0, 2,6,0,0,0, &
101 1,0,0,0,0, 2,0,0,0,0, 2,0,1,0,0, 2,0,2,0,0, 2,0,3,0,0, 2,0,4,0,0, 2,0,5,0,0, 2,0,6,0,0, &
102 2,0,7,0,0, 2,0,8,0,0, 2,0,9,0,0, 2,0,0,0,0, 2,1,0,0,0, 2,2,0,0,0, 2,3,0,0,0, 2,4,0,0,0, 2,5,0,0,0, 2,6,0,0,0, &
103 1,0,0,0,0, 2,0,0,0,0, 2,0,1,0,0, 2,0,2,0,0, 2,0,3,0,0, 2,0,4,0,0, 2,0,5,0,0, 2,0,6,0,0, &
104 2,0,7,0,0, 2,0,8,0,0, 2,0,9,0,0, 2,0,0,0,0, 2,1,0,0,0, 2,2,0,0,0, 2,3,0,0,0, 2,4,0,0,0, 2,5,0,0,0, 2,6,0,0,0, &
105 1,0,0,0,0, 2,0,0,0,0, 2,0,1,0,0, &
106 2,0,1,0,0, 2,0,1,0,0, 2,0,1,0,0, 2,0,1,0,0, 2,0,1,0,0, 2,0,1,0,0, 2,0,1,0,0, 2,0,1,0,0, &
107 2,0,1,0,0, 2,0,1,0,0, 2,0,1,0,0, 2,0,1,0,0, 2,0,1,0,0, 2,0,1,0,0, &
108 2,0,2,0,0, 2,0,3,0,0, 2,0,4,0,0, 2,0,5,0,0, 2,0,6,0,0, 2,0,7,0,0, 2,0,8,0,0, 2,0,9,0,0, &
109 2,0,0,0,0, 2,1,0,0,0, 2,2,0,0,0, 2,3,0,0,0, 2,4,0,0,0, 2,5,0,0,0, 2,6,0,0,0, &
110 0,0,0,0,0, 0,0,0,0,0, 0,0,0,0,0, &
111 0,0,0,0,0, 0,0,0,0,0, 0,0,0,0,0, 0,0,0,0,0, 0,0,0,0,0, 0,0,0,0,0, 0,0,0,0,0, 0,0,0,0,0, &
112 0,0,0,0,0, 0,0,0,0,0, 0,0,0,0,0, 0,0,0,0,0, 0,0,0,0,0, 0,0,0,0,0, &
113 0,0,0,0,0, 0,0,0,0,0, 0,0,0,0,0/), (/5, nelem+1/))
122 REAL(KIND=
dp),
DIMENSION(0:nelem), &
123 PARAMETER,
PRIVATE :: crad = (/0.00_dp, &
125 1.30_dp, 0.99_dp, 0.84_dp, 0.75_dp, 0.71_dp, 0.64_dp, 0.60_dp, 0.62_dp, &
126 1.60_dp, 1.40_dp, 1.24_dp, 1.14_dp, 1.09_dp, 1.04_dp, 1.00_dp, 1.01_dp, &
127 2.00_dp, 1.74_dp, 1.59_dp, 1.48_dp, 1.44_dp, 1.30_dp, 1.29_dp, 1.24_dp, &
128 1.18_dp, 1.17_dp, 1.22_dp, 1.20_dp, 1.23_dp, 1.20_dp, 1.20_dp, 1.18_dp, 1.17_dp, 1.16_dp, &
129 2.15_dp, 1.90_dp, 1.76_dp, 1.64_dp, 1.56_dp, 1.46_dp, 1.38_dp, 1.36_dp, &
130 1.34_dp, 1.30_dp, 1.36_dp, 1.40_dp, 1.42_dp, 1.40_dp, 1.40_dp, 1.37_dp, 1.36_dp, 1.36_dp, &
131 2.38_dp, 2.06_dp, 1.94_dp, &
132 1.84_dp, 1.90_dp, 1.88_dp, 1.86_dp, 1.85_dp, 1.83_dp, 1.82_dp, 1.81_dp, &
133 1.80_dp, 1.79_dp, 1.77_dp, 1.77_dp, 1.78_dp, 1.74_dp, &
134 1.64_dp, 1.58_dp, 1.50_dp, 1.41_dp, 1.36_dp, 1.32_dp, 1.30_dp, 1.30_dp, &
135 1.32_dp, 1.44_dp, 1.45_dp, 1.50_dp, 1.42_dp, 1.48_dp, 1.46_dp, &
136 2.42_dp, 2.11_dp, 2.01_dp, &
137 1.90_dp, 1.84_dp, 1.83_dp, 1.80_dp, 1.80_dp, 1.51_dp, 0.96_dp, 1.54_dp, &
138 1.83_dp, 1.50_dp, 1.50_dp, 1.50_dp, 1.50_dp, 1.50_dp, &
139 1.50_dp, 1.50_dp, 1.50_dp/)
144 REAL(KIND=
dp),
DIMENSION(0:nelem), &
145 PARAMETER,
PRIVATE :: clmt = (/0.00_dp, &
147 1.05_dp, 2.05_dp, 3.00_dp, 4.00_dp, 3.00_dp, 2.00_dp, 1.25_dp, 1.00_dp, &
148 1.05_dp, 2.05_dp, 3.00_dp, 4.00_dp, 3.00_dp, 2.00_dp, 1.25_dp, 1.00_dp, &
149 1.05_dp, 2.05_dp, 3.50_dp, 3.50_dp, 3.50_dp, 3.50_dp, 3.50_dp, 3.50_dp, &
150 3.50_dp, 3.50_dp, 3.50_dp, 2.50_dp, 2.50_dp, 3.50_dp, 3.50_dp, 3.50_dp, 1.25_dp, 1.00_dp, &
151 1.05_dp, 2.05_dp, 3.50_dp, 3.50_dp, 3.50_dp, 3.50_dp, 3.50_dp, 3.50_dp, &
152 3.50_dp, 3.50_dp, 3.50_dp, 2.50_dp, 2.50_dp, 3.50_dp, 3.50_dp, 3.50_dp, 1.25_dp, 1.00_dp, &
153 1.05_dp, 2.05_dp, 3.00_dp, &
154 3.00_dp, 3.00_dp, 3.00_dp, 3.00_dp, 3.00_dp, 3.00_dp, 3.00_dp, 3.00_dp, &
155 3.00_dp, 3.00_dp, 3.00_dp, 3.00_dp, 3.00_dp, 3.00_dp, &
156 3.50_dp, 3.50_dp, 3.50_dp, 3.50_dp, 3.50_dp, 3.50_dp, 3.50_dp, 3.50_dp, &
157 2.50_dp, 2.50_dp, 3.50_dp, 3.50_dp, 3.50_dp, 1.25_dp, 1.00_dp, &
158 1.05_dp, 2.05_dp, 3.00_dp, &
159 3.00_dp, 3.00_dp, 3.00_dp, 3.00_dp, 3.00_dp, 3.00_dp, 3.00_dp, 3.00_dp, &
160 3.00_dp, 3.00_dp, 3.00_dp, 3.00_dp, 3.00_dp, 3.00_dp, &
161 3.00_dp, 3.00_dp, 3.00_dp/)
166 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'xtb_parameters'
186 TYPE(xtb_atom_type),
POINTER :: param
187 CHARACTER(LEN=2),
INTENT(IN) :: element_symbol
188 CHARACTER(LEN=*),
INTENT(IN) :: parameter_file_path, parameter_file_name
189 TYPE(mp_para_env_type),
POINTER :: para_env
191 CHARACTER(len=2) :: enam, esym
192 CHARACTER(len=default_string_length) :: aname, filename
194 LOGICAL :: at_end, found
195 TYPE(cp_parser_type) :: parser
197 filename = adjustl(trim(parameter_file_path))//adjustl(trim(parameter_file_name))
204 CALL parser_get_object(parser, aname)
206 esym = element_symbol
209 IF (enam == esym)
THEN
211 CALL parser_get_object(parser, param%eta)
212 CALL parser_get_object(parser, param%xgamma)
213 CALL parser_get_object(parser, param%alpha)
214 CALL parser_get_object(parser, param%zneff)
216 CALL parser_get_object(parser, aname)
217 ia = ichar(aname(1:1))
218 IF (ia >= 49 .AND. ia <= 57)
THEN
219 CALL parser_get_object(parser, param%kpoly(i))
220 CALL parser_get_object(parser, param%kappa(i))
221 CALL parser_get_object(parser, param%hen(i))
222 CALL parser_get_object(parser, param%zeta(i))
224 param%nval(i) = ia - 48
225 SELECT CASE (aname(2:2))
235 cpabort(
"xTB PARAMETER ERROR")
247 param%typ =
"STANDARD"
248 param%symbol = element_symbol
249 param%defined = .true.
252 param%aname =
ptable(ia)%name
253 param%lmax = maxval(param%lval(1:param%nshell))
255 DO i = 1, param%nshell
257 param%natorb = param%natorb + (2*l + 1)
259 param%zeff = zval(ia)
261 esym = element_symbol
263 IF (
"X " == esym)
THEN
265 param%symbol = element_symbol
266 param%defined = .false.
274 param%defined = .false.
275 CALL cp_warn(__location__,
"xTB parameters for element "//element_symbol// &
276 " were not found in the parameter file "//adjustl(trim(filename)))
291 TYPE(xtb_atom_type),
POINTER :: param
292 CHARACTER(LEN=2),
INTENT(IN) :: element_symbol
293 TYPE(section_vals_type),
POINTER :: xtb_section
295 CHARACTER(LEN=2) :: label
296 CHARACTER(len=20*default_string_length) :: line_att
297 INTEGER :: i, ia, k, l, nshell
298 LOGICAL :: explicit, found, is_ok
299 TYPE(cp_sll_val_type),
POINTER ::
list
300 TYPE(section_vals_type),
POINTER :: ap_section
301 TYPE(val_type),
POINTER :: val
315 IF (.NOT. is_ok)
EXIT
316 CALL val_get(val, c_val=line_att)
318 READ (line_att, *) label
320 ia = ichar(label(1:1))
321 IF (ia >= 49 .AND. ia <= 57)
THEN
324 param%nval(k) = ia - 48
325 SELECT CASE (label(2:2))
335 cpabort(
"xTB PARAMETER ERROR")
338 READ (line_att, *) param%kpoly(k)
340 READ (line_att, *) param%kappa(k)
342 READ (line_att, *) param%hen(k)
344 READ (line_att, *) param%zeta(k)
350 READ (line_att, *) label
352 IF (label == element_symbol)
THEN
356 READ (line_att, *) param%eta
358 READ (line_att, *) param%xgamma
360 READ (line_att, *) param%alpha
362 READ (line_att, *) param%zneff
364 READ (line_att, *) label
366 ia = ichar(label(1:1))
367 cpassert((ia >= 49 .AND. ia <= 57))
368 param%nval(k) = ia - 48
369 SELECT CASE (label(2:2))
379 cpabort(
"xTB PARAMETER ERROR")
382 READ (line_att, *) param%kpoly(k)
384 READ (line_att, *) param%kappa(k)
386 READ (line_att, *) param%hen(k)
388 READ (line_att, *) param%zeta(k)
394 param%typ =
"STANDARD"
395 param%symbol = element_symbol
396 param%defined = .true.
399 param%aname =
ptable(ia)%name
400 param%lmax = maxval(param%lval(1:param%nshell))
402 param%nshell = nshell
403 DO i = 1, param%nshell
405 param%natorb = param%natorb + (2*l + 1)
407 param%zeff = zval(ia)
419 TYPE(xtb_atom_type),
POINTER :: param
421 INTEGER :: i, is, l, na
422 REAL(kind=
dp),
DIMENSION(5) :: kp
424 IF (param%defined)
THEN
428 DO is = 1, param%nshell
439 param%electronegativity = eneg(i)
441 param%rcov = crad(i)*
bohr
443 param%occupation(:) = occupation(:, i)
445 IF (abs(param%zeff - sum(param%occupation)) > 1.e-10_dp)
THEN
446 CALL cp_abort(__location__,
"Element <"//trim(param%aname)//
"> has inconsistent shell occupations")
449 param%hen = param%hen/
evolt
451 param%xgamma = 0.1_dp*param%xgamma
452 param%kpoly(:) = 0.01_dp*param%kpoly(:)
453 param%kappa(:) = 0.1_dp*param%kappa(:)
455 param%xgamma = -2.0_dp*param%xgamma
457 kp(:) = param%kappa(:)
458 param%kappa(:) = 0.0_dp
459 DO is = 1, param%nshell
461 IF (param%kappa(l + 1) == 0.0_dp)
THEN
462 param%kappa(l + 1) = kp(is)
464 cpassert(abs(param%kappa(l + 1) - kp(is)) < 1.e-10_dp)
468 IF (param%kx < -10._dp)
THEN
470 SELECT CASE (param%z)
474 param%kx = 0.1_dp*0.381742_dp
476 param%kx = 0.1_dp*0.321944_dp
478 param%kx = 0.1_dp*0.220000_dp
482 param%chmax = clmt(i)
495 TYPE(xtb_atom_type),
POINTER :: param
496 TYPE(gto_basis_set_type),
POINTER :: gto_basis_set
497 INTEGER,
INTENT(IN) :: ngauss
499 CHARACTER(LEN=6),
DIMENSION(:),
POINTER :: symbol
501 INTEGER,
DIMENSION(:),
POINTER :: lq, nq
502 REAL(kind=
dp),
DIMENSION(:),
POINTER :: zet
503 TYPE(sto_basis_set_type),
POINTER :: sto_basis_set
505 IF (
ASSOCIATED(param))
THEN
506 IF (param%defined)
THEN
507 NULLIFY (sto_basis_set)
509 nshell = param%nshell
511 ALLOCATE (symbol(1:nshell))
514 SELECT CASE (param%lval(i))
516 WRITE (symbol(i),
'(I1,A1)') param%nval(i),
"S"
518 WRITE (symbol(i),
'(I1,A1)') param%nval(i),
"P"
520 WRITE (symbol(i),
'(I1,A1)') param%nval(i),
"D"
522 WRITE (symbol(i),
'(I1,A1)') param%nval(i),
"F"
524 cpabort(
'BASIS SET OUT OF RANGE (lval)')
529 ALLOCATE (nq(nshell), lq(nshell), zet(nshell))
530 nq(1:nshell) = param%nval(1:nshell)
531 lq(1:nshell) = param%lval(1:nshell)
532 zet(1:nshell) = param%zeta(1:nshell)
533 CALL set_sto_basis_set(sto_basis_set, name=param%aname, nshell=nshell, symbol=symbol, &
534 nq=nq, lq=lq, zet=zet)
540 DEALLOCATE (symbol, nq, lq, zet)
544 cpabort(
"The pointer param is not associated")
558 INTEGER,
INTENT(IN) :: za, zb
559 TYPE(xtb_control_type),
INTENT(IN),
POINTER :: xtb_control
568 IF (xtb_control%kab_nval .GT. 0)
THEN
569 DO j = 1, xtb_control%kab_nval
570 IF ((za == xtb_control%kab_types(1, j) .AND. &
571 zb == xtb_control%kab_types(2, j)) .OR. &
572 (za == xtb_control%kab_types(2, j) .AND. &
573 zb == xtb_control%kab_types(1, j)))
THEN
575 kab = xtb_control%kab_vals(j)
581 IF (.NOT. custom)
THEN
582 IF (za == 1 .OR. zb == 1)
THEN
599 ELSEIF (za == 5 .OR. zb == 5)
THEN
606 ELSEIF (za == 7 .OR. zb == 7)
THEN
615 ELSEIF (za > 20 .AND. za < 30)
THEN
617 IF (zb > 20 .AND. zb < 30)
THEN
620 ELSEIF ((zb > 38 .AND. zb < 48) .OR. (zb > 56 .AND. zb < 80))
THEN
622 kab = 0.50_dp*(1.20_dp + 1.10_dp)
624 ELSEIF ((za > 38 .AND. za < 48) .OR. (za > 56 .AND. za < 80))
THEN
626 IF (zb > 20 .AND. zb < 30)
THEN
628 kab = 0.50_dp*(1.20_dp + 1.10_dp)
629 ELSEIF ((zb > 38 .AND. zb < 48) .OR. (zb > 56 .AND. zb < 80))
THEN
subroutine, public deallocate_sto_basis_set(sto_basis_set)
...
subroutine, public allocate_sto_basis_set(sto_basis_set)
...
subroutine, public create_gto_from_sto_basis(sto_basis_set, gto_basis_set, ngauss, ortho)
...
subroutine, public set_sto_basis_set(sto_basis_set, name, nshell, symbol, nq, lq, zet)
...
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
Utility routines to read data from files. Kept as close as possible to the old parser because.
subroutine, public parser_get_next_line(parser, nline, at_end)
Read the next input line and broadcast the input information. Skip (nline-1) lines and skip also all ...
Utility routines to read data from files. Kept as close as possible to the old parser because.
subroutine, public parser_release(parser)
releases the parser
subroutine, public parser_create(parser, file_name, unit_nr, para_env, end_section_label, separator_chars, comment_char, continuation_char, quote_char, section_char, parse_white_lines, initial_variables, apply_preprocessing)
Start a parser run. Initial variables allow to @SET stuff before opening the file.
Defines the basic variable types.
integer, parameter, public dp
integer, parameter, public default_string_length
An array-based list which grows on demand. When the internal array is full, a new array of twice the ...
Interface to the message passing library MPI.
Periodic Table related data definitions.
type(atom), dimension(0:nelem), public ptable
subroutine, public get_ptable_info(symbol, number, amass, ielement, covalent_radius, metallic_radius, vdw_radius, found)
Pass information about the kind given the element symbol.
Definition of physical constants:
real(kind=dp), parameter, public evolt
real(kind=dp), parameter, public bohr
Utilities for string manipulations.
subroutine, public remove_word(string)
remove a word from a string (words are separated by white spaces)
elemental subroutine, public uppercase(string)
Convert all lower case characters in a string to upper case.
subroutine, public xtb_parameters_set(param)
Read atom parameters for xTB Hamiltonian from input file.
subroutine, public xtb_parameters_init(param, element_symbol, parameter_file_path, parameter_file_name, para_env)
...
subroutine, public init_xtb_basis(param, gto_basis_set, ngauss)
...
subroutine, public xtb_parameters_read(param, element_symbol, xtb_section)
Read atom parameters for xTB Hamiltonian from input file.
real(kind=dp) function, public xtb_set_kab(za, zb, xtb_control)
...
Definition of the xTB parameter types.