Octopus
vdw_ts.F90
Go to the documentation of this file.
1!! Copyright (C) 2015 X. Andrade, R. A. DiStasio Jr., B. Santra
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
26module vdw_ts_oct_m
27 use atom_oct_m
28 use debug_oct_m
29 use global_oct_m
31 use io_oct_m
33 use ions_oct_m
34 use, intrinsic :: iso_fortran_env
38 use mesh_oct_m
39 use mpi_oct_m
41 use parser_oct_m
43 use ps_oct_m
44 use space_oct_m
47 use unit_oct_m
49
50 implicit none
51
52 private
53
54 public :: &
55 vdw_ts_t, &
57 vdw_ts_end, &
61
62 type vdw_ts_t
63 private
64 real(real64), allocatable :: c6free(:)
65 real(real64), allocatable :: dpfree(:)
66 real(real64), allocatable :: r0free(:)
67 real(real64), allocatable :: c6abfree(:, :)
68 real(real64), allocatable :: volfree(:)
69 real(real64), allocatable :: c6ab(:, :)
70 real(real64) :: cutoff
71 real(real64) :: damping
72 real(real64) :: sr
73
74 real(real64), allocatable :: derivative_coeff(:)
75 end type vdw_ts_t
76
77contains
78
79 subroutine vdw_ts_init(this, namespace, ions)
80 type(vdw_ts_t), intent(out) :: this
81 type(namespace_t), intent(in) :: namespace
82 type(ions_t), intent(in) :: ions
83
84 integer :: ispecies, jspecies
85 real(real64) :: num, den
86
87 push_sub(vdw_ts_init)
88
89 !%Variable VDW_TS_cutoff
90 !%Type float
91 !%Default 10.0
92 !%Section Hamiltonian::XC
93 !%Description
94 !% Set the value of the cutoff (unit of length) for the VDW correction in periodic system
95 !% in the Tkatchenko and Scheffler (vdw_ts) scheme only.
96 !%End
97 call parse_variable(namespace, 'VDW_TS_cutoff', 10.0_real64, this%cutoff, units_inp%length)
98
99
100 !%Variable VDW_TS_damping
101 !%Type float
102 !%Default 20.0
103 !%Section Hamiltonian::XC
104 !%Description
105 !% Set the value of the damping function (in unit of 1/length) steepness for the VDW correction in the
106 !% Tkatchenko-Scheffler scheme. See Equation (12) of Phys. Rev. Lett. 102 073005 (2009).
107 !%End
108 call parse_variable(namespace, 'VDW_TS_damping', 20.0_real64, this%damping, units_inp%length**(-1))
109
110 !%Variable VDW_TS_sr
111 !%Type float
112 !%Default 0.94
113 !%Section Hamiltonian::XC
114 !%Description
115 !% Set the value of the sr parameter in the damping function of the VDW correction in the
116 !% Tkatchenko-Scheffler scheme. See Equation (12) of Phys. Rev. Lett. 102 073005 (2009).
117 !% This parameter depends on the xc functional used.
118 !% The default value is 0.94, which holds for PBE. For PBE0, a value of 0.96 should be used.
119 !%End
120 call parse_variable(namespace, 'VDW_TS_sr', 0.94_real64, this%sr)
121
122
123 safe_allocate(this%c6free(1:ions%nspecies))
124 safe_allocate(this%dpfree(1:ions%nspecies))
125 safe_allocate(this%r0free(1:ions%nspecies))
126 safe_allocate(this%volfree(1:ions%nspecies))
127 safe_allocate(this%c6abfree(1:ions%nspecies, 1:ions%nspecies))
128 safe_allocate(this%c6ab(1:ions%natoms, 1:ions%natoms))
129 safe_allocate(this%derivative_coeff(1:ions%natoms))
130
131 do ispecies = 1, ions%nspecies
132 call get_vdw_param(namespace, ions%species(ispecies)%s%get_label(), &
133 this%c6free(ispecies), this%dpfree(ispecies), this%r0free(ispecies))
134 select type(spec=>ions%species(ispecies)%s)
135 type is(pseudopotential_t)
136 this%volfree(ispecies) = ps_density_volume(spec%ps, namespace)
137 class default
138 assert(.false.)
139 end select
140 end do
141
142 do ispecies = 1, ions%nspecies
143 do jspecies = 1, ions%nspecies
144 num = m_two*this%c6free(ispecies)*this%c6free(jspecies)
145 den = (this%dpfree(jspecies)/this%dpfree(ispecies))*this%c6free(ispecies) &
146 + (this%dpfree(ispecies)/this%dpfree(jspecies))*this%c6free(jspecies)
147 this%c6abfree(ispecies, jspecies) = num/den
148 end do
149 end do
150 pop_sub(vdw_ts_init)
151 end subroutine vdw_ts_init
152
153 !------------------------------------------
154
155 subroutine vdw_ts_end(this)
156 type(vdw_ts_t), intent(inout) :: this
158 push_sub(vdw_ts_end)
160 safe_deallocate_a(this%c6free)
161 safe_deallocate_a(this%dpfree)
162 safe_deallocate_a(this%r0free)
163 safe_deallocate_a(this%volfree)
164 safe_deallocate_a(this%c6abfree)
165 safe_deallocate_a(this%c6ab)
166 safe_deallocate_a(this%derivative_coeff)
168 pop_sub(vdw_ts_end)
169 end subroutine vdw_ts_end
170
171 !------------------------------------------
173 subroutine vdw_ts_calculate(this, namespace, space, latt, atom, natoms, pos, mesh, nspin, density, energy, potential, force)
174 type(vdw_ts_t), intent(inout) :: this
175 type(namespace_t), intent(in) :: namespace
176 class(space_t), intent(in) :: space
177 type(lattice_vectors_t), intent(in) :: latt
178 type(atom_t), intent(in) :: atom(:)
179 integer, intent(in) :: natoms
180 real(real64), intent(in) :: pos(1:space%dim,1:natoms)
181 class(mesh_t), intent(in) :: mesh
182 integer, intent(in) :: nspin
183 real(real64), intent(in) :: density(:, :)
184 real(real64), intent(out) :: energy
185 real(real64), contiguous, intent(out) :: potential(:)
186 real(real64), contiguous, intent(out) :: force(:, :)
187
188 interface
189 subroutine f90_vdw_calculate(natoms, dd, sr, zatom, coordinates, vol_ratio, &
190 energy, force, derivative_coeff)
191 use, intrinsic :: iso_fortran_env
192 implicit none
193 integer, intent(in) :: natoms
194 real(real64), intent(in) :: dd
195 real(real64), intent(in) :: sr
196 integer, intent(in) :: zatom
197 real(real64), intent(in) :: coordinates
198 real(real64), intent(in) :: vol_ratio
199 real(real64), intent(out) :: energy
200 real(real64), intent(out) :: force
201 real(real64), intent(out) :: derivative_coeff
202 end subroutine f90_vdw_calculate
203 end interface
204
205 type(lattice_iterator_t) :: latt_iter
206 integer :: iatom, jatom, ispecies, jspecies, jcopy, ip
207 real(real64) :: rr, rr6, dffdr0, ee, ff, dee, dffdrab, dffdvra, deabdvra
208 real(real64), allocatable :: coordinates(:,:), vol_ratio(:), dvadens(:), &
209 dr0dvra(:), r0ab(:,:)
210 type(hirshfeld_t) :: hirshfeld
211 integer, allocatable :: zatom(:)
212 real(real64) :: x_j(space%dim)
213
214 push_sub(vdw_ts_calculate)
215
216 safe_allocate(vol_ratio(1:natoms))
217 safe_allocate(dvadens(1:mesh%np))
218 safe_allocate(dr0dvra(1:natoms))
219
220 energy=m_zero
221 force(1:space%dim, 1:natoms) = m_zero
222 this%derivative_coeff(1:natoms) = m_zero
223 call hirshfeld_init(hirshfeld, mesh, namespace, space, latt, atom, natoms, pos, nspin)
224
225 do iatom = 1, natoms
226 call hirshfeld_volume_ratio(hirshfeld, iatom, density, vol_ratio(iatom))
227 end do
228
229 do iatom = 1, natoms
230 ispecies = atom(iatom)%species%get_index()
231 dr0dvra(iatom) = this%r0free(ispecies)/(m_three*(vol_ratio(iatom)**(m_two/m_three)))
232 do jatom = 1, natoms
233 jspecies = atom(jatom)%species%get_index()
234 this%c6ab(iatom,jatom) = vol_ratio(iatom)*vol_ratio(jatom)*this%c6abfree(ispecies,jspecies) !this operation is done again inside the .c part for the non periodic case
235 end do
236 end do
237
238 if (space%is_periodic()) then ! periodic case
239 safe_allocate(r0ab(1:natoms,1:natoms))
240
241 !Precomputing some quantities
242 do iatom = 1, natoms
243 ispecies = atom(iatom)%species%get_index()
244 do jatom = 1, natoms
245 jspecies = atom(jatom)%species%get_index()
246
247 r0ab(iatom,jatom) = (vol_ratio(iatom)**(m_one/m_three))*this%r0free(ispecies) &
248 + (vol_ratio(jatom)**(m_one/m_three))*this%r0free(jspecies)
249 end do
250 end do
251
252 latt_iter = lattice_iterator_t(latt, this%cutoff)
253 do jatom = 1, natoms
254 jspecies = atom(jatom)%species%get_index()
255
256 do jcopy = 1, latt_iter%n_cells ! one of the periodic copy is the initial atom
257 x_j = pos(:, jatom) + latt_iter%get(jcopy)
258
259 do iatom = 1, natoms
260 ispecies = atom(iatom)%species%get_index()
261 rr = norm2(x_j - pos(:, iatom))
262 rr6 = rr**6
263
264 if (rr < r_min_atom_dist) cycle !To avoid self interaction
265
266 ee = exp(-this%damping * ((rr / (this%sr*r0ab(iatom, jatom))) - m_one))
267 ff = m_one/(m_one + ee)
268 dee = ee*ff**2
269
270 !Calculate the derivative of the damping function with respect to the distance between atoms A and B.
271 dffdrab = (this%damping/(this%sr*r0ab(iatom, jatom)))*dee
272 !Calculate the derivative of the damping function with respect to the distance between the van der Waals radius.
273 dffdr0 = -this%damping*rr/(this%sr*r0ab(iatom, jatom)**2)*dee
274
275 energy = energy - m_half*ff*this%c6ab(iatom, jatom)/rr6
276
277 ! Derivative of the damping function with respecto to the volume ratio of atom A.
278 dffdvra = dffdr0*dr0dvra(iatom)
279
280 ! Calculation of the pair-wise partial energy derivative with respect to the volume ratio of atom A.
281 deabdvra = (dffdvra*this%c6ab(iatom, jatom) + ff*vol_ratio(jatom)*this%c6abfree(ispecies, jspecies))/rr6
282
283 this%derivative_coeff(iatom) = this%derivative_coeff(iatom) + deabdvra
284
285 end do
286 end do
287 end do
288
289 safe_deallocate_a(r0ab)
290 else ! Non periodic case
291 safe_allocate(coordinates(1:space%dim, 1:natoms))
292 safe_allocate(zatom(1:natoms))
293
294 do iatom = 1, natoms
295 coordinates(:, iatom) = pos(:, iatom)
296 zatom(iatom) = int(atom(iatom)%species%get_z())
297
298 end do
299
300 call f90_vdw_calculate(natoms, this%damping, this%sr, zatom(1), coordinates(1, 1), &
301 vol_ratio(1), energy, force(1, 1), this%derivative_coeff(1))
302
303
304 safe_deallocate_a(coordinates)
305 safe_deallocate_a(zatom)
306 end if
307
308 ! Calculate the potential
309 potential = m_zero
310 do iatom = 1, natoms
311 call hirshfeld_density_derivative(hirshfeld, iatom, dvadens)
312 call lalg_axpy(mesh%np, -this%derivative_coeff(iatom), dvadens, potential)
313 end do
314
315 if (debug%info) then
316 call dio_function_output(1_8, "./", "vvdw", namespace, space, mesh, potential, unit_one, ip)
317 end if
318
319 call hirshfeld_end(hirshfeld)
320
321 safe_deallocate_a(vol_ratio)
322 safe_deallocate_a(dvadens)
323 safe_deallocate_a(dr0dvra)
324
325 pop_sub(vdw_ts_calculate)
326 end subroutine vdw_ts_calculate
327
328
329
330 !------------------------------------------
331 subroutine vdw_ts_force_calculate(this, force_vdw, ions, mesh, nspin, density)
332 type(vdw_ts_t), intent(in) :: this
333 type(ions_t), intent(in) :: ions
334 real(real64), intent(inout) :: force_vdw(1:ions%space%dim, 1:ions%natoms)
335 class(mesh_t), intent(in) :: mesh
336 integer, intent(in) :: nspin
337 real(real64), intent(in) :: density(:, :)
338
339 type(hirshfeld_t) :: hirshfeld
340 type(lattice_iterator_t) :: latt_iter
341
342 integer :: iatom, jatom, ispecies, jspecies, jcopy
343 real(real64) :: rr, rr6, dffdr0, ee, ff, dee, dffdvra, deabdvra, deabdrab, x_j(ions%space%dim)
344 real(real64), allocatable :: vol_ratio(:), dvadrr(:), dr0dvra(:), r0ab(:,:), derivative_coeff(:), c6ab(:,:)
345
346 push_sub(vdw_ts_force_calculate)
347
348 call profiling_in("FORCE_VDW_TS")
349
350 safe_allocate(vol_ratio(1:ions%natoms))
351 safe_allocate(dvadrr(1:3))
352 safe_allocate(dr0dvra(1:ions%natoms))
353 safe_allocate(r0ab(1:ions%natoms,1:ions%natoms))
354 safe_allocate(derivative_coeff(1:ions%natoms))
355 safe_allocate(c6ab(1:ions%natoms,1:ions%natoms))
356
357
358 force_vdw(1:ions%space%dim, 1:ions%natoms) = m_zero
359 derivative_coeff(1:ions%natoms) = m_zero
360 c6ab(1:ions%natoms,1:ions%natoms) = m_zero
361 r0ab(1:ions%natoms,1:ions%natoms) = m_zero
362 dr0dvra(1:ions%natoms) = m_zero
363 dvadrr(1:ions%space%dim) = m_zero
364 vol_ratio(1:ions%natoms) = m_zero
365
366
367 call hirshfeld_init(hirshfeld, mesh, ions%namespace, ions%space, ions%latt, ions%atom, ions%natoms, ions%pos, nspin)
368
369
370 do iatom = 1, ions%natoms
371 call hirshfeld_volume_ratio(hirshfeld, iatom, density, vol_ratio(iatom))
372 end do
373
374 do iatom = 1, ions%natoms
375 ispecies = ions%atom(iatom)%species%get_index()
376 dr0dvra(iatom) = this%r0free(ispecies)/(m_three*(vol_ratio(iatom)**(m_two/m_three)))
377 do jatom = 1, ions%natoms
378 jspecies = ions%atom(jatom)%species%get_index()
379 c6ab(iatom, jatom) = vol_ratio(iatom)*vol_ratio(jatom)*this%c6abfree(ispecies, jspecies)
380 end do
381 end do
382
383 !Precomputing some quantities
384 do iatom = 1, ions%natoms
385 ispecies = ions%atom(iatom)%species%get_index()
386 do jatom = iatom, ions%natoms
387 jspecies = ions%atom(jatom)%species%get_index()
388
389 r0ab(iatom, jatom) = (vol_ratio(iatom)**(m_one/m_three))*this%r0free(ispecies) &
390 + (vol_ratio(jatom)**(m_one/m_three))*this%r0free(jspecies)
391 if (iatom /= jatom) r0ab(jatom, iatom) = r0ab(iatom, jatom)
392 end do
393 end do
394
395 latt_iter = lattice_iterator_t(ions%latt, this%cutoff)
396 do jatom = 1, ions%natoms
397 jspecies = ions%atom(jatom)%species%get_index()
398
399 do jcopy = 1, latt_iter%n_cells ! one of the periodic copy is the initial atom
400 x_j = ions%pos(:, jatom) + latt_iter%get(jcopy)
401 do iatom = 1, ions%natoms
402 ispecies = ions%atom(iatom)%species%get_index()
403 rr = norm2(x_j - ions%pos(:, iatom))
404 rr6 = rr**6
405
406 if (rr < r_min_atom_dist) cycle !To avoid self interaction
407
408 ee = exp(-this%damping*(rr/(this%sr*r0ab(iatom, jatom)) - m_one))
409 ff = m_one/(m_one + ee)
410 dee = ee*ff**2
411 !Calculate the derivative of the damping function with respect to the van der Waals radius.
412 dffdr0 = -this%damping*rr/( this%sr*r0ab(iatom, jatom)**2)*dee
413 ! Calculation of the pair-wise partial energy derivative with respect to the distance between atoms A and B.
414 deabdrab = c6ab(iatom,jatom)*(this%damping/(this%sr*r0ab(iatom, jatom))*dee - 6.0_real64*ff/rr)/rr6
415 ! Derivative of the damping function with respecto to the volume ratio of atom A.
416 dffdvra = dffdr0*dr0dvra(iatom)
417 ! Calculation of the pair-wise partial energy derivative with respect to the volume ratio of atom A.
418 deabdvra = (dffdvra*c6ab(iatom, jatom) + ff*vol_ratio(jatom)*this%c6abfree(ispecies, jspecies))/rr6
419 !Summing for using later
420 derivative_coeff(iatom) = derivative_coeff(iatom) + deabdvra
421
422 ! Calculation of the pair-wise partial energy derivative with respect to the distance between atoms A and B.
423 force_vdw(:, iatom) = force_vdw(:, iatom) + m_half*deabdrab*(ions%pos(:, iatom) - x_j)/rr
424 end do
425 end do
426 end do
427
428 do iatom = 1, ions%natoms
429 do jatom = 1, ions%natoms
430 call hirshfeld_position_derivative(hirshfeld, ions%space, iatom, jatom, density, dvadrr) !dvadrr_ij = \frac{\delta V_i}{\delta \vec{x_j}}
431 force_vdw(:, jatom)= force_vdw(:, jatom) + derivative_coeff(iatom)*dvadrr ! ions%atom(jatom)%f_vdw = sum_i coeff_i * dvadrr_ij
432 end do
433 end do
434
435 call hirshfeld_end(hirshfeld)
436
437 safe_deallocate_a(vol_ratio)
438 safe_deallocate_a(dvadrr)
439 safe_deallocate_a(dr0dvra)
440 safe_deallocate_a(r0ab)
441 safe_deallocate_a(derivative_coeff)
442 safe_deallocate_a(c6ab)
443
444 call profiling_out("FORCE_VDW_TS")
445
447 end subroutine vdw_ts_force_calculate
448
449 !------------------------------------------
450 subroutine vdw_ts_write_c6ab(this, ions, dir, fname, namespace)
451 type(vdw_ts_t), intent(in) :: this
452 type(ions_t), intent(in) :: ions
453 character(len=*), intent(in) :: dir, fname
454 type(namespace_t), intent(in) :: namespace
455
456 integer :: iunit, iatom, jatom
457
458 push_sub(vdw_ts_write_c6ab)
459
460 if (mpi_grp_is_root(mpi_world)) then
461 call io_mkdir(dir, namespace)
462 iunit = io_open(trim(dir) // "/" // trim(fname), namespace, action='write')
463 write(iunit, '(a)') ' # Atom1 Atom2 C6_{12}^{eff}'
464
465
466 do iatom = 1, ions%natoms
467 do jatom = 1, ions%natoms
468 write(iunit, '(3x, i5, i5, e15.6)') iatom, jatom, this%c6ab(iatom, jatom)
469 end do
470 end do
471 call io_close(iunit)
472 end if
473 pop_sub(vdw_ts_write_c6ab)
474
475 end subroutine vdw_ts_write_c6ab
476
477
478 !------------------------------------------
479
480 subroutine get_vdw_param(namespace, atom, c6, alpha, r0)
481 type(namespace_t), intent(in) :: namespace
482 character(len=*), intent(in) :: atom
483 real(real64), intent(out) :: c6
484 real(real64), intent(out) :: alpha
485 real(real64), intent(out) :: r0
486
487 push_sub(get_vdw_param)
488
489 select case (trim(atom))
490
491 case ('H')
492 alpha = 4.500000_real64
493 c6 = 6.500000_real64
494 r0 = 3.100000_real64
495
496 case ('He')
497 alpha = 1.380000_real64
498 c6 = 1.460000_real64
499 r0 = 2.650000_real64
500
501 case ('Li')
502 alpha = 164.200000_real64
503 c6 = 1387.000000_real64
504 r0 = 4.160000_real64
505
506 case ('Be')
507 alpha = 38.000000_real64
508 c6 = 214.000000_real64
509 r0 = 4.170000_real64
510
511 case ('B')
512 alpha = 21.000000_real64
513 c6 = 99.500000_real64
514 r0 = 3.890000_real64
515
516 case ('C')
517 alpha = 12.000000_real64
518 c6 = 46.600000_real64
519 r0 = 3.590000_real64
520
521 case ('N')
522 alpha = 7.400000_real64
523 c6 = 24.200000_real64
524 r0 = 3.340000_real64
525
526 case ('O')
527 alpha = 5.400000_real64
528 c6 = 15.600000_real64
529 r0 = 3.190000_real64
530
531 case ('F')
532 alpha = 3.800000_real64
533 c6 = 9.520000_real64
534 r0 = 3.040000_real64
535
536 case ('Ne')
537 alpha = 2.670000_real64
538 c6 = 6.380000_real64
539 r0 = 2.910000_real64
540
541 case ('Na')
542 alpha = 162.700000_real64
543 c6 = 1556.000000_real64
544 r0 = 3.730000_real64
545
546 case ('Mg')
547 alpha = 71.000000_real64
548 c6 = 627.000000_real64
549 r0 = 4.270000_real64
550
551 case ('Al')
552 alpha = 60.000000_real64
553 c6 = 528.000000_real64
554 r0 = 4.330000_real64
555
556 case ('Si')
557 alpha = 37.000000_real64
558 c6 = 305.000000_real64
559 r0 = 4.200000_real64
560
561 case ('P')
562 alpha = 25.000000_real64
563 c6 = 185.000000_real64
564 r0 = 4.010000_real64
565
566 case ('S')
567 alpha = 19.600000_real64
568 c6 = 134.000000_real64
569 r0 = 3.860000_real64
570
571 case ('Cl')
572 alpha = 15.000000_real64
573 c6 = 94.600000_real64
574 r0 = 3.710000_real64
575
576 case ('Ar')
577 alpha = 11.100000_real64
578 c6 = 64.300000_real64
579 r0 = 3.550000_real64
580
581 case ('K')
582 alpha = 292.900000_real64
583 c6 = 3897.000000_real64
584 r0 = 3.710000_real64
585
586 case ('Ca')
587 alpha = 160.000000_real64
588 c6 = 2221.000000_real64
589 r0 = 4.650000_real64
590
591 case ('Sc')
592 alpha = 120.000000_real64
593 c6 = 1383.000000_real64
594 r0 = 4.590000_real64
595
596 case ('Ti')
597 alpha = 98.000000_real64
598 c6 = 1044.000000_real64
599 r0 = 4.510000_real64
600
601 case ('V')
602 alpha = 84.000000_real64
603 c6 = 832.000000_real64
604 r0 = 4.440000_real64
605
606 case ('Cr')
607 alpha = 78.000000_real64
608 c6 = 602.000000_real64
609 r0 = 3.990000_real64
610
611 case ('Mn')
612 alpha = 63.000000_real64
613 c6 = 552.000000_real64
614 r0 = 3.970000_real64
615
616 case ('Fe')
617 alpha = 56.000000_real64
618 c6 = 482.000000_real64
619 r0 = 4.230000_real64
620
621 case ('Co')
622 alpha = 50.000000_real64
623 c6 = 408.000000_real64
624 r0 = 4.180000_real64
625
626 case ('Ni')
627 alpha = 48.000000_real64
628 c6 = 373.000000_real64
629 r0 = 3.820000_real64
630
631 case ('Cu')
632 alpha = 42.000000_real64
633 c6 = 253.000000_real64
634 r0 = 3.760000_real64
635
636 case ('Zn')
637 alpha = 40.000000_real64
638 c6 = 284.000000_real64
639 r0 = 4.020000_real64
640
641 case ('Ga')
642 alpha = 60.000000_real64
643 c6 = 498.000000_real64
644 r0 = 4.190000_real64
645
646 case ('Ge')
647 alpha = 41.000000_real64
648 c6 = 354.000000_real64
649 r0 = 4.200000_real64
650
651 case ('As')
652 alpha = 29.000000_real64
653 c6 = 246.000000_real64
654 r0 = 4.110000_real64
655
656 case ('Se')
657 alpha = 25.000000_real64
658 c6 = 210.000000_real64
659 r0 = 4.040000_real64
660
661 case ('Br')
662 alpha = 20.000000_real64
663 c6 = 162.000000_real64
664 r0 = 3.930000_real64
665
666 case ('Kr')
667 alpha = 16.800000_real64
668 c6 = 129.600000_real64
669 r0 = 3.820000_real64
670
671 case ('Rb')
672 alpha = 319.200000_real64
673 c6 = 4691.000000_real64
674 r0 = 3.720000_real64
675
676 case ('Sr')
677 alpha = 199.000000_real64
678 c6 = 3170.000000_real64
679 r0 = 4.540000_real64
680
681 case ('Y')
682 alpha = 126.7370_real64
683 c6 = 1968.580_real64
684 r0 = 4.81510_real64
685
686 case ('Zr')
687 alpha = 119.97_real64
688 c6 = 1677.91_real64
689 r0 = 4.53_real64
690
691 case ('Nb')
692 alpha = 101.603_real64
693 c6 = 1263.61_real64
694 r0 = 4.2365_real64
695
696 case ('Mo')
697 alpha = 88.4225785_real64
698 c6 = 1028.73_real64
699 r0 = 4.099_real64
700
701 case ('Tc')
702 alpha = 80.083_real64
703 c6 = 1390.87_real64
704 r0 = 4.076_real64
705
706 case ('Ru')
707 alpha = 65.8950_real64
708 c6 = 609.754_real64
709 r0 = 3.99530_real64
710
711
712 case ('Rh')
713 alpha = 56.1_real64
714 c6 = 469.0_real64
715 r0 = 3.95_real64
716
717 case ('Pd')
718 alpha = 23.680000_real64
719 c6 = 157.500000_real64
720 r0 = 3.66000_real64
721
722 case ('Ag')
723 alpha = 50.600000_real64
724 c6 = 339.000000_real64
725 r0 = 3.820000_real64
726
727 case ('Cd')
728 alpha = 39.7_real64
729 c6 = 452.0_real64
730 r0 = 3.99_real64
731 case ('Sn')
732 alpha = 55.9500_real64
733 c6 = 587.41700_real64
734 r0 = 4.303000_real64
735 !
736 case ('Sb')
737 alpha = 43.671970_real64
738 c6 = 459.322_real64
739 r0 = 4.2760_real64
740
741
742 case ('Te')
743 alpha = 37.65_real64
744 c6 = 396.0_real64
745 r0 = 4.22_real64
746
747 case ('I')
748 alpha = 35.000000_real64
749 c6 = 385.000000_real64
750 r0 = 4.170000_real64
751
752 case ('Xe')
753 alpha = 27.300000_real64
754 c6 = 285.900000_real64
755 r0 = 4.080000_real64
756
757 case ('Cs')
758 alpha = 427.12_real64
759 c6 = 6582.08_real64
760 r0 = 3.78_real64
761
762 case ('Ba')
763 alpha = 275.0_real64
764 c6 = 5727.0_real64
765 r0 = 4.77_real64
766
767 case ('Hf')
768 alpha = 99.52_real64
769 c6 = 1274.8_real64
770 r0 = 4.21_real64
771
772 case ('Ta')
773 alpha = 82.53_real64
774 c6 = 1019.92_real64
775 r0 = 4.15_real64
776
777 case ('W')
778 alpha = 71.041_real64
779 c6 = 847.93_real64
780 r0 = 4.08_real64
781
782 case ('Re')
783 alpha = 63.04_real64
784 c6 = 710.2_real64
785 r0 = 4.02_real64
786
787 case ('Os')
788 alpha = 55.055_real64
789 c6 = 596.67_real64
790 r0 = 3.84_real64
791
792 case ('Ir')
793 alpha = 42.51_real64
794 c6 = 359.1_real64
795 r0 = 4.00_real64
796
797 case ('Pt')
798 alpha = 39.68_real64
799 c6 = 347.1_real64
800 r0 = 3.92_real64
801
802 case ('Au')
803 alpha = 36.5_real64
804 c6 = 298.0_real64
805 r0 = 3.86_real64
806
807 case ('Hg')
808 alpha = 33.9_real64
809 c6 = 392.0_real64
810 r0 = 3.98_real64
811
812 case ('Tl')
813 alpha = 69.92_real64
814 c6 = 717.44_real64
815 r0 = 3.91_real64
816
817 case ('Pb')
818 alpha = 61.8_real64
819 c6 = 697.0_real64
820 r0 = 4.31_real64
821
822 case ('Bi')
823 alpha = 49.02_real64
824 c6 = 571.0_real64
825 r0 = 4.32_real64
826
827 case ('Po')
828 alpha = 45.013_real64
829 c6 = 530.92_real64
830 r0 = 4.097_real64
831
832 case ('At')
833 alpha = 38.93_real64
834 c6 = 457.53_real64
835 r0 = 4.07_real64
836
837 case ('Rn')
838 alpha = 33.54_real64
839 c6 = 390.63_real64
840 r0 = 4.23_real64
841
842
843 case default
844
845 call messages_write('vdw ts: reference free atom parameters not available for species '//trim(atom))
846 call messages_fatal(namespace=namespace)
847
848 end select
849
850 pop_sub(get_vdw_param)
851 end subroutine get_vdw_param
852
853end module vdw_ts_oct_m
854
855!! Local Variables:
856!! mode: f90
857!! coding: utf-8
858!! End:
constant times a vector plus a vector
Definition: lalg_basic.F90:171
double exp(double __x) __attribute__((__nothrow__
type(debug_t), save, public debug
Definition: debug.F90:156
real(real64), parameter, public m_two
Definition: global.F90:190
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public r_min_atom_dist
Minimal distance between two distinguishable atoms.
Definition: global.F90:183
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_one
Definition: global.F90:189
real(real64), parameter, public m_three
Definition: global.F90:191
subroutine, public hirshfeld_init(this, mesh, namespace, space, latt, atom, natoms, pos, nspin)
Definition: hirshfeld.F90:166
subroutine, public hirshfeld_end(this)
Definition: hirshfeld.F90:258
subroutine, public hirshfeld_density_derivative(this, iatom, ddensity)
Definition: hirshfeld.F90:360
subroutine, public hirshfeld_position_derivative(this, space, iatom, jatom, density, dposition)
Definition: hirshfeld.F90:392
subroutine, public hirshfeld_volume_ratio(this, iatom, density, volume_ratio)
Definition: hirshfeld.F90:322
subroutine, public dio_function_output(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
Top-level IO routine for functions defined on the mesh.
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:418
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:311
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:352
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:415
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
Definition: mpi.F90:430
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
Definition: ps.F90:114
real(real64) function, public ps_density_volume(ps, namespace)
Definition: ps.F90:1654
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
This module defines the unit system, used for input and output.
type(unit_system_t), public units_inp
the units systems for reading and writing
type(unit_t), public unit_one
some special units required for particular quantities
Tkatchenko-Scheffler pairwise method for van der Waals (vdW, dispersion) interactions.
Definition: vdw_ts.F90:119
subroutine, public vdw_ts_init(this, namespace, ions)
Definition: vdw_ts.F90:173
subroutine, public vdw_ts_calculate(this, namespace, space, latt, atom, natoms, pos, mesh, nspin, density, energy, potential, force)
Calculate the potential for the Tatchenko-Scheffler vdW correction as described in Phys....
Definition: vdw_ts.F90:267
subroutine, public vdw_ts_write_c6ab(this, ions, dir, fname, namespace)
Definition: vdw_ts.F90:544
subroutine get_vdw_param(namespace, atom, c6, alpha, r0)
Definition: vdw_ts.F90:574
subroutine, public vdw_ts_force_calculate(this, force_vdw, ions, mesh, nspin, density)
Definition: vdw_ts.F90:425
subroutine, public vdw_ts_end(this)
Definition: vdw_ts.F90:249
The following class implements a lattice iterator. It allows one to loop over all cells that are with...
Describes mesh distribution to nodes.
Definition: mesh.F90:186