40 use,
intrinsic :: iso_fortran_env
98 class(space_t),
intent(in) :: space
99 type(grid_t),
intent(in) :: gr
100 type(ions_t),
intent(in) :: ions
101 type(epot_t),
intent(in) :: ep
102 type(states_elec_t),
intent(in) :: st
103 type(phase_t),
intent(in) :: phase
104 real(real64),
intent(inout) :: x(:)
105 integer,
intent(in) :: lda_u
125 type(grid_t),
intent(in) :: gr
126 type(namespace_t),
intent(in) :: namespace
127 type(ions_t),
intent(inout) :: ions
128 type(hamiltonian_elec_t),
intent(in) :: hm
129 type(states_elec_t),
intent(in) :: psi
130 type(states_elec_t),
intent(in) :: chi
131 real(real64),
intent(inout) :: ff(:, :)
132 real(real64),
intent(in) :: qq(:, :)
134 integer :: jatom, idim, jdim
135 integer,
target :: j, ist, ik, iatom
136 real(real64) :: rr, w2r, w1r, xx(ions%space%dim), dq, pdot3p, pdot3m, pdot3p2, pdot3m2, dforce1, dforce2
137 complex(real64),
allocatable :: zpsi(:, :), derpsi(:, :, :)
138 real(real64),
allocatable :: forceks1p(:), forceks1m(:), forceks1p2(:), forceks1m2(:), dforceks1(:)
146 do iatom = 1, ions%natoms
147 do jatom = 1, ions%natoms
148 if (jatom == iatom) cycle
149 xx(1:ions%space%dim) = ions%pos(:, jatom) - ions%pos(:, iatom)
150 rr = norm2(xx(1:ions%space%dim))
151 w1r = - ions%charge(iatom) * ions%charge(jatom) / rr**2
152 w2r =
m_two * ions%charge(iatom) * ions%charge(jatom) / rr**3
153 do idim = 1, ions%space%dim
154 do jdim = 1, ions%space%dim
155 ff(iatom, idim) = ff(iatom, idim) + (qq(jatom, jdim) - qq(iatom, jdim)) * w2r * (
m_one/rr**2) * xx(idim) * xx(jdim)
156 ff(iatom, idim) = ff(iatom, idim) - (qq(jatom, jdim) - qq(iatom, jdim)) * w1r * (
m_one/rr**3) * xx(idim) * xx(jdim)
157 if (jdim == idim)
then
158 ff(iatom, idim) = ff(iatom, idim) + (qq(jatom, jdim) - qq(iatom, jdim)) * w1r * (
m_one/rr)
165 safe_allocate(derpsi(1:gr%np_part, 1:ions%space%dim, 1:psi%d%dim))
169 safe_allocate(forceks1p(1:ions%space%dim))
170 safe_allocate(forceks1m(1:ions%space%dim))
171 safe_allocate(forceks1p2(1:ions%space%dim))
172 safe_allocate(forceks1m2(1:ions%space%dim))
173 safe_allocate(dforceks1(1:ions%space%dim))
174 safe_allocate(zpsi(1:gr%np_part, 1:psi%d%dim))
181 do iatom = 1, ions%natoms
182 do j = 1, ions%space%dim
183 call force1(ions%pos(j, iatom) + dq, forceks1p, pdot3p)
184 call force1(ions%pos(j, iatom) - dq, forceks1m, pdot3m)
185 call force1(ions%pos(j, iatom) + dq/
m_two, forceks1p2, pdot3p2)
186 call force1(ions%pos(j, iatom) - dq/
m_two, forceks1m2, pdot3m2)
187 dforceks1 = ((
m_four/
m_three) * (forceks1p2 - forceks1m2) - (
m_one / 6.0_real64) * (forceks1p - forceks1m)) / dq
188 dforce1 = sum(qq(iatom, :) * dforceks1(:))
189 dforce2 = ((
m_four/
m_three) * (pdot3p2 - pdot3m2) - (
m_one / 6.0_real64) * (pdot3p - pdot3m)) / dq
190 ff(iatom, j) = ff(iatom, j) -
m_two * psi%occ(ist, ik) * dforce1 +
m_two * dforce2
196 safe_deallocate_a(zpsi)
197 safe_deallocate_a(forceks1p)
198 safe_deallocate_a(forceks1m)
199 safe_deallocate_a(forceks1p2)
200 safe_deallocate_a(forceks1m2)
201 safe_deallocate_a(dforceks1)
202 safe_deallocate_a(derpsi)
209 subroutine force1(qq, res, pdot3)
210 real(real64),
intent(in) :: qq
211 real(real64),
contiguous,
intent(inout) :: res(:)
212 real(real64),
intent(inout) :: pdot3
216 complex(real64),
allocatable :: viapsi(:, :), zpsi(:, :)
218 qold = ions%pos(j, iatom)
219 ions%pos(j, iatom) = qq
220 safe_allocate(viapsi(1:gr%np_part, 1:psi%d%dim))
221 safe_allocate(zpsi(1:gr%np_part, 1:psi%d%dim))
225 ions%pos(:, iatom), iatom, gr, zpsi, viapsi)
228 do m = 1, ubound(res, 1)
229 res(m) = real(
zmf_dotp(gr, viapsi(:, 1), derpsi(:, m, 1), reduce = .false.), real64)
231 call gr%allreduce(res)
234 pdot3 = real(
m_zi *
zmf_dotp(gr, zpsi(:, 1), viapsi(:, 1)), real64)
235 ions%pos(j, iatom) = qold
237 safe_deallocate_a(viapsi)
245 subroutine forces_calculate(gr, namespace, ions, hm, ext_partners, st, ks, vhxc_old, t, dt)
246 type(
grid_t),
intent(in) :: gr
248 type(
ions_t),
intent(inout) :: ions
252 type(
v_ks_t),
intent(in) :: ks
253 real(real64),
optional,
intent(in) :: vhxc_old(:,:)
254 real(real64),
optional,
intent(in) :: t
255 real(real64),
optional,
intent(in) :: dt
257 integer :: j, iatom, idir
258 real(real64) :: xx(ions%space%dim), time, global_force(ions%space%dim)
259 real(real64),
allocatable :: force(:, :), force_loc(:, :), force_nl(:, :), force_u(:, :)
260 real(real64),
allocatable :: force_nlcc(: ,:)
261 real(real64),
allocatable :: force_scf(:, :)
268 if (
present(t)) time = t
271 do iatom = 1, ions%natoms
272 ions%atom(iatom)%f_ii(1:ions%space%dim) =
m_zero
273 ions%atom(iatom)%f_vdw(1:ions%space%dim) =
m_zero
274 ions%atom(iatom)%f_loc(1:ions%space%dim) =
m_zero
275 ions%atom(iatom)%f_nl(1:ions%space%dim) =
m_zero
276 ions%atom(iatom)%f_u(1:ions%space%dim) =
m_zero
277 ions%atom(iatom)%f_fields(1:ions%space%dim) =
m_zero
278 ions%atom(iatom)%f_nlcc(1:ions%space%dim) =
m_zero
279 ions%atom(iatom)%f_scf(1:ions%space%dim) =
m_zero
280 ions%atom(iatom)%f_photons(1:ions%space%dim) =
m_zero
285 if (ks%vdw%vdw_correction == option__vdwcorrection__vdw_ts)
then
289 do iatom = 1, ions%natoms
290 ions%tot_force(:, iatom) = hm%ep%fii(1:ions%space%dim, iatom) + hm%ep%vdw_forces(1:ions%space%dim, iatom)
291 if (ks%has_photons)
then
292 ions%tot_force(:, iatom) = ions%tot_force(:, iatom) &
293 -
p_proton_charge*ions%charge(iatom)*hm%ep%photon_forces(1:ions%space%dim)
295 ions%atom(iatom)%f_ii(1:ions%space%dim) = hm%ep%fii(1:ions%space%dim, iatom)
296 ions%atom(iatom)%f_vdw(1:ions%space%dim) = hm%ep%vdw_forces(1:ions%space%dim, iatom)
297 ions%atom(iatom)%f_photons(1:ions%space%dim) = -
p_proton_charge*ions%charge(iatom)*hm%ep%photon_forces(1:ions%space%dim)
301 global_force = ions%global_force(time)
304 do iatom = 1, ions%natoms
305 ions%tot_force(:, iatom) = ions%tot_force(:, iatom) + global_force
306 ions%atom(iatom)%f_ii(1:ions%space%dim) = ions%atom(iatom)%f_ii(1:ions%space%dim) + global_force
310 safe_allocate(force(1:ions%space%dim, 1:ions%natoms))
311 safe_allocate(force_loc(1:ions%space%dim, 1:ions%natoms))
312 safe_allocate(force_nl(1:ions%space%dim, 1:ions%natoms))
313 safe_allocate(force_u(1:ions%space%dim, 1:ions%natoms))
314 safe_allocate(force_scf(1:ions%space%dim, 1:ions%natoms))
315 safe_allocate(force_nlcc(1:ions%space%dim, 1:ions%natoms))
324 .and. ks%theory_level /=
hartree .and. ks%theory_level /=
rdmft)
then
325 call forces_from_nlcc(gr, ions, st%d%spin_channels, hm%ks_pot%vxc, force_nlcc)
329 if (
present(vhxc_old))
then
330 call forces_from_scf(gr, ions, st%d%spin_channels, hm%ks_pot%vhxc, vhxc_old, force_scf)
335 if (ions%force_total_enforce)
then
344 do iatom = 1, ions%natoms
345 do idir = 1, ions%space%dim
346 ions%tot_force(idir, iatom) = ions%tot_force(idir, iatom) + force(idir, iatom) &
347 + force_scf(idir, iatom) + force_nlcc(idir, iatom)
348 ions%atom(iatom)%f_loc(idir) = force_loc(idir, iatom)
349 ions%atom(iatom)%f_nl(idir) = force_nl(idir, iatom)
350 ions%atom(iatom)%f_u(idir) = force_u(idir, iatom)
351 ions%atom(iatom)%f_nlcc(idir) = force_nlcc(idir, iatom)
352 ions%atom(iatom)%f_scf(idir) = force_scf(idir, iatom)
356 safe_deallocate_a(force)
357 safe_deallocate_a(force_loc)
358 safe_deallocate_a(force_nl)
359 safe_deallocate_a(force_u)
360 safe_deallocate_a(force_nlcc)
361 safe_deallocate_a(force_scf)
365 if (
present(t) .and.
associated(lasers))
then
366 do j = 1, lasers%no_lasers
371 do iatom = 1, ions%natoms
373 ions%tot_force(:, iatom) = ions%tot_force(:, iatom) + ions%charge(iatom)*xx(:)
374 ions%atom(iatom)%f_fields(1:ions%space%dim) = ions%atom(iatom)%f_fields(1:ions%space%dim) &
375 + ions%charge(iatom)*xx(:)
387 do iatom = 1, ions%natoms
389 ions%tot_force(:, iatom) = ions%tot_force(:, iatom) + ions%charge(iatom)*xx(:)
390 ions%atom(iatom)%f_fields(1:ions%space%dim) = ions%atom(iatom)%f_fields(1:ions%space%dim) &
391 + ions%charge(iatom)*xx(:)
395 write(
message(1),
'(a)')
'The forces are currently not supported for nondipole '
396 write(
message(2),
'(a)')
'strong-field approximation Hamiltonian approach.'
401 write(
message(1),
'(a)')
'The forces are currently not properly calculated if time-dependent'
402 write(
message(2),
'(a)')
'magnetic fields are present.'
408 if (
allocated(hm%ep%e_field))
then
409 do iatom = 1, ions%natoms
411 ions%tot_force(:, iatom) = ions%tot_force(:, iatom) + ions%charge(iatom)*hm%ep%e_field(1:ions%space%dim)
412 ions%atom(iatom)%f_fields(1:ions%space%dim) = ions%atom(iatom)%f_fields(1:ions%space%dim) &
413 + ions%charge(iatom)*hm%ep%e_field(1:ions%space%dim)
417 if (
allocated(hm%ep%b_field) .or.
allocated(hm%ep%a_static))
then
418 write(
message(1),
'(a)')
'The forces are currently not properly calculated if static'
419 write(
message(2),
'(a)')
'magnetic fields or static vector potentials are present.'
424 write(
message(1),
'(a)')
'The forces are currently not properly calculated if gauge-field'
425 write(
message(2),
'(a)')
'is applied.'
434 if (hm%kpoints%use_symmetries .or. st%symmetrize_density)
then
439 if (ions%space%is_periodic())
then
451 type(
ions_t),
intent(in) :: ions
452 real(real64),
intent(inout) :: force(:, :)
454 real(real64),
allocatable :: total_force(:)
459 safe_allocate(total_force(1:ions%space%dim))
461 total_force(1:ions%space%dim) = sum(force, dim=2) / ions%natoms
463 do iatom = 1, ions%natoms
464 force(1:ions%space%dim, iatom) = force(1:ions%space%dim, iatom) - total_force(1:ions%space%dim)
467 safe_deallocate_a(total_force)
474 type(
ions_t),
intent(in) :: ions
475 real(real64),
intent(inout) :: total_torque(:)
477 real(real64) :: center_of_mass(ions%space%dim), rr(3), ff(3)
482 center_of_mass = ions%center_of_mass()
487 do iatom = 1, ions%natoms
488 rr(1:ions%space%dim) = ions%pos(1:ions%space%dim, iatom) - center_of_mass
489 ff(1:ions%space%dim) = ions%tot_force(1:ions%space%dim, iatom)
500 integer,
intent(in) :: iunit
501 type(
ions_t),
intent(in) :: ions
502 character(len=*),
intent(in) :: dir
505 integer :: iatom, idir, ii, iunit2
506 real(real64) :: torque(1:3)
513 write(iunit,
'(a,10x,99(14x,a))')
' Ion', (
index2axis(idir), idir = 1, ions%space%dim)
514 do iatom = 1, ions%natoms
515 write(iunit,
'(i4,a10,10es17.8)') iatom, trim(ions%atom(iatom)%species%get_label()), &
518 write(iunit,
'(1x,100a1)') (
"-", ii = 1, 13 + ions%space%dim * 15)
519 write(iunit,
'(a14, 10es17.8)')
" Max abs force", &
521 write(iunit,
'(a14, 10es17.8)')
" Total force", &
524 if ((ions%space%dim == 2 .or. ions%space%dim == 3) .and. .not. ions%space%is_periodic())
then
526 write(iunit,
'(a14, 10es17.8)')
' Total torque', &
532 iunit2 =
io_open(trim(dir)//
'/forces', namespace, action=
'write', position=
'asis')
533 write(iunit2,
'(a)') &
534 ' # Total force (x,y,z) Ion-Ion (x,y,z) VdW (x,y,z) Local (x,y,z) NL (x,y,z)' // &
535 ' Fields (x,y,z) Hubbard(x,y,z) SCF(x,y,z) NLCC(x,y,z) Phot (x,y,z)'
536 do iatom = 1, ions%natoms
537 write(iunit2,
'(i4,a10,30es17.8)') iatom, trim(ions%atom(iatom)%species%get_label()), &
561 type(
grid_t),
intent(in) :: mesh
562 type(
ions_t),
intent(inout) :: ions
563 integer,
intent(in) :: spin_channels
564 real(real64),
intent(in) :: vxc(:,:)
565 real(real64),
contiguous,
intent(out) :: force_nlcc(:, :)
567 integer :: is, iatom, idir
568 real(real64),
allocatable :: pot(:), grad_vxc(:, :)
574 safe_allocate(pot(1:mesh%np_part))
575 safe_allocate(grad_vxc(1:mesh%np, 1:ions%space%dim))
577 pot(1:mesh%np) = vxc(1:mesh%np,1)
578 do is = 2, spin_channels
579 pot(1:mesh%np) = pot(1:mesh%np) + vxc(1:mesh%np, is)
586 do iatom = ions%atoms_dist%start, ions%atoms_dist%end
587 call species_get_nlcc(ions%atom(iatom)%species, ions%space, ions%latt, ions%pos(:, iatom), mesh, pot)
589 do idir = 1, ions%space%dim
590 force_nlcc(idir, iatom) = force_nlcc(idir, iatom) &
591 -
dmf_dotp(mesh, pot(:), grad_vxc(:, idir), reduce = .false.)
595 safe_deallocate_a(pot)
596 safe_deallocate_a(grad_vxc)
598 if (ions%atoms_dist%parallel)
call dforces_gather(ions, force_nlcc)
601 call mesh%allreduce(force_nlcc)
614 subroutine forces_from_scf(mesh, ions, spin_channels, vhxc, vhxc_old, force_scf)
615 class(
mesh_t),
intent(in) :: mesh
616 type(
ions_t),
intent(inout) :: ions
617 integer,
intent(in) :: spin_channels
618 real(real64),
intent(in) :: vhxc(:,:)
619 real(real64),
intent(in) :: vhxc_old(:,:)
620 real(real64),
contiguous,
intent(out) :: force_scf(:, :)
622 integer :: is, iatom, idir
623 real(real64),
allocatable :: dvhxc(:,:), drho(:,:,:)
629 safe_allocate(dvhxc(1:mesh%np, 1:spin_channels))
630 safe_allocate(drho(1:mesh%np, 1:spin_channels, 1:ions%space%dim))
633 do is = 1, spin_channels
634 dvhxc(1:mesh%np, is) = vhxc(1:mesh%np, is) - vhxc_old(1:mesh%np, is)
639 do iatom = ions%atoms_dist%start, ions%atoms_dist%end
640 select type(spec=>ions%atom(iatom)%species)
646 ions%pos(:, iatom), mesh, spin_channels, drho)
648 do idir = 1, ions%space%dim
649 do is = 1, spin_channels
650 force_scf(idir, iatom) = force_scf(idir, iatom) &
651 -
dmf_dotp(mesh, drho(:,is,idir), dvhxc(:,is), reduce = .false.)
658 safe_deallocate_a(dvhxc)
659 safe_deallocate_a(drho)
661 if (ions%atoms_dist%parallel)
call dforces_gather(ions, force_scf)
664 call mesh%allreduce(force_scf)
674 class(
mesh_t),
intent(in) :: mesh
675 class(
space_t),
intent(in) :: space
676 real(real64),
intent(in) :: vpsl(:)
677 real(real64),
intent(in) :: gdensity(:, :)
678 real(real64),
intent(inout) :: force(:)
681 real(real64) :: force_tmp(1:space%dim)
685 do idir = 1, space%dim
686 force_tmp(idir) =
dmf_dotp(mesh, vpsl(1:mesh%np), gdensity(:, idir), reduce = .false.)
689 call mesh%allreduce(force_tmp)
691 force(1:space%dim) = force(1:space%dim) + force_tmp(1:space%dim)
699 type(
ions_t),
intent(in) :: ions
700 real(real64),
intent(inout) :: force(:,:)
702 integer :: iatom, iop, iatom_symm
703 real(real64) :: symmetrized_force(ions%space%dim), force_tmp(ions%space%dim)
704 real(real64),
allocatable :: force_sym(:,:)
708 safe_allocate(force_sym(1:ions%space%dim,1:ions%natoms))
709 do iatom = 1, ions%natoms
710 symmetrized_force =
m_zero
711 do iop = 1, ions%symm%nops
712 iatom_symm = ions%inv_map_symm_atoms(iatom, iop)
714 symmetrized_force = symmetrized_force + force_tmp
717 do iop = 1, ions%symm%nops_nonsymmorphic
718 iatom_symm = ions%inv_map_symm_atoms(iatom, iop + ions%symm%nops)
720 symmetrized_force = symmetrized_force + force_tmp
723 force_sym(:, iatom) = symmetrized_force / (ions%symm%nops + ions%symm%nops_nonsymmorphic)
727 safe_deallocate_a(force_sym)
735#include "forces_inc.F90"
738#include "complex.F90"
739#include "forces_inc.F90"
subroutine force1(qq, res, pdot3)
scales a vector by a constant
This module implements batches of mesh functions.
This module implements common operations on batches of mesh functions.
Module implementing boundary conditions in Octopus.
This module implements a calculator for the density and defines related functions.
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public dderivatives_grad(der, ff, op_ff, ghost_update, set_bc, to_cartesian)
apply the gradient to a mesh function
subroutine, public zderivatives_grad(der, ff, op_ff, ghost_update, set_bc, to_cartesian)
apply the gradient to a mesh function
logical function, public list_has_gauge_field(partners)
type(lasers_t) function, pointer, public list_get_lasers(partners)
subroutine, public zforces_from_potential(gr, namespace, space, ions, hm, st, force, force_loc, force_nl, force_u)
Ref: Kikuji Hirose, Tomoya Ono, Yoshitaka Fujimoto, and Shigeru Tsukamoto, First-principles calculati...
subroutine, public forces_costate_calculate(gr, namespace, ions, hm, psi, chi, ff, qq)
subroutine forces_from_scf(mesh, ions, spin_channels, vhxc, vhxc_old, force_scf)
subroutine symmetrize_force(ions, force)
Given the forces on all atoms, this symmetrizes them using symmorphic and non-symmorphic operations.
subroutine, public zforces_born_charges(gr, namespace, space, ions, ep, st, kpoints, lr, lr2, born_charges, lda_u_level)
lr, lr2 are wfns from electric perturbation; lr is for +omega, lr2 is for -omega. for each atom,...
subroutine forces_compute_total_torque(ions, total_torque)
Computes the total torque acting on the system.
subroutine, public forces_write_info(iunit, ions, dir, namespace)
subroutine forces_set_total_to_zero(ions, force)
subroutine, public dforces_from_potential(gr, namespace, space, ions, hm, st, force, force_loc, force_nl, force_u)
Ref: Kikuji Hirose, Tomoya Ono, Yoshitaka Fujimoto, and Shigeru Tsukamoto, First-principles calculati...
subroutine, public zforces_derivative(gr, namespace, space, ions, ep, st, kpoints, lr, lr2, force_deriv, lda_u_level)
subroutine dtotal_force_from_potential(space, gr, ions, ep, st, phase, x, lda_u_level)
subroutine, public forces_calculate(gr, namespace, ions, hm, ext_partners, st, ks, vhxc_old, t, dt)
subroutine, public total_force_calculate(space, gr, ions, ep, st, phase, x, lda_u)
This computes the total forces on the ions created by the electrons (it excludes the force due to pos...
subroutine, public dforces_derivative(gr, namespace, space, ions, ep, st, kpoints, lr, lr2, force_deriv, lda_u_level)
subroutine ztotal_force_from_potential(space, gr, ions, ep, st, phase, x, lda_u_level)
subroutine, public dforces_born_charges(gr, namespace, space, ions, ep, st, kpoints, lr, lr2, born_charges, lda_u_level)
lr, lr2 are wfns from electric perturbation; lr is for +omega, lr2 is for -omega. for each atom,...
subroutine total_force_from_local_potential(mesh, space, vpsl, gdensity, force)
subroutine dforces_gather(ions, force)
subroutine forces_from_nlcc(mesh, ions, spin_channels, vxc, force_nlcc)
real(real64), parameter, public m_two
real(real64), parameter, public p_proton_charge
real(real64), parameter, public m_zero
real(real64), parameter, public m_four
complex(real64), parameter, public m_z0
complex(real64), parameter, public m_zi
real(real64), parameter, public m_one
real(real64), parameter, public m_three
This module implements the underlying real-space grid.
subroutine, public zhamiltonian_elec_apply_atom(hm, namespace, space, latt, species, pos, ia, mesh, psi, vpsi)
This module defines classes and functions for interaction partners.
subroutine, public io_close(iunit, grp)
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
A module to handle KS potential, without the external potential.
integer, parameter, public rdmft
integer, parameter, public hartree
integer, parameter, public independent_particles
logical function, public lasers_with_nondipole_field(lasers)
Check if a nondipole SFA correction should be computed for the given laser.
integer, parameter, public e_field_electric
integer, parameter, public e_field_vector_potential
subroutine, public laser_electric_field(laser, field, time, dt)
Returns a vector with the electric field, no matter whether the laser is described directly as an ele...
integer, parameter, public e_field_scalar_potential
integer pure elemental function, public laser_kind(laser)
subroutine, public laser_field(laser, field, time)
Retrieves the value of either the electric or the magnetic field. If the laser is given by a scalar p...
integer, parameter, public e_field_magnetic
This modules implements the routines for doing constrain DFT for noncollinear magnetism.
integer, parameter, public constrain_none
This module is intended to contain "only mathematical" functions and procedures.
pure real(real64) function, dimension(1:3), public dcross_product(a, b)
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
subroutine, public messages_not_implemented(feature, namespace)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
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.
pure logical function, public ps_has_density(ps)
subroutine, public species_atom_density_grad(species, namespace, space, latt, pos, mesh, spin_channels, drho)
subroutine, public species_get_nlcc(species, space, latt, pos, mesh, rho_core, accumulate)
pure logical function, public states_are_real(st)
This module handles spin dimensions of the states and the k-point distribution.
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
character(len=20) pure function, public units_abbrev(this)
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
This module is intended to contain simple general-purpose utility functions and procedures.
character pure function, public index2axis(idir)
Tkatchenko-Scheffler pairwise method for van der Waals (vdW, dispersion) interactions.
subroutine, public vdw_ts_force_calculate(this, force_vdw, ions, mesh, nspin, density)
Description of the grid, containing information on derivatives, stencil, and symmetries.
Describes mesh distribution to nodes.
The states_elec_t class contains all electronic wave functions.