43 use,
intrinsic :: iso_fortran_env
93 type(namespace_t),
intent(in) :: namespace
94 type(grid_t),
intent(inout) :: gr
95 type(hamiltonian_elec_t),
intent(inout) :: hm
96 type(states_elec_t),
target,
intent(inout) :: st
97 type(ions_t),
intent(inout) :: ions
98 type(v_ks_t),
intent(in) :: ks
99 type(partner_list_t),
intent(in) :: ext_partners
101 real(real64),
allocatable :: rho_total(:)
102 real(real64) :: stress(3,3)
103 real(real64) :: stress_kin(3,3), stress_Hartree(3,3), stress_xc(3,3), stress_xc_nlcc(3,3)
104 real(real64) :: stress_ps(3,3), stress_ps_nl(3,3), stress_ps_local(3,3), stress_ii(3,3)
105 real(real64) :: stress_hubbard(3,3)
107 real(real64),
allocatable :: vh(:)
108 real(real64),
allocatable :: grad_vh(:,:)
109 real(real64) :: ehartree
110 real(real64),
contiguous,
pointer :: rho(:)
116 write(
message(1),
'(a)')
'The stress tensors for real wavefunctions has not been implemented!'
118 if (hm%kpoints%full%npoints == 1)
then
119 write(
message(2),
'(a)')
'For testing this feature, you can add ForceComplex=yes to the input file'
126 if (ions%space%periodic_dim == 1)
then
130 if (.not. ions%space%is_periodic())
then
131 write(
message(1),
'(a)')
'The stress tensor cannot be computed for isolated systems'
135 if (ks%vdw%vdw_correction /= option__vdwcorrection__none .and. .not. any(ks%vdw%vdw_correction ==
d3_lib_options))
then
136 write(
message(1),
'(a)')
'The stress tensor is currently only implemented with DFT-D3 vdW correction'
140 if (hm%pcm%run_pcm)
then
144 if (
allocated(hm%v_static))
then
148 if (ks%has_photons)
then
152 if (.not. hm%vnl%apply_projector_matrices)
then
166 if ( .not.
in_family(hm%xc%family, [xc_family_lda, xc_family_gga]))
then
167 write(
message(1),
'(a)')
'The stress tensor computation is currently only possible at the Kohn-Sham DFT level'
168 write(
message(2),
'(a)')
'with LDA and GGA functionals or for independent particles.'
172 if (
in_family(hm%xc%family, [xc_family_gga]) .and. st%d%ispin ==
spinors)
then
183 safe_allocate(rho_total(1:gr%np_part))
185 rho_total(ip) = sum(st%rho(ip, 1:st%d%nspin))
194 safe_allocate(vh(1:gr%np_part))
195 safe_allocate(grad_vh(1:gr%np, 1:gr%der%dim))
197 call lalg_copy(gr%np, hm%ks_pot%vhartree, vh)
199 if (hm%d%spin_channels > 1)
then
200 safe_allocate(rho(1:gr%np_part))
208 if (hm%d%spin_channels > 1)
then
209 safe_deallocate_p(rho)
214 ehartree = hm%energy%hartree
221 call stress_from_kinetic(gr, ions%space, hm, st, gr%symm, ions%latt%rcell_volume, stress_kin)
222 stress = stress + stress_kin
229 call stress_from_hartree(gr, ions%space, ions%latt%rcell_volume, vh, grad_vh, ehartree, stress_hartree)
230 stress = stress + stress_hartree
232 call stress_from_xc(hm%energy, ions%latt%rcell_volume, ions%space%periodic_dim, stress_xc)
235 if (
allocated(st%rho_core))
then
236 call stress_from_xc_nlcc(ions%latt%rcell_volume, gr, st, ions, hm%ks_pot%vxc, stress_xc_nlcc)
241 stress_xc = stress_xc + ks%stress_xc_gga / ions%latt%rcell_volume
242 stress = stress + stress_xc + stress_xc_nlcc
246 stress_ps = stress_ps_local
247 stress = stress + stress_ps_local
249 safe_deallocate_a(vh)
250 safe_deallocate_a(grad_vh)
253 stress_ps = stress_ps + stress_ps_nl
254 stress = stress + stress_ps_nl
256 call stress_from_hubbard(namespace, gr, st, hm, ions%space, ions%latt%rcell_volume, stress_hubbard)
257 stress = stress + stress_hubbard
259 call ion_interaction_stress(ions%ion_interaction, ions%space, ions%latt, ions%atom, ions%natoms, ions%pos, stress_ii)
260 stress = stress + stress_ii
266 st%stress_tensors%kinetic = -stress_kin
267 st%stress_tensors%Hartree = -stress_hartree
268 st%stress_tensors%xc = -stress_xc
269 st%stress_tensors%xc_nlcc = -stress_xc_nlcc
270 st%stress_tensors%ps_local = -stress_ps_local
271 st%stress_tensors%ps_nl = -stress_ps_nl
272 st%stress_tensors%hubbard = -stress_hubbard
273 st%stress_tensors%ion_ion = -stress_ii
276 if (ks%vdw%vdw_correction /= option__vdwcorrection__none)
then
277 st%stress_tensors%vdw = hm%ep%vdw_stress
279 st%stress_tensors%vdw =
m_zero
281 stress = stress + st%stress_tensors%vdw
284 if (hm%kpoints%use_symmetries)
then
290 st%stress_tensors%total = stress
294 st%stress_tensors%kinetic_sumrule =
m_zero
296 st%stress_tensors%Hartree_sumrule =
m_zero
297 if(ions%space%periodic_dim == 3)
then
298 st%stress_tensors%kinetic_sumrule = (stress_kin(1,1) + stress_kin(2,2) + stress_kin(3,3))*ions%latt%rcell_volume
299 st%stress_tensors%kinetic_sumrule = st%stress_tensors%kinetic_sumrule -
m_two * hm%energy%kinetic
301 st%stress_tensors%hartree_sumrule = (stress_hartree(1,1) + stress_hartree(2,2) + stress_hartree(3,3))*ions%latt%rcell_volume
302 st%stress_tensors%hartree_sumrule = st%stress_tensors%hartree_sumrule - hm%energy%hartree
305 safe_deallocate_a(rho_total)
327 type(
grid_t),
intent(in) :: gr
328 class(
space_t),
intent(in) :: space
332 real(real64),
intent(in) :: rcell_volume
333 real(real64),
intent(out) :: stress_kin(3, 3)
335 integer :: ik, ist, idir, jdir, ib, minst, maxst
336 complex(real64),
allocatable :: stress_l_block(:)
344 safe_allocate(stress_l_block(1:st%block_size))
346 do ik = st%d%kpt%start, st%d%kpt%end
349 do ib = st%group%block_start, st%group%block_end
359 do idir = 1, space%periodic_dim
360 do jdir = idir, space%periodic_dim
363 do ist = minst, maxst
364 stress_kin(idir,jdir) = stress_kin(idir,jdir) &
365 + st%kweights(ik) * st%occ(ist, ik) &
366 * real(stress_l_block(ist - minst + 1), real64)
371 do idir = 1, space%dim
372 call gpsib(idir)%end()
379 if (st%parallel_in_states .or. st%d%kpt%parallel)
then
388 if (hm%kpoints%use_symmetries)
then
392 stress_kin = stress_kin / rcell_volume
416 type(
grid_t),
intent(in) :: gr
417 class(
space_t),
intent(in) :: space
418 real(real64),
intent(in) :: volume
419 real(real64),
intent(in) :: vh(:)
420 real(real64),
intent(in) :: grad_vh(:,:)
421 real(real64),
intent(in) :: ehartree
422 real(real64),
intent(out) :: stress_Hartree(3, 3)
424 integer :: idir, jdir
429 stress_hartree(:,:) =
m_zero
431 do idir = 1, space%periodic_dim
432 do jdir = idir, space%periodic_dim
433 stress_hartree(idir, jdir) = -
dmf_dotp(gr, grad_vh(:,idir), grad_vh(:, jdir))/
m_four/
m_pi
435 stress_hartree(idir, idir) = stress_hartree(idir, idir) + ehartree
440 stress_hartree = stress_hartree/volume
464 subroutine stress_from_xc(energy, rcell_volume, periodic_dim, stress_xc)
465 type(
energy_t),
intent(in) :: energy
466 real(real64),
intent(in) :: rcell_volume
467 integer,
intent(in) :: periodic_dim
468 real(real64),
intent(out) :: stress_xc(3, 3)
476 do idir = 1, periodic_dim
477 stress_xc(idir, idir) = - energy%exchange - energy%correlation + energy%intnvxc
479 stress_xc(:,:) = stress_xc(:,:) / rcell_volume
496 real(real64),
intent(in) :: rcell_volume
497 type(
grid_t),
intent(in) :: gr
499 type(
ions_t),
intent(in) :: ions
500 real(real64),
intent(in) :: vxc(:,:)
501 real(real64),
intent(out) :: stress_xc_nlcc(3, 3)
503 integer :: idir, jdir, iat, ip
504 real(real64),
allocatable :: nlcc(:), gvxc(:,:), nlcc_x(:,:), vxc_tot(:)
509 assert(
allocated(st%rho_core))
513 safe_allocate(vxc_tot(1:gr%np_part))
514 safe_allocate(nlcc(1:gr%np))
515 safe_allocate(nlcc_x(1:gr%np, 1:gr%der%dim))
518 call lalg_copy(gr%np, vxc(:, 1), vxc_tot)
519 if(st%d%nspin > 1)
call lalg_axpy(gr%np,
m_one, vxc(:, 2), vxc_tot)
525 do iat = ions%atoms_dist%start, ions%atoms_dist%end
527 ions%pos(:,iat), gr, nlcc)
529 do idir = 1, ions%space%periodic_dim
532 nlcc_x(ip, idir) = nlcc_x(ip, idir) + (gr%x(ip, idir) - ions%pos(idir,iat)) * nlcc(ip)
534 stress_xc_nlcc(idir, idir) = stress_xc_nlcc(idir, idir) +
dmf_dotp(gr, nlcc, vxc_tot)
537 safe_deallocate_a(nlcc)
539 if (ions%atoms_dist%parallel)
then
544 safe_allocate(gvxc(1:gr%np, 1:gr%der%dim))
547 do idir = 1, ions%space%periodic_dim
548 do jdir = idir, ions%space%periodic_dim
549 stress_xc_nlcc(idir, jdir) = stress_xc_nlcc(idir, jdir) +
dmf_dotp(gr, gvxc(:,idir), nlcc_x(:,jdir))
552 safe_deallocate_a(vxc_tot)
553 safe_deallocate_a(gvxc)
554 safe_deallocate_a(nlcc_x)
558 stress_xc_nlcc(:,:) = stress_xc_nlcc(:,:) / rcell_volume
584 type(
grid_t),
target,
intent(in) :: gr
587 type(
ions_t),
intent(in) :: ions
588 real(real64),
intent(out) :: stress_ps_nl(3, 3)
590 integer :: ik, ist, idir, jdir
591 integer :: ib, minst, maxst
592 type(
wfs_elec_t) :: psib, rvnl_psib(3), gpsib(3)
593 complex(real64),
allocatable :: stress_tmp(:)
600 safe_allocate(stress_tmp(1:st%block_size))
604 do ik = st%d%kpt%start, st%d%kpt%end
608 do ib = st%group%block_start, st%group%block_end
619 do idir = 1, gr%der%dim
620 call psib%copy_to(rvnl_psib(idir))
623 call hm%vnl%zr_vn_local(gr, st%d, gr%der%boundaries%spiral, psib, rvnl_psib)
625 do idir = 1, ions%space%periodic_dim
626 do jdir = idir, ions%space%periodic_dim
629 do ist = minst, maxst
630 stress_ps_nl(idir, jdir) = stress_ps_nl(idir, jdir) &
631 +
m_two * st%kweights(ik) * st%occ(ist, ik) * real(stress_tmp(ist-minst+1), real64)
637 do idir = 1, gr%der%dim
638 call rvnl_psib(idir)%end()
639 call gpsib(idir)%end()
645 safe_deallocate_a(stress_tmp)
647 if (st%parallel_in_states .or. st%d%kpt%parallel)
then
655 if (hm%kpoints%use_symmetries)
then
660 do idir = 1, ions%space%periodic_dim
661 stress_ps_nl(idir, idir) = stress_ps_nl(idir, idir) + hm%energy%extern_non_local
664 stress_ps_nl = stress_ps_nl/ions%latt%rcell_volume
690 type(
grid_t),
target,
intent(in) :: gr
693 type(
ions_t),
intent(in) :: ions
694 real(real64),
contiguous,
intent(inout) :: rho_total(:)
695 real(real64),
intent(in) :: vh(:)
696 real(real64),
intent(in) :: grad_vh(:,:)
697 real(real64),
intent(out) :: stress_ps_local(3, 3)
700 real(real64) :: stress_SR(3, 3), stress_LR(3, 3)
701 real(real64) :: energy_ps_SR, charge, zi
702 real(real64),
allocatable :: vloc(:), rvloc(:,:), rho_local_lr(:), rho_lr(:)
703 real(real64),
allocatable :: grad_rho(:,:), rho_lr_x(:,:), vlr(:), grad_vlr(:,:)
704 integer :: idir, jdir, iatom
705 type(
ps_t),
pointer :: spec_ps
713 safe_allocate(vloc(1:gr%np))
715 safe_allocate(rvloc(1:gr%np, 1:gr%der%dim))
717 do iatom = 1, ions%natoms
720 safe_deallocate_a(vloc)
722 safe_allocate(grad_rho(1:gr%np,1:gr%der%dim))
725 energy_ps_sr = hm%energy%extern_local
726 do idir = 1, ions%space%periodic_dim
727 do jdir = idir, ions%space%periodic_dim
728 stress_sr(idir, jdir) = stress_sr(idir, jdir) &
729 +
dmf_dotp(gr, rvloc(:, jdir), grad_rho(:, idir))
731 stress_sr(idir,idir) = stress_sr(idir,idir) + energy_ps_sr
736 stress_sr = stress_sr/ions%latt%rcell_volume
738 safe_deallocate_a(rvloc)
739 safe_deallocate_a(grad_rho)
747 safe_allocate(rho_lr(1:gr%np_part))
748 safe_allocate(rho_lr_x(1:gr%np, 1:gr%der%dim))
751 safe_allocate(rho_local_lr(1:gr%np))
752 do iatom = ions%atoms_dist%start, ions%atoms_dist%end
753 assert(ions%atom(iatom)%species%is_ps())
755 ions%pos(:, iatom), gr, rho_local_lr, nlr_x=rho_lr_x)
759 safe_deallocate_a(rho_local_lr)
761 if (ions%atoms_dist%parallel)
then
766 do idir = 1, ions%space%periodic_dim
767 do jdir = idir, ions%space%periodic_dim
768 stress_lr(idir, jdir) = stress_lr(idir, jdir) +
dmf_dotp(gr, rho_lr_x(:,jdir), grad_vh(:, idir))
771 safe_deallocate_a(rho_lr_x)
773 safe_allocate(vlr(1:gr%np_part))
775 safe_deallocate_a(rho_lr)
777 safe_allocate(grad_vlr(1:gr%np, 1:gr%der%dim))
779 safe_deallocate_a(vlr)
781 do idir = 1, ions%space%periodic_dim
782 do jdir = idir, ions%space%periodic_dim
783 stress_lr(idir, jdir) = stress_lr(idir, jdir) -
dmf_dotp(gr, grad_vh(:,idir), grad_vlr(:, jdir))/
m_two/
m_pi
789 safe_deallocate_a(grad_vlr)
793 if (ions%space%periodic_dim == 3)
then
795 do iatom = 1, ions%natoms
796 charge = charge + ions%atom(iatom)%species%get_zval()
799 do iatom = 1, ions%natoms
800 select type(spec => ions%atom(iatom)%species)
805 do idir = 1, ions%space%periodic_dim
806 stress_lr(idir, idir) = stress_lr(idir, idir) &
807 +
m_two*
m_pi*spec_ps%sigma_erf**2*charge*zi /ions%latt%rcell_volume
813 stress_lr = stress_lr/ions%latt%rcell_volume
815 stress_ps_local = stress_sr + stress_lr
824 class(
mesh_t),
intent(in) :: mesh
825 type(
ions_t),
intent(in) :: ions
826 integer,
intent(in) :: iatom
827 real(real64),
intent(inout) :: vpsl(:)
828 real(real64),
intent(inout) :: rvpsl(:,:)
831 real(real64) :: radius, vl_ip
833 type(
ps_t),
pointer :: ps
837 if (.not. ions%atom(iatom)%species%is_ps())
then
844 select type(spec=>ions%atom(iatom)%species)
849 radius = ps%vl%x_threshold*1.05_real64
851 call submesh_init(sphere, ions%space, mesh, ions%latt, ions%pos(:, iatom), radius)
857 vpsl(sphere%map(ip)) = vpsl(sphere%map(ip)) + vl_ip
858 rvpsl(sphere%map(ip), 1:ions%space%periodic_dim) = rvpsl(sphere%map(ip), 1:ions%space%periodic_dim) &
859 + sphere%rel_x(1:ions%space%periodic_dim, ip) * vl_ip
889 type(
grid_t),
target,
intent(in) :: gr
892 type(
space_t),
intent(in) :: space
893 real(real64),
intent(in) :: rcell_volume
894 real(real64),
intent(out) :: stress_hubbard(3, 3)
896 integer :: ik, ist, idir, jdir
897 integer :: ib, minst, maxst
898 type(
wfs_elec_t) :: psib, rvu_psib(3), gpsib(3)
899 complex(real64),
allocatable :: stress_tmp(:)
910 safe_allocate(stress_tmp(1:st%block_size))
914 do ik = st%d%kpt%start, st%d%kpt%end
918 do ib = st%group%block_start, st%group%block_end
928 do idir = 1, gr%der%dim
929 call psib%copy_to(rvu_psib(idir))
933 call zlda_u_rvu(hm%lda_u, gr, space, hm%d, namespace, psib, rvu_psib)
939 do ist = minst, maxst
940 stress_hubbard(idir, jdir) = stress_hubbard(idir, jdir) &
941 +
m_two * st%kweights(ik) * st%occ(ist, ik) * real(stress_tmp(ist-minst+1), real64)
947 do idir = 1, gr%der%dim
948 call rvu_psib(idir)%end()
949 call gpsib(idir)%end()
955 safe_deallocate_a(stress_tmp)
957 if (st%parallel_in_states .or. st%d%kpt%parallel)
then
965 if (hm%kpoints%use_symmetries)
then
971 stress_hubbard(idir, idir) = stress_hubbard(idir, idir) + hm%energy%int_dft_u
974 stress_hubbard = stress_hubbard/rcell_volume
981 subroutine output_stress(iunit, space_dim, stress_tensors, all_terms)
982 integer,
intent(in) :: iunit
983 integer,
intent(in) :: space_dim
984 type(
stress_t),
intent(in) :: stress_tensors
985 logical,
optional,
intent(in) :: all_terms
987 logical :: write_all_terms
988 character(len=16) :: stress_unit
997 if (write_all_terms)
then
998 write(iunit,
'(3a)')
'Kinetic stress tensor [', trim(stress_unit),
'] ='
1000 if (space_dim == 3)
then
1001 write(iunit,
'(a, es15.6, 3a)')
'Kinetic pressure sumrule violation: ', &
1008 write(iunit,
'(3a)')
'Hartree stress tensor [', trim(stress_unit),
'] ='
1010 if (space_dim == 3)
then
1011 write(iunit,
'(a, es15.6, 3a)')
'Hartree pressure sumrule violation: ', &
1017 write(iunit,
'(3a)')
'XC stress tensor [', trim(stress_unit),
'] ='
1020 write(iunit,
'(3a)')
'XC NLCC stress tensor [', trim(stress_unit),
'] ='
1023 write(iunit,
'(3a)')
'Local pseudo. stress tensor [', trim(stress_unit),
'] ='
1026 write(iunit,
'(3a)')
'Nonlocal pseudo. stress tensor [', trim(stress_unit),
'] ='
1029 write(iunit,
'(3a)')
'Ion-ion stress tensor [', trim(stress_unit),
'] ='
1032 write(iunit,
'(3a)')
'vdW stress tensor [', trim(stress_unit),
'] ='
1035 write(iunit,
'(3a)')
'Hubbard stress tensor [', trim(stress_unit),
'] ='
1039 write(iunit,
'(3a)')
'Total stress tensor [', trim(stress_unit),
'] ='
1047 integer,
intent(in) :: iunit
1048 integer,
intent(in) :: space_dim
1049 real(real64),
intent(in) :: total_stress_tensor(3,3)
1052 real(real64) :: pressure
1053 character(len=16) :: stress_unit
1059 do idim = 1, space_dim
1060 pressure = pressure - total_stress_tensor(idim, idim) / real(space_dim, real64)
1063 write(iunit,
'(3a,es16.8)', advance=
"no")
'Pressure [', trim(stress_unit),
'] = ', &
1065 if (space_dim == 3)
then
1074 integer,
intent(in) :: ounit
1075 integer,
intent(in) :: space_dim
1076 real(real64),
intent(in) :: tensor(3,3)
1078 real(real64) :: tensor_with_unit(3,3)
1079 integer :: idim, jdim
1083 write(ounit,
'(a9,2x)', advance=
"no")
"T_{ij}"
1084 do jdim = 1, space_dim
1085 write(ounit,
'(i18)', advance=
"no") jdim
1088 do idim = 1, space_dim
1089 write(ounit,
'(i9,2x)', advance=
"no") idim
1090 do jdim = 1, space_dim
1091 write(ounit,
'(es18.9)', advance=
"no") tensor_with_unit(idim, jdim)
constant times a vector plus a vector
Copies a vector x, to a vector y.
This module implements common operations on batches of mesh functions.
subroutine, public batch_set_zero(this, np, async)
fill all mesh functions of the batch with zero
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_batch_grad(der, ffb, opffb, ghost_update, set_bc, to_cartesian, metric, factor)
apply the gradient to a batch of mesh functions
integer, parameter, public spinors
subroutine, public energy_calc_total(namespace, space, hm, gr, st, ext_partners, iunit, full)
This subroutine calculates the total energy of the system. Basically, it adds up the KS eigenvalues,...
integer, parameter, public scalar_relativistic_zora
integer, parameter, public fully_relativistic_zora
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
real(real64), parameter, public m_four
real(real64), parameter, public m_pi
some mathematical constants
real(real64), parameter, public m_epsilon
real(real64), parameter, public m_one
This module implements the underlying real-space grid.
subroutine, public hamiltonian_elec_copy_and_set_phase(hm, gr, kpt, psib, psib_with_phase)
Copy a batch to another batch and apply the Bloch phase to it.
This module defines classes and functions for interaction partners.
subroutine, public ion_interaction_stress(this, space, latt, atom, natoms, pos, stress_ii)
Computes the contribution to the stress tensor the ion-ion energy.
A module to handle KS potential, without the external potential.
integer, parameter, public independent_particles
integer, parameter, public generalized_kohn_sham_dft
integer, parameter, public kohn_sham_dft
integer, parameter, public dft_u_none
subroutine, public zlda_u_rvu(this, mesh, space, d, namespace, psib, gpsib)
This routine computes .
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.
subroutine, public dsymmetrize_matrix(nn, aa)
This module defines functions over batches of mesh functions.
subroutine, public zmesh_batch_dotp_vector(mesh, aa, bb, dot, reduce, cproduct)
calculate the vector of dot-products of mesh functions between two batches
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 dpoisson_solve(this, namespace, pot, rho, all_nodes, kernel, reset)
Calculates the Poisson equation. Given the density returns the corresponding potential.
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.
subroutine, public species_get_long_range_density(species, namespace, space, latt, pos, mesh, rho, sphere_inout, nlr_x)
subroutine, public species_get_nlcc(species, space, latt, pos, mesh, rho_core, accumulate)
real(real64) function, public spline_eval(spl, x)
This module handles spin dimensions of the states and the k-point distribution.
integer pure function, public states_elec_block_max(st, ib)
return index of last state in block ib
integer pure function, public states_elec_block_min(st, ib)
return index of first state in block ib
This module implements the calculation of the stress tensor.
subroutine stress_from_kinetic(gr, space, hm, st, symm, rcell_volume, stress_kin)
Computes the contribution to the stress tensor from the kinetic energy.
subroutine stress_from_xc(energy, rcell_volume, periodic_dim, stress_xc)
Computes the contribution to the stress tensor from the xc energy.
subroutine print_stress_tensor(ounit, space_dim, tensor)
subroutine, public output_pressure(iunit, space_dim, total_stress_tensor)
subroutine epot_local_pseudopotential_sr(mesh, ions, iatom, vpsl, rvpsl)
subroutine, public stress_calculate(namespace, gr, hm, st, ions, ks, ext_partners)
This computes the total stress on the lattice.
subroutine stress_from_hubbard(namespace, gr, st, hm, space, rcell_volume, stress_hubbard)
Computes the contribution to the stress tensor from the Hubbard energy.
subroutine stress_from_xc_nlcc(rcell_volume, gr, st, ions, vxc, stress_xc_nlcc)
Computes the NLCC contribution to the stress tensor from the xc energy.
subroutine stress_from_pseudo_nonloc(gr, st, hm, ions, stress_ps_nl)
Computes the contribution to the stress tensor from the nonlocal part of the pseudopotentials.
subroutine stress_from_hartree(gr, space, volume, vh, grad_vh, ehartree, stress_Hartree)
Computes the contribution to the stress tensor from the Hartree energy.
subroutine, public output_stress(iunit, space_dim, stress_tensors, all_terms)
subroutine stress_from_pseudo_local(gr, st, hm, ions, rho_total, vh, grad_vh, stress_ps_local)
Computes the contribution from the local part of the pseudopotential.
subroutine, public submesh_end(this)
subroutine, public submesh_init(this, space, mesh, latt, center, rc)
subroutine, public dsymmetrize_tensor_cart(symm, tensor, use_non_symmorphic)
Symmetric a rank-2 tensor defined in Cartesian space.
type(type_t), parameter, public type_cmplx
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
type(unit_t), public unit_gpa
For output pressure in GPa.
logical pure function, public xc_is_energy_functional(xcs)
Is one of the x or c functional is not an energy functional.
pure logical function, public in_family(family, xc_families)
A module that takes care of xc contribution from vdW interactions.
integer(int64), dimension(5), parameter, public d3_lib_options
VDWCORRECTION options that correspond to the DFT-D3 library.
Description of the grid, containing information on derivatives, stencil, and symmetries.
Describes mesh distribution to nodes.
A type storing the information and data about a pseudopotential.
The states_elec_t class contains all electronic wave functions.
A submesh is a type of mesh, used for the projectors in the pseudopotentials It contains points on a ...
batches of electronic states