34 use,
intrinsic :: iso_fortran_env
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
74 real(real64),
allocatable :: derivative_coeff(:)
80 type(vdw_ts_t),
intent(out) :: this
81 type(namespace_t),
intent(in) :: namespace
82 type(ions_t),
intent(in) :: ions
84 integer :: ispecies, jspecies
85 real(real64) :: num, den
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))
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)
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
156 type(vdw_ts_t),
intent(inout) :: this
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)
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
176 class(
space_t),
intent(in) :: space
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(:, :)
189 subroutine f90_vdw_calculate(natoms, dd, sr, zatom, coordinates, vol_ratio, &
190 energy, force, derivative_coeff)
191 use,
intrinsic :: iso_fortran_env
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
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(:,:)
211 integer,
allocatable :: zatom(:)
212 real(real64) :: x_j(space%dim)
216 safe_allocate(vol_ratio(1:natoms))
217 safe_allocate(dvadens(1:mesh%np))
218 safe_allocate(dr0dvra(1:natoms))
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)
230 ispecies = atom(iatom)%species%get_index()
233 jspecies = atom(jatom)%species%get_index()
234 this%c6ab(iatom,jatom) = vol_ratio(iatom)*vol_ratio(jatom)*this%c6abfree(ispecies,jspecies)
238 if (space%is_periodic())
then
239 safe_allocate(r0ab(1:natoms,1:natoms))
243 ispecies = atom(iatom)%species%get_index()
245 jspecies = atom(jatom)%species%get_index()
247 r0ab(iatom,jatom) = (vol_ratio(iatom)**(
m_one/
m_three))*this%r0free(ispecies) &
254 jspecies = atom(jatom)%species%get_index()
256 do jcopy = 1, latt_iter%n_cells
257 x_j = pos(:, jatom) + latt_iter%get(jcopy)
260 ispecies = atom(iatom)%species%get_index()
261 rr = norm2(x_j - pos(:, iatom))
266 ee =
exp(-this%damping * ((rr / (this%sr*r0ab(iatom, jatom))) -
m_one))
271 dffdrab = (this%damping/(this%sr*r0ab(iatom, jatom)))*dee
273 dffdr0 = -this%damping*rr/(this%sr*r0ab(iatom, jatom)**2)*dee
275 energy = energy -
m_half*ff*this%c6ab(iatom, jatom)/rr6
278 dffdvra = dffdr0*dr0dvra(iatom)
281 deabdvra = (dffdvra*this%c6ab(iatom, jatom) + ff*vol_ratio(jatom)*this%c6abfree(ispecies, jspecies))/rr6
283 this%derivative_coeff(iatom) = this%derivative_coeff(iatom) + deabdvra
289 safe_deallocate_a(r0ab)
291 safe_allocate(coordinates(1:space%dim, 1:natoms))
292 safe_allocate(zatom(1:natoms))
295 coordinates(:, iatom) = pos(:, iatom)
296 zatom(iatom) = int(atom(iatom)%species%get_z())
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))
304 safe_deallocate_a(coordinates)
305 safe_deallocate_a(zatom)
312 call lalg_axpy(mesh%np, -this%derivative_coeff(iatom), dvadens, potential)
321 safe_deallocate_a(vol_ratio)
322 safe_deallocate_a(dvadens)
323 safe_deallocate_a(dr0dvra)
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(:, :)
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(:,:)
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))
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
367 call hirshfeld_init(hirshfeld, mesh, ions%namespace, ions%space, ions%latt, ions%atom, ions%natoms, ions%pos, nspin)
370 do iatom = 1, ions%natoms
374 do iatom = 1, ions%natoms
375 ispecies = ions%atom(iatom)%species%get_index()
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)
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()
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)
396 do jatom = 1, ions%natoms
397 jspecies = ions%atom(jatom)%species%get_index()
399 do jcopy = 1, latt_iter%n_cells
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))
408 ee =
exp(-this%damping*(rr/(this%sr*r0ab(iatom, jatom)) -
m_one))
412 dffdr0 = -this%damping*rr/( this%sr*r0ab(iatom, jatom)**2)*dee
414 deabdrab = c6ab(iatom,jatom)*(this%damping/(this%sr*r0ab(iatom, jatom))*dee - 6.0_real64*ff/rr)/rr6
416 dffdvra = dffdr0*dr0dvra(iatom)
418 deabdvra = (dffdvra*c6ab(iatom, jatom) + ff*vol_ratio(jatom)*this%c6abfree(ispecies, jspecies))/rr6
420 derivative_coeff(iatom) = derivative_coeff(iatom) + deabdvra
423 force_vdw(:, iatom) = force_vdw(:, iatom) +
m_half*deabdrab*(ions%pos(:, iatom) - x_j)/rr
428 do iatom = 1, ions%natoms
429 do jatom = 1, ions%natoms
431 force_vdw(:, jatom)= force_vdw(:, jatom) + derivative_coeff(iatom)*dvadrr
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)
452 type(
ions_t),
intent(in) :: ions
453 character(len=*),
intent(in) :: dir, fname
456 integer :: iunit, iatom, jatom
462 iunit =
io_open(trim(dir) //
"/" // trim(fname), namespace, action=
'write')
463 write(iunit,
'(a)')
' # Atom1 Atom2 C6_{12}^{eff}'
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)
482 character(len=*),
intent(in) :: atom
483 real(real64),
intent(out) :: c6
484 real(real64),
intent(out) :: alpha
485 real(real64),
intent(out) :: r0
489 select case (trim(atom))
492 alpha = 4.500000_real64
497 alpha = 1.380000_real64
502 alpha = 164.200000_real64
503 c6 = 1387.000000_real64
507 alpha = 38.000000_real64
508 c6 = 214.000000_real64
512 alpha = 21.000000_real64
513 c6 = 99.500000_real64
517 alpha = 12.000000_real64
518 c6 = 46.600000_real64
522 alpha = 7.400000_real64
523 c6 = 24.200000_real64
527 alpha = 5.400000_real64
528 c6 = 15.600000_real64
532 alpha = 3.800000_real64
537 alpha = 2.670000_real64
542 alpha = 162.700000_real64
543 c6 = 1556.000000_real64
547 alpha = 71.000000_real64
548 c6 = 627.000000_real64
552 alpha = 60.000000_real64
553 c6 = 528.000000_real64
557 alpha = 37.000000_real64
558 c6 = 305.000000_real64
562 alpha = 25.000000_real64
563 c6 = 185.000000_real64
567 alpha = 19.600000_real64
568 c6 = 134.000000_real64
572 alpha = 15.000000_real64
573 c6 = 94.600000_real64
577 alpha = 11.100000_real64
578 c6 = 64.300000_real64
582 alpha = 292.900000_real64
583 c6 = 3897.000000_real64
587 alpha = 160.000000_real64
588 c6 = 2221.000000_real64
592 alpha = 120.000000_real64
593 c6 = 1383.000000_real64
597 alpha = 98.000000_real64
598 c6 = 1044.000000_real64
602 alpha = 84.000000_real64
603 c6 = 832.000000_real64
607 alpha = 78.000000_real64
608 c6 = 602.000000_real64
612 alpha = 63.000000_real64
613 c6 = 552.000000_real64
617 alpha = 56.000000_real64
618 c6 = 482.000000_real64
622 alpha = 50.000000_real64
623 c6 = 408.000000_real64
627 alpha = 48.000000_real64
628 c6 = 373.000000_real64
632 alpha = 42.000000_real64
633 c6 = 253.000000_real64
637 alpha = 40.000000_real64
638 c6 = 284.000000_real64
642 alpha = 60.000000_real64
643 c6 = 498.000000_real64
647 alpha = 41.000000_real64
648 c6 = 354.000000_real64
652 alpha = 29.000000_real64
653 c6 = 246.000000_real64
657 alpha = 25.000000_real64
658 c6 = 210.000000_real64
662 alpha = 20.000000_real64
663 c6 = 162.000000_real64
667 alpha = 16.800000_real64
668 c6 = 129.600000_real64
672 alpha = 319.200000_real64
673 c6 = 4691.000000_real64
677 alpha = 199.000000_real64
678 c6 = 3170.000000_real64
682 alpha = 126.7370_real64
687 alpha = 119.97_real64
692 alpha = 101.603_real64
697 alpha = 88.4225785_real64
702 alpha = 80.083_real64
707 alpha = 65.8950_real64
718 alpha = 23.680000_real64
719 c6 = 157.500000_real64
723 alpha = 50.600000_real64
724 c6 = 339.000000_real64
732 alpha = 55.9500_real64
733 c6 = 587.41700_real64
737 alpha = 43.671970_real64
748 alpha = 35.000000_real64
749 c6 = 385.000000_real64
753 alpha = 27.300000_real64
754 c6 = 285.900000_real64
758 alpha = 427.12_real64
778 alpha = 71.041_real64
788 alpha = 55.055_real64
828 alpha = 45.013_real64
845 call messages_write(
'vdw ts: reference free atom parameters not available for species '//trim(atom))
constant times a vector plus a vector
double exp(double __x) __attribute__((__nothrow__
type(debug_t), save, public debug
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
real(real64), parameter, public r_min_atom_dist
Minimal distance between two distinguishable atoms.
real(real64), parameter, public m_half
real(real64), parameter, public m_one
real(real64), parameter, public m_three
subroutine, public hirshfeld_init(this, mesh, namespace, space, latt, atom, natoms, pos, nspin)
subroutine, public hirshfeld_end(this)
subroutine, public hirshfeld_density_derivative(this, iatom, ddensity)
subroutine, public hirshfeld_position_derivative(this, space, iatom, jatom, density, dposition)
subroutine, public hirshfeld_volume_ratio(this, iatom, density, volume_ratio)
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.
subroutine, public io_close(iunit, grp)
subroutine, public io_mkdir(fname, namespace, parents)
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
This module defines the meshes, which are used in Octopus.
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
type(mpi_grp_t), public mpi_world
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
real(real64) function, public ps_density_volume(ps, namespace)
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
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.
subroutine, public vdw_ts_init(this, namespace, ions)
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....
subroutine, public vdw_ts_write_c6ab(this, ions, dir, fname, namespace)
subroutine get_vdw_param(namespace, atom, c6, alpha, r0)
subroutine, public vdw_ts_force_calculate(this, force_vdw, ions, mesh, nspin, density)
subroutine, public vdw_ts_end(this)
The following class implements a lattice iterator. It allows one to loop over all cells that are with...
Describes mesh distribution to nodes.