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 (ks%vdw%vdw_correction /= option__vdwcorrection__none .and. .not. any(ks%vdw%vdw_correction ==
d3_lib_options))
then
131 write(
message(1),
'(a)')
'The stress tensor is currently only implemented with DFT-D3 vdW correction'
135 if (hm%pcm%run_pcm)
then
139 if (
allocated(hm%v_static))
then
143 if (ks%has_photons)
then
147 if (.not. hm%vnl%apply_projector_matrices)
then
157 if (.not. xc_is_energy_functional(hm%xc))
then
161 if ( .not. in_family(hm%xc%family, [xc_family_lda, xc_family_gga]))
then
162 write(
message(1),
'(a)')
'The stress tensor computation is currently only possible at the Kohn-Sham DFT level'
163 write(
message(2),
'(a)')
'with LDA and GGA functionals or for independent particles.'
167 if (in_family(hm%xc%family, [xc_family_gga]) .and. st%d%ispin ==
spinors)
then
178 safe_allocate(rho_total(1:gr%np_part))
180 rho_total(ip) = sum(st%rho(ip, 1:st%d%nspin))
189 safe_allocate(vh(1:gr%np_part))
190 safe_allocate(grad_vh(1:gr%np, 1:gr%der%dim))
192 call lalg_copy(gr%np, hm%ks_pot%vhartree, vh)
194 if (hm%d%spin_channels > 1)
then
195 safe_allocate(rho(1:gr%np_part))
203 if (hm%d%spin_channels > 1)
then
204 safe_deallocate_p(rho)
209 ehartree = hm%energy%hartree
216 call stress_from_kinetic(gr, ions%space, hm, st, gr%symm, ions%latt%rcell_volume, stress_kin)
217 stress = stress + stress_kin
223 call stress_from_hartree(gr, ions%space, ions%latt%rcell_volume, vh, grad_vh, ehartree, stress_hartree)
224 stress = stress + stress_hartree
226 call stress_from_xc(hm%energy, ions%latt%rcell_volume, ions%space%periodic_dim, stress_xc)
229 if (
allocated(st%rho_core))
then
230 call stress_from_xc_nlcc(ions%latt%rcell_volume, gr, st, ions, hm%ks_pot%vxc, stress_xc_nlcc)
231 stress_xc = stress_xc + stress_xc_nlcc
234 stress_xc = stress_xc + ks%stress_xc_gga / ions%latt%rcell_volume
235 stress = stress + stress_xc
239 stress_ps = stress_ps_local
240 stress = stress + stress_ps_local
242 safe_deallocate_a(vh)
243 safe_deallocate_a(grad_vh)
246 stress_ps = stress_ps + stress_ps_nl
247 stress = stress + stress_ps_nl
249 call stress_from_hubbard(namespace, gr, st, hm, ions%space, ions%latt%rcell_volume, stress_hubbard)
250 stress = stress + stress_hubbard
252 call ion_interaction_stress(ions%ion_interaction, ions%space, ions%latt, ions%atom, ions%natoms, ions%pos, stress_ii)
253 stress = stress + stress_ii
260 st%stress_tensors%kinetic = stress_kin
261 st%stress_tensors%Hartree = stress_hartree
262 st%stress_tensors%xc = stress_xc
263 st%stress_tensors%ps_local = stress_ps_local
264 st%stress_tensors%ps_nl = stress_ps_nl
265 st%stress_tensors%hubbard = stress_hubbard
266 st%stress_tensors%ion_ion = stress_ii
269 if (ks%vdw%vdw_correction /= option__vdwcorrection__none)
then
270 st%stress_tensors%vdw = hm%ep%vdw_stress
272 st%stress_tensors%vdw =
m_zero
274 stress = stress + st%stress_tensors%vdw
277 if (hm%kpoints%use_symmetries)
then
283 st%stress_tensors%total = stress
287 st%stress_tensors%kinetic_sumrule =
m_zero
289 st%stress_tensors%Hartree_sumrule =
m_zero
290 if(ions%space%periodic_dim == 3)
then
291 st%stress_tensors%kinetic_sumrule = (stress_kin(1,1) + stress_kin(2,2) + stress_kin(3,3))*ions%latt%rcell_volume
292 st%stress_tensors%kinetic_sumrule = st%stress_tensors%kinetic_sumrule -
m_two * hm%energy%kinetic
294 st%stress_tensors%hartree_sumrule = (stress_hartree(1,1) + stress_hartree(2,2) + stress_hartree(3,3))*ions%latt%rcell_volume
295 st%stress_tensors%hartree_sumrule = st%stress_tensors%hartree_sumrule - hm%energy%hartree
298 safe_deallocate_a(rho_total)
320 type(
grid_t),
intent(in) :: gr
321 class(
space_t),
intent(in) :: space
325 real(real64),
intent(in) :: rcell_volume
326 real(real64),
intent(out) :: stress_kin(3, 3)
328 integer :: ik, ist, idir, jdir, ib, minst, maxst
329 complex(real64),
allocatable :: stress_l_block(:)
337 safe_allocate(stress_l_block(1:st%block_size))
339 do ik = st%d%kpt%start, st%d%kpt%end
342 do ib = st%group%block_start, st%group%block_end
352 do idir = 1, space%periodic_dim
353 do jdir = idir, space%periodic_dim
356 do ist = minst, maxst
357 stress_kin(idir,jdir) = stress_kin(idir,jdir) &
358 + st%kweights(ik) * st%occ(ist, ik) &
359 * real(stress_l_block(ist - minst + 1), real64)
364 do idir = 1, space%dim
365 call gpsib(idir)%end()
372 if (st%parallel_in_states .or. st%d%kpt%parallel)
then
381 if (hm%kpoints%use_symmetries)
then
385 stress_kin = stress_kin / rcell_volume
409 type(
grid_t),
intent(in) :: gr
410 class(
space_t),
intent(in) :: space
411 real(real64),
intent(in) :: volume
412 real(real64),
intent(in) :: vh(:)
413 real(real64),
intent(in) :: grad_vh(:,:)
414 real(real64),
intent(in) :: ehartree
415 real(real64),
intent(out) :: stress_Hartree(3, 3)
417 integer :: idir, jdir
422 stress_hartree(:,:) =
m_zero
424 do idir = 1, space%periodic_dim
425 do jdir = idir, space%periodic_dim
426 stress_hartree(idir, jdir) = -
dmf_dotp(gr, grad_vh(:,idir), grad_vh(:, jdir))/
m_four/
m_pi
428 stress_hartree(idir, idir) = stress_hartree(idir, idir) + ehartree
433 stress_hartree = stress_hartree/volume
457 subroutine stress_from_xc(energy, rcell_volume, periodic_dim, stress_xc)
458 type(
energy_t),
intent(in) :: energy
459 real(real64),
intent(in) :: rcell_volume
460 integer,
intent(in) :: periodic_dim
461 real(real64),
intent(out) :: stress_xc(3, 3)
469 do idir = 1, periodic_dim
470 stress_xc(idir, idir) = - energy%exchange - energy%correlation + energy%intnvxc
472 stress_xc(:,:) = stress_xc(:,:) / rcell_volume
487 real(real64),
intent(in) :: rcell_volume
488 type(
grid_t),
intent(in) :: gr
490 type(
ions_t),
intent(in) :: ions
491 real(real64),
intent(in) :: vxc(:,:)
492 real(real64),
intent(out) :: stress_xc_nlcc(3, 3)
494 integer :: idir, jdir, iat
495 real(real64),
allocatable :: gnlcc(:,:), gnlcc_x(:,:,:), vxc_tot(:)
500 assert(
allocated(st%rho_core))
505 safe_allocate(gnlcc(gr%np, gr%der%dim))
506 safe_allocate(gnlcc_x(gr%np, gr%der%dim, gr%der%dim))
508 do iat = ions%atoms_dist%start, ions%atoms_dist%end
509 assert(ions%atom(iat)%species%is_ps())
511 ions%pos(:,iat), gr, gnlcc, gnlcc_x)
513 safe_deallocate_a(gnlcc)
515 if (ions%atoms_dist%parallel)
then
520 safe_allocate(vxc_tot(1:gr%np))
521 vxc_tot = vxc(1:gr%np, 1)
522 if(st%d%nspin > 1) vxc_tot = vxc_tot + vxc(1:gr%np, 2)
524 do idir = 1, ions%space%periodic_dim
525 do jdir = idir, ions%space%periodic_dim
526 stress_xc_nlcc(idir, jdir) =
dmf_dotp(gr, vxc_tot, gnlcc_x(:,idir, jdir))
529 safe_deallocate_a(vxc_tot)
530 safe_deallocate_a(gnlcc_x)
534 stress_xc_nlcc(:,:) = stress_xc_nlcc(:,:) / rcell_volume
560 type(
grid_t),
target,
intent(in) :: gr
563 type(
ions_t),
intent(in) :: ions
564 real(real64),
intent(out) :: stress_ps_nl(3, 3)
566 integer :: ik, ist, idir, jdir
567 integer :: ib, minst, maxst
568 type(
wfs_elec_t) :: psib, rvnl_psib(3), gpsib(3)
569 complex(real64),
allocatable :: stress_tmp(:)
576 safe_allocate(stress_tmp(1:st%block_size))
580 do ik = st%d%kpt%start, st%d%kpt%end
584 do ib = st%group%block_start, st%group%block_end
595 do idir = 1, gr%der%dim
596 call psib%copy_to(rvnl_psib(idir))
599 call hm%vnl%zr_vn_local(gr, st%d, gr%der%boundaries%spiral, psib, rvnl_psib)
601 do idir = 1, ions%space%periodic_dim
602 do jdir = idir, ions%space%periodic_dim
605 do ist = minst, maxst
606 stress_ps_nl(idir, jdir) = stress_ps_nl(idir, jdir) &
607 +
m_two * st%kweights(ik) * st%occ(ist, ik) * real(stress_tmp(ist-minst+1), real64)
613 do idir = 1, gr%der%dim
614 call rvnl_psib(idir)%end()
615 call gpsib(idir)%end()
621 safe_deallocate_a(stress_tmp)
623 if (st%parallel_in_states .or. st%d%kpt%parallel)
then
631 if (hm%kpoints%use_symmetries)
then
636 do idir = 1, ions%space%periodic_dim
637 stress_ps_nl(idir, idir) = stress_ps_nl(idir, idir) + hm%energy%extern_non_local
640 stress_ps_nl = stress_ps_nl/ions%latt%rcell_volume
666 type(
grid_t),
target,
intent(in) :: gr
669 type(
ions_t),
intent(in) :: ions
670 real(real64),
contiguous,
intent(inout) :: rho_total(:)
671 real(real64),
intent(in) :: vh(:)
672 real(real64),
intent(in) :: grad_vh(:,:)
673 real(real64),
intent(out) :: stress_ps_local(3, 3)
676 real(real64) :: stress_SR(3, 3), stress_LR(3, 3)
677 real(real64) :: energy_ps_SR, charge, zi
678 real(real64),
allocatable :: vloc(:), rvloc(:,:), rho_local_lr(:), rho_lr(:)
679 real(real64),
allocatable :: grad_rho(:,:), rho_lr_x(:,:), vlr(:), grad_vlr(:,:)
680 integer :: idir, jdir, iatom
681 type(
ps_t),
pointer :: spec_ps
689 safe_allocate(vloc(1:gr%np))
691 safe_allocate(rvloc(1:gr%np, 1:gr%der%dim))
693 do iatom = 1, ions%natoms
696 safe_deallocate_a(vloc)
698 safe_allocate(grad_rho(1:gr%np,1:gr%der%dim))
701 energy_ps_sr = hm%energy%extern_local
702 do idir = 1, ions%space%periodic_dim
703 do jdir = idir, ions%space%periodic_dim
704 stress_sr(idir, jdir) = stress_sr(idir, jdir) &
705 +
dmf_dotp(gr, rvloc(:, jdir), grad_rho(:, idir))
707 stress_sr(idir,idir) = stress_sr(idir,idir) + energy_ps_sr
712 stress_sr = stress_sr/ions%latt%rcell_volume
714 safe_deallocate_a(rvloc)
715 safe_deallocate_a(grad_rho)
723 safe_allocate(rho_lr(1:gr%np_part))
724 safe_allocate(rho_lr_x(1:gr%np, 1:gr%der%dim))
727 safe_allocate(rho_local_lr(1:gr%np))
728 do iatom = ions%atoms_dist%start, ions%atoms_dist%end
729 assert(ions%atom(iatom)%species%is_ps())
731 ions%pos(:, iatom), gr, rho_local_lr, nlr_x=rho_lr_x)
735 safe_deallocate_a(rho_local_lr)
737 if (ions%atoms_dist%parallel)
then
742 do idir = 1, ions%space%periodic_dim
743 do jdir = idir, ions%space%periodic_dim
744 stress_lr(idir, jdir) = stress_lr(idir, jdir) +
dmf_dotp(gr, rho_lr_x(:,jdir), grad_vh(:, idir))
747 safe_deallocate_a(rho_lr_x)
749 safe_allocate(vlr(1:gr%np_part))
751 safe_deallocate_a(rho_lr)
753 safe_allocate(grad_vlr(1:gr%np, 1:gr%der%dim))
755 safe_deallocate_a(vlr)
757 do idir = 1, ions%space%periodic_dim
758 do jdir = idir, ions%space%periodic_dim
759 stress_lr(idir, jdir) = stress_lr(idir, jdir) -
dmf_dotp(gr, grad_vh(:,idir), grad_vlr(:, jdir))/
m_two/
m_pi
765 safe_deallocate_a(grad_vlr)
769 if (ions%space%periodic_dim == 3)
then
771 do iatom = 1, ions%natoms
772 charge = charge + ions%atom(iatom)%species%get_zval()
775 do iatom = 1, ions%natoms
776 select type(spec => ions%atom(iatom)%species)
781 do idir = 1, ions%space%periodic_dim
782 stress_lr(idir, idir) = stress_lr(idir, idir) &
783 +
m_two*
m_pi*spec_ps%sigma_erf**2*charge*zi /ions%latt%rcell_volume
789 stress_lr = stress_lr/ions%latt%rcell_volume
791 stress_ps_local = stress_sr + stress_lr
800 class(
mesh_t),
intent(in) :: mesh
801 type(
ions_t),
intent(in) :: ions
802 integer,
intent(in) :: iatom
803 real(real64),
intent(inout) :: vpsl(:)
804 real(real64),
intent(inout) :: rvpsl(:,:)
807 real(real64) :: radius, vl_ip
809 type(
ps_t),
pointer :: ps
813 if (.not. ions%atom(iatom)%species%is_ps())
then
820 select type(spec=>ions%atom(iatom)%species)
825 radius = ps%vl%x_threshold*1.05_real64
827 call submesh_init(sphere, ions%space, mesh, ions%latt, ions%pos(:, iatom), radius)
833 vpsl(sphere%map(ip)) = vpsl(sphere%map(ip)) + vl_ip
834 rvpsl(sphere%map(ip), 1:ions%space%periodic_dim) = rvpsl(sphere%map(ip), 1:ions%space%periodic_dim) &
835 + sphere%rel_x(1:ions%space%periodic_dim, ip) * vl_ip
865 type(
grid_t),
target,
intent(in) :: gr
868 type(
space_t),
intent(in) :: space
869 real(real64),
intent(in) :: rcell_volume
870 real(real64),
intent(out) :: stress_hubbard(3, 3)
872 integer :: ik, ist, idir, jdir
873 integer :: ib, minst, maxst
874 type(
wfs_elec_t) :: psib, rvu_psib(3), gpsib(3)
875 complex(real64),
allocatable :: stress_tmp(:)
886 safe_allocate(stress_tmp(1:st%block_size))
890 do ik = st%d%kpt%start, st%d%kpt%end
894 do ib = st%group%block_start, st%group%block_end
904 do idir = 1, gr%der%dim
905 call psib%copy_to(rvu_psib(idir))
909 call zlda_u_rvu(hm%lda_u, gr, space, hm%d, namespace, psib, rvu_psib)
915 do ist = minst, maxst
916 stress_hubbard(idir, jdir) = stress_hubbard(idir, jdir) &
917 +
m_two * st%kweights(ik) * st%occ(ist, ik) * real(stress_tmp(ist-minst+1), real64)
923 do idir = 1, gr%der%dim
924 call rvu_psib(idir)%end()
925 call gpsib(idir)%end()
931 safe_deallocate_a(stress_tmp)
933 if (st%parallel_in_states .or. st%d%kpt%parallel)
then
941 if (hm%kpoints%use_symmetries)
then
947 stress_hubbard(idir, idir) = stress_hubbard(idir, idir) + hm%energy%int_dft_u
950 stress_hubbard = stress_hubbard/rcell_volume
957 subroutine output_stress(iunit, space_dim, stress_tensors, all_terms)
958 integer,
intent(in) :: iunit
959 integer,
intent(in) :: space_dim
960 type(
stress_t),
intent(in) :: stress_tensors
961 logical,
optional,
intent(in) :: all_terms
963 logical :: write_all_terms
964 character(len=16) :: stress_unit
973 if (write_all_terms)
then
974 write(iunit,
'(3a)')
'Kinetic stress tensor [', trim(stress_unit),
'] ='
976 if (space_dim == 3)
then
977 write(iunit,
'(a, es15.6, 3a)')
'Kinetic pressure sumrule violation: ', &
984 write(iunit,
'(3a)')
'Hartree stress tensor [', trim(stress_unit),
'] ='
986 if (space_dim == 3)
then
987 write(iunit,
'(a, es15.6, 3a)')
'Hartree pressure sumrule violation: ', &
993 write(iunit,
'(3a)')
'XC stress tensor [', trim(stress_unit),
'] ='
996 write(iunit,
'(3a)')
'Local pseudo. stress tensor [', trim(stress_unit),
'] ='
999 write(iunit,
'(3a)')
'Nonlocal pseudo. stress tensor [', trim(stress_unit),
'] ='
1002 write(iunit,
'(3a)')
'Ion-ion stress tensor [', trim(stress_unit),
'] ='
1005 write(iunit,
'(3a)')
'vdW stress tensor [', trim(stress_unit),
'] ='
1008 write(iunit,
'(3a)')
'Hubbard stress tensor [', trim(stress_unit),
'] ='
1012 write(iunit,
'(3a)')
'Total stress tensor [', trim(stress_unit),
'] ='
1020 integer,
intent(in) :: iunit
1021 integer,
intent(in) :: space_dim
1022 real(real64),
intent(in) :: total_stress_tensor(3,3)
1025 real(real64),
parameter :: au_to_GPa = 29421.02648438959_real64
1028 real(real64) :: pressure =
m_zero
1029 character(len=16) :: stress_unit
1034 do idim = 1, space_dim
1035 pressure = pressure - total_stress_tensor(idim, idim) / real(space_dim, real64)
1038 write(iunit,
'(3a,es16.8)', advance=
"no")
'Pressure [', trim(stress_unit),
'] = ', &
1040 if (space_dim == 3)
then
1041 write(iunit,
'(2x,a,f16.8)')
'Pressure [GPa] = ', pressure * au_to_gpa
1049 integer,
intent(in) :: ounit
1050 integer,
intent(in) :: space_dim
1051 real(real64),
intent(in) :: tensor(3,3)
1053 real(real64) :: tensor_with_unit(3,3)
1054 integer :: idim, jdim
1058 write(ounit,
'(a9,2x)', advance=
"no")
"T_{ij}"
1059 do jdim = 1, space_dim
1060 write(ounit,
'(i18)', advance=
"no") jdim
1063 do idim = 1, space_dim
1064 write(ounit,
'(i9,2x)', advance=
"no") idim
1065 do jdim = 1, space_dim
1066 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_grad(species, space, latt, pos, mesh, rho_core_grad, gnlcc_x)
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), 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
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