Line data Source code
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 Definition of the xTB parameter types.
10 : !> \author JGH (10.2018)
11 : ! **************************************************************************************************
12 : ! To be done:
13 : ! 1) Ewald defaults options for GMAX, ALPHA, RCUT
14 : ! 2) QM/MM debugging of forces -- done
15 : ! 3) Periodic displacement field (debugging)
16 : ! 4) Check for RTP and EMD
17 : ! 5) Wannier localization
18 : ! 6) Charge Mixing methods: Broyden/Pulay (more debugging needed, also add to DFTB)
19 : ! **************************************************************************************************
20 : MODULE xtb_types
21 :
22 : USE cp_log_handling, ONLY: cp_get_default_logger,&
23 : cp_logger_type
24 : USE cp_output_handling, ONLY: cp_p_file,&
25 : cp_print_key_finished_output,&
26 : cp_print_key_should_output,&
27 : cp_print_key_unit_nr
28 : USE input_section_types, ONLY: section_vals_type
29 : USE kinds, ONLY: default_string_length,&
30 : dp
31 : #include "./base/base_uses.f90"
32 :
33 : IMPLICIT NONE
34 :
35 : PRIVATE
36 :
37 : ! *** Global parameters ***
38 :
39 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'xtb_types'
40 :
41 : ! **************************************************************************************************
42 : TYPE xtb_atom_type
43 : ! PRIVATE
44 : CHARACTER(LEN=default_string_length) :: typ = ""
45 : CHARACTER(LEN=default_string_length) :: aname = ""
46 : CHARACTER(LEN=2) :: symbol = ""
47 : LOGICAL :: defined = .FALSE.
48 : INTEGER :: z = -1 !atomic number
49 : REAL(KIND=dp) :: zeff = -1.0_dp !effective core charge
50 : INTEGER :: natorb = -1 !number of orbitals
51 : INTEGER :: lmax = -1 !max angular momentum
52 : !
53 : REAL(KIND=dp) :: rcut = -1.0_dp !cutoff radius for sr-Coulomb
54 : REAL(KIND=dp) :: rcov = -1.0_dp !covalent radius
55 : REAL(KIND=dp) :: electronegativity = -1.0_dp !electronegativity
56 : !
57 : REAL(KIND=dp) :: kx = -1.0_dp !scaling for halogen term
58 : !
59 : REAL(KIND=dp) :: eta = -1.0_dp !Atomic Hubbard parameter
60 : REAL(KIND=dp) :: xgamma = -1.0_dp !charge derivative of eta
61 : REAL(KIND=dp) :: alpha = -1.0_dp !exponential scaling parameter for repulsion potential
62 : REAL(KIND=dp) :: zneff = -1.0_dp !effective core charge for repulsion potential
63 : ! shell specific parameters
64 : INTEGER :: nshell = -1 !number of orbital shells
65 : INTEGER, DIMENSION(5) :: nval = -1 ! n-quantum number of shell i
66 : INTEGER, DIMENSION(5) :: lval = -1 ! l-quantum number of shell i
67 : INTEGER, DIMENSION(5) :: occupation = -1 ! occupation of shell i
68 : REAL(KIND=dp), DIMENSION(5) :: kpoly = -1.0_dp
69 : REAL(KIND=dp), DIMENSION(5) :: kappa = -1.0_dp
70 : REAL(KIND=dp), DIMENSION(5) :: hen = -1.0_dp
71 : REAL(KIND=dp), DIMENSION(5) :: zeta = -1.0_dp
72 : ! gfn0 params
73 : REAL(KIND=dp) :: en = -1.0_dp
74 : REAL(KIND=dp) :: kqat2 = -1.0_dp
75 : REAL(KIND=dp), DIMENSION(5) :: kq = -1.0_dp
76 : REAL(KIND=dp), DIMENSION(5) :: kcn = -1.0_dp
77 : ! charge equilibration parameter gfn0
78 : REAL(KIND=dp) :: xi = -1.0_dp
79 : REAL(KIND=dp) :: kappa0 = -1.0_dp
80 : REAL(KIND=dp) :: alpg = -1.0_dp
81 : ! AO to shell pointer
82 : INTEGER, DIMENSION(25) :: nao = -1, lao = -1
83 : ! Upper limit of Mulliken charge
84 : REAL(KIND=dp) :: chmax = -1.0_dp
85 : END TYPE xtb_atom_type
86 :
87 : ! *** Public data types ***
88 :
89 : PUBLIC :: xtb_atom_type, get_xtb_atom_param, set_xtb_atom_param, write_xtb_atom_param
90 : PUBLIC :: allocate_xtb_atom_param, deallocate_xtb_atom_param
91 :
92 : CONTAINS
93 :
94 : ! **************************************************************************************************
95 : !> \brief ...
96 : !> \param xtb_parameter ...
97 : ! **************************************************************************************************
98 2000 : SUBROUTINE allocate_xtb_atom_param(xtb_parameter)
99 :
100 : TYPE(xtb_atom_type), POINTER :: xtb_parameter
101 :
102 2000 : IF (ASSOCIATED(xtb_parameter)) &
103 0 : CALL deallocate_xtb_atom_param(xtb_parameter)
104 :
105 212000 : ALLOCATE (xtb_parameter)
106 :
107 : xtb_parameter%defined = .FALSE.
108 2000 : xtb_parameter%aname = ""
109 2000 : xtb_parameter%symbol = ""
110 2000 : xtb_parameter%typ = "NONE"
111 : xtb_parameter%z = -1
112 : xtb_parameter%zeff = -1.0_dp
113 2000 : xtb_parameter%natorb = 0
114 : xtb_parameter%lmax = -1
115 2000 : xtb_parameter%rcut = 0.0_dp
116 2000 : xtb_parameter%rcov = 0.0_dp
117 2000 : xtb_parameter%electronegativity = 0.0_dp
118 2000 : xtb_parameter%kx = -100.0_dp
119 2000 : xtb_parameter%eta = 0.0_dp
120 2000 : xtb_parameter%xgamma = 0.0_dp
121 2000 : xtb_parameter%alpha = 0.0_dp
122 2000 : xtb_parameter%zneff = 0.0_dp
123 2000 : xtb_parameter%nshell = 0
124 12000 : xtb_parameter%nval = 0
125 12000 : xtb_parameter%lval = 0
126 12000 : xtb_parameter%occupation = 0
127 12000 : xtb_parameter%kpoly = 0.0_dp
128 12000 : xtb_parameter%kappa = 0.0_dp
129 12000 : xtb_parameter%hen = 0.0_dp
130 12000 : xtb_parameter%zeta = 0.0_dp
131 2000 : xtb_parameter%en = 0.0_dp
132 2000 : xtb_parameter%kqat2 = 0.0_dp
133 12000 : xtb_parameter%kq = 0.0_dp
134 12000 : xtb_parameter%kcn = 0.0_dp
135 2000 : xtb_parameter%xi = 0.0_dp
136 2000 : xtb_parameter%kappa0 = 0.0_dp
137 2000 : xtb_parameter%alpg = 0.0_dp
138 52000 : xtb_parameter%nao = 0
139 52000 : xtb_parameter%lao = 0
140 2000 : xtb_parameter%chmax = 0.0_dp
141 :
142 2000 : END SUBROUTINE allocate_xtb_atom_param
143 :
144 : ! **************************************************************************************************
145 : !> \brief ...
146 : !> \param xtb_parameter ...
147 : ! **************************************************************************************************
148 2000 : SUBROUTINE deallocate_xtb_atom_param(xtb_parameter)
149 :
150 : TYPE(xtb_atom_type), POINTER :: xtb_parameter
151 :
152 2000 : CPASSERT(ASSOCIATED(xtb_parameter))
153 2000 : DEALLOCATE (xtb_parameter)
154 :
155 2000 : END SUBROUTINE deallocate_xtb_atom_param
156 :
157 : ! **************************************************************************************************
158 : !> \brief ...
159 : !> \param xtb_parameter ...
160 : !> \param symbol ...
161 : !> \param aname ...
162 : !> \param typ ...
163 : !> \param defined ...
164 : !> \param z ...
165 : !> \param zeff ...
166 : !> \param natorb ...
167 : !> \param lmax ...
168 : !> \param nao ...
169 : !> \param lao ...
170 : !> \param rcut ...
171 : !> \param rcov ...
172 : !> \param kx ...
173 : !> \param eta ...
174 : !> \param xgamma ...
175 : !> \param alpha ...
176 : !> \param zneff ...
177 : !> \param nshell ...
178 : !> \param nval ...
179 : !> \param lval ...
180 : !> \param kpoly ...
181 : !> \param kappa ...
182 : !> \param hen ...
183 : !> \param zeta ...
184 : !> \param xi ...
185 : !> \param kappa0 ...
186 : !> \param alpg ...
187 : !> \param occupation ...
188 : !> \param electronegativity ...
189 : !> \param chmax ...
190 : !> \param en ...
191 : !> \param kqat2 ...
192 : !> \param kcn ...
193 : !> \param kq ...
194 : ! **************************************************************************************************
195 43926367 : SUBROUTINE get_xtb_atom_param(xtb_parameter, symbol, aname, typ, defined, z, zeff, natorb, lmax, nao, lao, &
196 : rcut, rcov, kx, eta, xgamma, alpha, zneff, nshell, nval, lval, kpoly, kappa, &
197 : hen, zeta, xi, kappa0, alpg, occupation, electronegativity, chmax, &
198 : en, kqat2, kcn, kq)
199 :
200 : TYPE(xtb_atom_type), POINTER :: xtb_parameter
201 : CHARACTER(LEN=2), INTENT(OUT), OPTIONAL :: symbol
202 : CHARACTER(LEN=default_string_length), &
203 : INTENT(OUT), OPTIONAL :: aname, typ
204 : LOGICAL, INTENT(OUT), OPTIONAL :: defined
205 : INTEGER, INTENT(OUT), OPTIONAL :: z
206 : REAL(KIND=dp), INTENT(OUT), OPTIONAL :: zeff
207 : INTEGER, INTENT(OUT), OPTIONAL :: natorb, lmax
208 : INTEGER, DIMENSION(25), INTENT(OUT), OPTIONAL :: nao, lao
209 : REAL(KIND=dp), INTENT(OUT), OPTIONAL :: rcut, rcov, kx, eta, xgamma, alpha, zneff
210 : INTEGER, INTENT(OUT), OPTIONAL :: nshell
211 : INTEGER, DIMENSION(5), INTENT(OUT), OPTIONAL :: nval, lval
212 : REAL(KIND=dp), DIMENSION(5), INTENT(OUT), OPTIONAL :: kpoly, kappa, hen, zeta
213 : REAL(KIND=dp), INTENT(OUT), OPTIONAL :: xi, kappa0, alpg
214 : INTEGER, DIMENSION(5), INTENT(OUT), OPTIONAL :: occupation
215 : REAL(KIND=dp), INTENT(OUT), OPTIONAL :: electronegativity, chmax, en, kqat2
216 : REAL(KIND=dp), DIMENSION(5), INTENT(OUT), OPTIONAL :: kcn, kq
217 :
218 43926367 : CPASSERT(ASSOCIATED(xtb_parameter))
219 :
220 43926367 : IF (PRESENT(symbol)) symbol = xtb_parameter%symbol
221 43926367 : IF (PRESENT(aname)) aname = xtb_parameter%aname
222 43926367 : IF (PRESENT(typ)) typ = xtb_parameter%typ
223 43926367 : IF (PRESENT(defined)) defined = xtb_parameter%defined
224 43926367 : IF (PRESENT(z)) z = xtb_parameter%z
225 43926367 : IF (PRESENT(zeff)) zeff = xtb_parameter%zeff
226 43926367 : IF (PRESENT(natorb)) natorb = xtb_parameter%natorb
227 43926367 : IF (PRESENT(lmax)) lmax = xtb_parameter%lmax
228 113527017 : IF (PRESENT(nao)) nao = xtb_parameter%nao
229 283769717 : IF (PRESENT(lao)) lao = xtb_parameter%lao
230 : !
231 43926367 : IF (PRESENT(rcut)) rcut = xtb_parameter%rcut
232 43926367 : IF (PRESENT(rcov)) rcov = xtb_parameter%rcov
233 43926367 : IF (PRESENT(kx)) kx = xtb_parameter%kx
234 43926367 : IF (PRESENT(electronegativity)) electronegativity = xtb_parameter%electronegativity
235 43926367 : IF (PRESENT(eta)) eta = xtb_parameter%eta
236 43926367 : IF (PRESENT(xgamma)) xgamma = xtb_parameter%xgamma
237 43926367 : IF (PRESENT(alpha)) alpha = xtb_parameter%alpha
238 43926367 : IF (PRESENT(zneff)) zneff = xtb_parameter%zneff
239 43926367 : IF (PRESENT(nshell)) nshell = xtb_parameter%nshell
240 43926367 : IF (PRESENT(nval)) nval = xtb_parameter%nval
241 44103527 : IF (PRESENT(lval)) lval = xtb_parameter%lval
242 44212827 : IF (PRESENT(occupation)) occupation = xtb_parameter%occupation
243 57634687 : IF (PRESENT(kpoly)) kpoly = xtb_parameter%kpoly
244 120020827 : IF (PRESENT(kappa)) kappa = xtb_parameter%kappa
245 48276317 : IF (PRESENT(hen)) hen = xtb_parameter%hen
246 43987347 : IF (PRESENT(zeta)) zeta = xtb_parameter%zeta
247 43926367 : IF (PRESENT(chmax)) chmax = xtb_parameter%chmax
248 43926367 : IF (PRESENT(xi)) xi = xtb_parameter%xi
249 43926367 : IF (PRESENT(kappa0)) kappa0 = xtb_parameter%kappa0
250 43926367 : IF (PRESENT(alpg)) alpg = xtb_parameter%alpg
251 43926367 : IF (PRESENT(en)) en = xtb_parameter%en
252 43926367 : IF (PRESENT(kqat2)) kqat2 = xtb_parameter%kqat2
253 43960917 : IF (PRESENT(kcn)) kcn = xtb_parameter%kcn
254 43960917 : IF (PRESENT(kq)) kq = xtb_parameter%kq
255 :
256 43926367 : END SUBROUTINE get_xtb_atom_param
257 :
258 : ! **************************************************************************************************
259 : !> \brief ...
260 : !> \param xtb_parameter ...
261 : !> \param aname ...
262 : !> \param typ ...
263 : !> \param defined ...
264 : !> \param z ...
265 : !> \param zeff ...
266 : !> \param natorb ...
267 : !> \param lmax ...
268 : !> \param nao ...
269 : !> \param lao ...
270 : !> \param rcut ...
271 : !> \param rcov ...
272 : !> \param kx ...
273 : !> \param eta ...
274 : !> \param xgamma ...
275 : !> \param alpha ...
276 : !> \param zneff ...
277 : !> \param nshell ...
278 : !> \param nval ...
279 : !> \param lval ...
280 : !> \param kpoly ...
281 : !> \param kappa ...
282 : !> \param hen ...
283 : !> \param zeta ...
284 : !> \param xi ...
285 : !> \param kappa0 ...
286 : !> \param alpg ...
287 : !> \param electronegativity ...
288 : !> \param occupation ...
289 : !> \param chmax ...
290 : !> \param en ...
291 : !> \param kqat2 ...
292 : !> \param kcn ...
293 : !> \param kq ...
294 : ! **************************************************************************************************
295 0 : SUBROUTINE set_xtb_atom_param(xtb_parameter, aname, typ, defined, z, zeff, natorb, lmax, nao, lao, &
296 : rcut, rcov, kx, eta, xgamma, alpha, zneff, nshell, nval, lval, kpoly, kappa, &
297 : hen, zeta, xi, kappa0, alpg, electronegativity, occupation, chmax, &
298 : en, kqat2, kcn, kq)
299 :
300 : TYPE(xtb_atom_type), POINTER :: xtb_parameter
301 : CHARACTER(LEN=default_string_length), INTENT(IN), &
302 : OPTIONAL :: aname, typ
303 : LOGICAL, INTENT(IN), OPTIONAL :: defined
304 : INTEGER, INTENT(IN), OPTIONAL :: z
305 : REAL(KIND=dp), INTENT(IN), OPTIONAL :: zeff
306 : INTEGER, INTENT(IN), OPTIONAL :: natorb, lmax
307 : INTEGER, DIMENSION(25), INTENT(IN), OPTIONAL :: nao, lao
308 : REAL(KIND=dp), INTENT(IN), OPTIONAL :: rcut, rcov, kx, eta, xgamma, alpha, zneff
309 : INTEGER, INTENT(IN), OPTIONAL :: nshell
310 : INTEGER, DIMENSION(5), INTENT(IN), OPTIONAL :: nval, lval
311 : REAL(KIND=dp), DIMENSION(5), INTENT(IN), OPTIONAL :: kpoly, kappa, hen, zeta
312 : REAL(KIND=dp), INTENT(IN), OPTIONAL :: xi, kappa0, alpg, electronegativity
313 : INTEGER, DIMENSION(5), INTENT(IN), OPTIONAL :: occupation
314 : REAL(KIND=dp), INTENT(IN), OPTIONAL :: chmax, en, kqat2
315 : REAL(KIND=dp), DIMENSION(5), INTENT(IN), OPTIONAL :: kcn, kq
316 :
317 0 : CPASSERT(ASSOCIATED(xtb_parameter))
318 :
319 0 : IF (PRESENT(aname)) xtb_parameter%aname = aname
320 0 : IF (PRESENT(typ)) xtb_parameter%typ = typ
321 0 : IF (PRESENT(defined)) xtb_parameter%defined = defined
322 0 : IF (PRESENT(z)) xtb_parameter%z = z
323 0 : IF (PRESENT(zeff)) xtb_parameter%zeff = zeff
324 0 : IF (PRESENT(natorb)) xtb_parameter%natorb = natorb
325 0 : IF (PRESENT(lmax)) xtb_parameter%lmax = lmax
326 0 : IF (PRESENT(nao)) xtb_parameter%nao = nao
327 0 : IF (PRESENT(lao)) xtb_parameter%lao = lao
328 : !
329 0 : IF (PRESENT(rcut)) xtb_parameter%rcut = rcut
330 0 : IF (PRESENT(rcov)) xtb_parameter%rcov = rcov
331 0 : IF (PRESENT(kx)) xtb_parameter%kx = kx
332 0 : IF (PRESENT(electronegativity)) xtb_parameter%electronegativity = electronegativity
333 0 : IF (PRESENT(eta)) xtb_parameter%eta = eta
334 0 : IF (PRESENT(xgamma)) xtb_parameter%xgamma = xgamma
335 0 : IF (PRESENT(alpha)) xtb_parameter%alpha = alpha
336 0 : IF (PRESENT(zneff)) xtb_parameter%zneff = zneff
337 0 : IF (PRESENT(nshell)) xtb_parameter%nshell = nshell
338 0 : IF (PRESENT(nval)) xtb_parameter%nval = nval
339 0 : IF (PRESENT(lval)) xtb_parameter%lval = lval
340 0 : IF (PRESENT(occupation)) xtb_parameter%occupation = occupation
341 0 : IF (PRESENT(kpoly)) xtb_parameter%kpoly = kpoly
342 0 : IF (PRESENT(kappa)) xtb_parameter%kappa = kappa
343 0 : IF (PRESENT(hen)) xtb_parameter%hen = hen
344 0 : IF (PRESENT(zeta)) xtb_parameter%zeta = zeta
345 0 : IF (PRESENT(chmax)) xtb_parameter%chmax = chmax
346 : !
347 0 : IF (PRESENT(xi)) xtb_parameter%xi = xi
348 0 : IF (PRESENT(kappa0)) xtb_parameter%kappa0 = kappa0
349 0 : IF (PRESENT(alpg)) xtb_parameter%alpg = alpg
350 0 : IF (PRESENT(en)) xtb_parameter%en = en
351 0 : IF (PRESENT(kqat2)) xtb_parameter%kqat2 = kqat2
352 0 : IF (PRESENT(kcn)) xtb_parameter%kcn = kcn
353 0 : IF (PRESENT(kq)) xtb_parameter%kq = kq
354 :
355 0 : END SUBROUTINE set_xtb_atom_param
356 :
357 : ! **************************************************************************************************
358 : !> \brief ...
359 : !> \param xtb_parameter ...
360 : !> \param gfn_type ...
361 : !> \param subsys_section ...
362 : ! **************************************************************************************************
363 2000 : SUBROUTINE write_xtb_atom_param(xtb_parameter, gfn_type, subsys_section)
364 :
365 : TYPE(xtb_atom_type), POINTER :: xtb_parameter
366 : INTEGER, INTENT(IN) :: gfn_type
367 : TYPE(section_vals_type), POINTER :: subsys_section
368 :
369 : CHARACTER(LEN=default_string_length) :: aname, bb
370 : INTEGER :: i, io_unit, m, natorb, nshell
371 : INTEGER, DIMENSION(5) :: lval, nval, occupation
372 : LOGICAL :: defined
373 : REAL(dp) :: zeff
374 : REAL(KIND=dp) :: alpha, en, eta, xgamma, zneff
375 : REAL(KIND=dp), DIMENSION(5) :: hen, kappa, kpoly, zeta
376 : TYPE(cp_logger_type), POINTER :: logger
377 :
378 2000 : NULLIFY (logger)
379 2000 : logger => cp_get_default_logger()
380 2000 : IF (ASSOCIATED(xtb_parameter) .AND. &
381 : BTEST(cp_print_key_should_output(logger%iter_info, subsys_section, &
382 : "PRINT%KINDS/POTENTIAL"), cp_p_file)) THEN
383 :
384 : io_unit = cp_print_key_unit_nr(logger, subsys_section, "PRINT%KINDS", &
385 0 : extension=".Log")
386 :
387 0 : IF (io_unit > 0) THEN
388 0 : SELECT CASE (gfn_type)
389 : CASE (0)
390 0 : CPABORT("gfn_type = 0 missing code")
391 : CASE (1)
392 0 : CALL get_xtb_atom_param(xtb_parameter, aname=aname, defined=defined, zeff=zeff, natorb=natorb)
393 0 : CALL get_xtb_atom_param(xtb_parameter, nshell=nshell, lval=lval, nval=nval, occupation=occupation)
394 0 : CALL get_xtb_atom_param(xtb_parameter, kpoly=kpoly, kappa=kappa, hen=hen, zeta=zeta)
395 0 : CALL get_xtb_atom_param(xtb_parameter, electronegativity=en, xgamma=xgamma, eta=eta, alpha=alpha, zneff=zneff)
396 :
397 0 : bb = " "
398 0 : WRITE (UNIT=io_unit, FMT="(/,A,T67,A14)") " xTB parameters: ", TRIM(aname)
399 0 : IF (defined) THEN
400 0 : m = 5 - nshell
401 0 : WRITE (UNIT=io_unit, FMT="(T16,A,T71,F10.2)") "Effective core charge:", zeff
402 0 : WRITE (UNIT=io_unit, FMT="(T16,A,T71,I10)") "Number of orbitals:", natorb
403 0 : WRITE (UNIT=io_unit, FMT="(T16,A,T41,A,5(A4,I1,I2,A1))") "Basis set [nl]", bb(1:8*m), &
404 0 : (" [", nval(i), lval(i), "]", i=1, nshell)
405 0 : WRITE (UNIT=io_unit, FMT="(T16,A,T41,A,5F8.3)") "Slater Exponent", bb(1:8*m), (zeta(i), i=1, nshell)
406 0 : WRITE (UNIT=io_unit, FMT="(T16,A,T41,A,5I8)") "Ref. occupation", bb(1:8*m), (occupation(i), i=1, nshell)
407 0 : WRITE (UNIT=io_unit, FMT="(T16,A,T41,A,5F8.3)") "Energy levels [au]", bb(1:8*m), (hen(i), i=1, nshell)
408 0 : WRITE (UNIT=io_unit, FMT="(T16,A,T41,A,5F8.3)") "Kpoly", bb(1:8*m), (kpoly(i), i=1, nshell)
409 0 : WRITE (UNIT=io_unit, FMT="(T16,A,T71,F10.3)") "Electronegativity", en
410 0 : WRITE (UNIT=io_unit, FMT="(T16,A,T71,F10.3)") "Mataga-Nishimoto constant (eta)", eta
411 0 : WRITE (UNIT=io_unit, FMT="(T16,A,T41,A,5F8.3)") "Mataga-Nishimoto scaling kappa", bb(1:8*m), &
412 0 : (kappa(i), i=1, nshell)
413 0 : WRITE (UNIT=io_unit, FMT="(T16,A,T71,F10.3)") "3rd Order constant", xgamma
414 0 : WRITE (UNIT=io_unit, FMT="(T16,A,T61,2F10.3)") "Repulsion potential [Z,alpha]", zneff, alpha
415 : ELSE
416 0 : WRITE (UNIT=io_unit, FMT="(T55,A)") "Parameters are not defined"
417 : END IF
418 : CASE (2)
419 0 : CPABORT("gfn_type = 2 not yet defined")
420 : END SELECT
421 : END IF
422 0 : CALL cp_print_key_finished_output(io_unit, logger, subsys_section, "PRINT%KINDS")
423 : END IF
424 :
425 2000 : END SUBROUTINE write_xtb_atom_param
426 :
427 0 : END MODULE xtb_types
428 :
|