Octopus
stress.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2016 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2023 N. Tancogne-Dejean
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
22! ---------------------------------------------------------
25module stress_oct_m
28 use comm_oct_m
29 use debug_oct_m
33 use energy_oct_m
35 use epot_oct_m
36 use global_oct_m
37 use grid_oct_m
39 use io_oct_m
41 use ions_oct_m
43 use, intrinsic :: iso_fortran_env
46 use lda_u_oct_m
48 use math_oct_m
49 use mesh_oct_m
53 use mpi_oct_m
58 use ps_oct_m
60 use space_oct_m
70 use types_oct_m
71 use unit_oct_m
73 use v_ks_oct_m
75 use xc_oct_m
76 use xc_f03_lib_m
78 implicit none
79
80 private
81 public :: &
85
86contains
87
88 ! ---------------------------------------------------------
90 subroutine stress_calculate(namespace, gr, hm, st, ions, ks, ext_partners)
91 type(namespace_t), intent(in) :: namespace
92 type(grid_t), intent(inout) :: gr
93 type(hamiltonian_elec_t), intent(inout) :: hm
94 type(states_elec_t), intent(inout) :: st
95 type(ions_t), intent(inout) :: ions
96 type(v_ks_t), intent(in) :: ks
97 type(partner_list_t), intent(in) :: ext_partners
98
99 real(real64), allocatable :: rho_total(:)
100 real(real64) :: stress(3,3) ! stress tensor in Cartesian coordinate
101 real(real64) :: stress_kin(3,3), stress_Hartree(3,3), stress_xc(3,3), stress_xc_nlcc(3,3)
102 real(real64) :: stress_ps(3,3), stress_ps_nl(3,3), stress_ps_local(3,3), stress_ii(3,3)
103 real(real64) :: stress_hubbard(3,3)
104 integer :: ip
105 real(real64), allocatable :: vh(:)
106 real(real64), allocatable :: grad_vh(:,:)
107 real(real64) :: ehartree
108
109 call profiling_in("STRESS_CALCULATE")
110 push_sub(stress_calculate)
111
112 if (st%wfs_type /= type_cmplx) then
113 write(message(1),'(a)') 'The stress tensors for real wavefunctions has not been implemented!'
114
115 if (hm%kpoints%full%npoints == 1) then
116 write(message(2),'(a)') 'For testing this feature, you can add ForceComplex=yes to the input file'
117 call messages_fatal(2, namespace=namespace)
118 end if
119
120 call messages_fatal(1, namespace=namespace)
121 end if
122
123 if (ions%space%periodic_dim == 1) then
124 call messages_not_implemented("Stress tensor for 1D periodic systems", namespace=namespace)
125 end if
126
127 if ( .not. (ks%theory_level == kohn_sham_dft .and. in_family(hm%xc%family, [xc_family_lda, xc_family_gga])) &
128 .and. .not. (ks%theory_level == independent_particles)) then
129 write(message(1),'(a)') 'The stress tensor computation is currently only possible at the Kohn-Sham DFT level'
130 write(message(2),'(a)') 'with LDA and GGA functionals or for independent particles.'
131 call messages_fatal(2, namespace=namespace)
132 end if
133
134 if (in_family(hm%xc%family, [xc_family_gga]) .and. st%d%ispin == spinors) then
135 call messages_not_implemented("Stress tensor for GGAs with spinors", namespace=namespace)
136 end if
137
138 if (ks%vdw%vdw_correction /= option__vdwcorrection__none .and. .not. any(ks%vdw%vdw_correction == d3_lib_options)) then
139 write(message(1),'(a)') 'The stress tensor is currently only implemented with DFT-D3 vdW correction'
140 call messages_fatal(1, namespace=namespace)
141 end if
142
143 if (hm%pcm%run_pcm) then
144 call messages_not_implemented('Stress tensor with PCM')
145 end if
146
147 if (allocated(hm%v_static)) then
148 call messages_not_implemented('Stress tensor with static electric fields')
149 end if
150
151 if (ks%has_photons) then
152 call messages_not_implemented('Stress tensor with photon modes')
153 end if
154
155 if (.not. hm%vnl%apply_projector_matrices) then
156 call messages_not_implemented('Stress tensor with relativistic Kleinman-Bylander pseudopotential')
157 end if
158
159 if (hm%ep%reltype == scalar_relativistic_zora .or. hm%ep%reltype == fully_relativistic_zora) then
160 call messages_not_implemented('Stress tensor with ZORA')
161 end if
162
163 if (.not. xc_is_energy_functional(hm%xc)) then
164 call messages_not_implemented("Stress tensor with xc functionals that are not energy functionals")
165 end if
166
167 stress(:,:) = m_zero
168
169 safe_allocate(rho_total(1:gr%np_part))
170 do ip = 1, gr%np
171 rho_total(ip) = sum(st%rho(ip, 1:st%d%nspin))
172 end do
173
174 ! As we rely on some of the full energy components, we need to recompute it first
175 ! TODO: We should restrict the components of the energy needed to be computed
176 call energy_calc_total(namespace, ions%space, hm, gr, st, ext_partners, iunit = 0, full = .true.)
177
178 ! In order to get the electrostatic part (Hartree and local pseudopotential part),
179 ! we need to get the Hartree potential and its gradient
180 safe_allocate(vh(1:gr%np_part))
181 safe_allocate(grad_vh(1:gr%np, 1:gr%der%dim))
182 call lalg_copy(gr%np, hm%vhartree, vh)
183 ehartree = hm%energy%hartree
184
185 ! We also compute the gradient here
186 call dderivatives_grad(gr%der, vh, grad_vh)
187
188 ! We now compute the various contributions to the stress tensor
189
190 ! Stress from kinetic energy of electrons
191 call stress_from_kinetic(gr, ions%space, hm, st, gr%symm, ions%latt%rcell_volume, stress_kin)
192 stress = stress + stress_kin
193
194 if (ks%theory_level == independent_particles) then
195 stress_hartree = m_zero
196 stress_xc = m_zero
197 else
198 call stress_from_hartree(gr, ions%space, ions%latt%rcell_volume, vh, grad_vh, ehartree, stress_hartree)
199 stress = stress + stress_hartree
200
201 call stress_from_xc(hm%energy, ions%latt%rcell_volume, ions%space%periodic_dim, stress_xc)
202
203 ! Nonlinear core correction contribution
204 if (allocated(st%rho_core)) then
205 call stress_from_xc_nlcc(ions%latt%rcell_volume, gr, st, ions, hm%vxc, stress_xc_nlcc)
206 stress_xc = stress_xc + stress_xc_nlcc
207 end if
208 ! Adds the beyond LDA contribution to the stress tensor
209 stress_xc = stress_xc + ks%stress_xc_gga / ions%latt%rcell_volume
210 stress = stress + stress_xc
211
212 end if
213
214 call stress_from_pseudo_local(gr, st, hm, ions, rho_total, vh, grad_vh, stress_ps_local)
215 stress_ps = stress_ps_local
216 stress = stress + stress_ps_local
217
218 safe_deallocate_a(vh)
219 safe_deallocate_a(grad_vh)
220
221 call stress_from_pseudo_nonloc(gr, st, hm, ions, stress_ps_nl)
222 stress_ps = stress_ps + stress_ps_nl
223 stress = stress + stress_ps_nl
224
225 call stress_from_hubbard(namespace, gr, st, hm, ions%space, ions%latt%rcell_volume, stress_hubbard)
226 stress = stress + stress_hubbard
227
228 call ion_interaction_stress(ions%ion_interaction, ions%space, ions%latt, ions%atom, ions%natoms, ions%pos, stress_ii)
229 stress = stress + stress_ii
230 ! Stress from kinetic energy of ion
231 ! Stress from ion-field interaction
232
233 ! Sign changed to fit conventional definition
234 stress = -stress
235
236 st%stress_tensors%kinetic = stress_kin
237 st%stress_tensors%Hartree = stress_hartree
238 st%stress_tensors%xc = stress_xc
239 st%stress_tensors%ps_local = stress_ps_local
240 st%stress_tensors%ps_nl = stress_ps_nl
241 st%stress_tensors%hubbard = stress_hubbard
242 st%stress_tensors%ion_ion = stress_ii
243
244 ! Stress contribution from vdW D3
245 if (ks%vdw%vdw_correction /= option__vdwcorrection__none) then
246 st%stress_tensors%vdw = hm%ep%vdw_stress
247 else
248 st%stress_tensors%vdw = m_zero
249 end if
250 stress = stress + st%stress_tensors%vdw
251
252 ! Symmetrize the stress tensor if we use k-point symmetries
253 if (hm%kpoints%use_symmetries) then
254 call dsymmetrize_tensor_cart(gr%symm, stress, use_non_symmorphic=.true.)
255 ! We guarantee that the matrix is truely symmetric. There could be small numerical assymetries after symmetrization
256 call dsymmetrize_matrix(ions%space%periodic_dim, stress)
257 end if
258
259 st%stress_tensors%total = stress
260
261 ! Some sumrule for validation
262 ! Sumrule is -3P_{kin}\Omega = 2 E_{kin}
263 st%stress_tensors%kinetic_sumrule = m_zero
264 ! Sumrule is -3P_{Hartree}\Omega = E_{Hartree}
265 st%stress_tensors%Hartree_sumrule = m_zero
266 if(ions%space%periodic_dim == 3) then
267 st%stress_tensors%kinetic_sumrule = (stress_kin(1,1) + stress_kin(2,2) + stress_kin(3,3))*ions%latt%rcell_volume
268 st%stress_tensors%kinetic_sumrule = st%stress_tensors%kinetic_sumrule - m_two * hm%energy%kinetic
269
270
271 st%stress_tensors%hartree_sumrule = (stress_hartree(1,1) + stress_hartree(2,2) + stress_hartree(3,3))*ions%latt%rcell_volume
272 st%stress_tensors%hartree_sumrule = st%stress_tensors%hartree_sumrule - hm%energy%hartree
273 end if
274
275 safe_deallocate_a(rho_total)
276
277 pop_sub(stress_calculate)
278 call profiling_out("STRESS_CALCULATE")
279 end subroutine stress_calculate
280
281 ! -------------------------------------------------------
296 subroutine stress_from_kinetic(gr, space, hm, st, symm, rcell_volume, stress_kin)
297 type(grid_t), intent(in) :: gr
298 class(space_t), intent(in) :: space
299 type(hamiltonian_elec_t), intent(in) :: hm
300 type(states_elec_t), intent(inout) :: st
301 type(symmetries_t), intent(in) :: symm
302 real(real64), intent(in) :: rcell_volume
303 real(real64), intent(out) :: stress_kin(3, 3)
304
305 integer :: ik, ist, idir, jdir, ib, minst, maxst
306 complex(real64), allocatable :: stress_l_block(:)
307 type(wfs_elec_t) :: psib, gpsib(space%dim)
308
309 call profiling_in("STRESS_FROM_KINETIC")
310 push_sub(stress_from_kinetic)
311
312 stress_kin(:,:) = m_zero
313
314 safe_allocate(stress_l_block(1:st%block_size))
315
316 do ik = st%d%kpt%start, st%d%kpt%end
317 if (st%kweights(ik) <= m_epsilon) cycle
318
319 do ib = st%group%block_start, st%group%block_end
320 minst = states_elec_block_min(st, ib)
321 maxst = states_elec_block_max(st, ib)
322
323 call hamiltonian_elec_copy_and_set_phase(hm, gr, st%d%kpt, st%group%psib(ib, ik), psib)
324
325 ! calculate the gradient
326 call zderivatives_batch_grad(gr%der, psib, gpsib, set_bc=.false.)
327
328 ! Accumulate the result
329 do idir = 1, space%periodic_dim
330 do jdir = idir, space%periodic_dim
331 call zmesh_batch_dotp_vector(gr, gpsib(idir), gpsib(jdir), stress_l_block)
332
333 do ist = minst, maxst
334 stress_kin(idir,jdir) = stress_kin(idir,jdir) &
335 + st%kweights(ik) * st%occ(ist, ik) &
336 * real(stress_l_block(ist - minst + 1), real64)
337 end do
338 end do
339 end do
340
341 do idir = 1, space%dim
342 call gpsib(idir)%end()
343 end do
344 call psib%end()
345
346 end do
347 end do
348
349 if (st%parallel_in_states .or. st%d%kpt%parallel) then
350 call comm_allreduce(st%st_kpt_mpi_grp, stress_kin)
351 end if
352
353
354 ! Symmetrize the kinetic stress tensor
355 call upper_triangular_to_hermitian(space%periodic_dim, stress_kin)
356
357 ! Symmetrize the stress tensor if we use k-point symmetries
358 if (hm%kpoints%use_symmetries) then
359 call dsymmetrize_tensor_cart(symm, stress_kin, use_non_symmorphic=.true.)
360 end if
361
362 stress_kin = stress_kin / rcell_volume
363
364 call profiling_out("STRESS_FROM_KINETIC")
365 pop_sub(stress_from_kinetic)
366 end subroutine stress_from_kinetic
367
368 ! -------------------------------------------------------
385 subroutine stress_from_hartree(gr, space, volume, vh, grad_vh, ehartree, stress_Hartree)
386 type(grid_t), intent(in) :: gr
387 class(space_t), intent(in) :: space
388 real(real64), intent(in) :: volume
389 real(real64), intent(in) :: vh(:)
390 real(real64), intent(in) :: grad_vh(:,:)
391 real(real64), intent(in) :: ehartree
392 real(real64), intent(out) :: stress_Hartree(3, 3)
393
394 integer :: idir, jdir
395
396 call profiling_in("STRESS_FROM_HARTREE")
397 push_sub(stress_from_hartree)
398
399 stress_hartree(:,:) = m_zero
400
401 do idir = 1, space%periodic_dim
402 do jdir = idir, space%periodic_dim
403 stress_hartree(idir, jdir) = -dmf_dotp(gr, grad_vh(:,idir), grad_vh(:, jdir))/m_four/m_pi
404 end do
405 stress_hartree(idir, idir) = stress_hartree(idir, idir) + ehartree
406 end do
407
408 call upper_triangular_to_hermitian(space%periodic_dim, stress_hartree)
409
410 stress_hartree = stress_hartree/volume
411
412 call profiling_out("STRESS_FROM_HARTREE")
413 pop_sub(stress_from_hartree)
414 end subroutine stress_from_hartree
415
416
417 ! -------------------------------------------------------
431 !
432 ! Note: We assume hm%energy%echange, correlation, and intnvxc
433 ! have already been calculated somewhere else.
434 subroutine stress_from_xc(energy, rcell_volume, periodic_dim, stress_xc)
435 type(energy_t), intent(in) :: energy
436 real(real64), intent(in) :: rcell_volume
437 integer, intent(in) :: periodic_dim
438 real(real64), intent(out) :: stress_xc(3, 3)
439
440 integer :: idir
441
442 call profiling_in("STRESS_FROM_XC")
443 push_sub(stress_from_xc)
444
445 stress_xc = m_zero
446 do idir = 1, periodic_dim
447 stress_xc(idir, idir) = - energy%exchange - energy%correlation + energy%intnvxc
448 end do
449 stress_xc(:,:) = stress_xc(:,:) / rcell_volume
450
451 call profiling_out("STRESS_FROM_XC")
452 pop_sub(stress_from_xc)
453 end subroutine stress_from_xc
454
455
456 ! -------------------------------------------------------
463 subroutine stress_from_xc_nlcc(rcell_volume, gr, st, ions, vxc, stress_xc_nlcc)
464 real(real64), intent(in) :: rcell_volume
465 type(grid_t), intent(in) :: gr
466 type(states_elec_t), intent(in) :: st
467 type(ions_t), intent(in) :: ions
468 real(real64), intent(in) :: vxc(:,:)
469 real(real64), intent(out) :: stress_xc_nlcc(3, 3)
470
471 integer :: idir, jdir, iat
472 real(real64), allocatable :: gnlcc(:,:), gnlcc_x(:,:,:), vxc_tot(:)
473
474 call profiling_in("STRESS_FROM_XC_NLCC")
475 push_sub(stress_from_xc_nlcc)
476
477 assert(allocated(st%rho_core))
479 stress_xc_nlcc = m_zero
480
481 ! We first accumulate the contribution from all the pseudo-ions
482 safe_allocate(gnlcc(gr%np, gr%der%dim))
483 safe_allocate(gnlcc_x(gr%np, gr%der%dim, gr%der%dim))
484 gnlcc_x = m_zero
485 do iat = ions%atoms_dist%start, ions%atoms_dist%end
486 assert(ions%atom(iat)%species%is_ps())
487 call species_get_nlcc_grad(ions%atom(iat)%species, ions%space, ions%latt, &
488 ions%pos(:,iat), gr, gnlcc, gnlcc_x)
489 end do
490 safe_deallocate_a(gnlcc)
491
492 if (ions%atoms_dist%parallel) then
493 call comm_allreduce(ions%atoms_dist%mpi_grp, gnlcc_x)
494 end if
495
496 ! Sum over spin of the xc potential
497 safe_allocate(vxc_tot(1:gr%np))
498 vxc_tot = vxc(1:gr%np, 1)
499 if(st%d%nspin > 1) vxc_tot = vxc_tot + vxc(1:gr%np, 2)
500
501 do idir = 1, ions%space%periodic_dim
502 do jdir = idir, ions%space%periodic_dim
503 stress_xc_nlcc(idir, jdir) = dmf_dotp(gr, vxc_tot, gnlcc_x(:,idir, jdir))
504 end do
505 end do
506 safe_deallocate_a(vxc_tot)
507 safe_deallocate_a(gnlcc_x)
508
509 call upper_triangular_to_hermitian(ions%space%periodic_dim, stress_xc_nlcc)
510
511 stress_xc_nlcc(:,:) = stress_xc_nlcc(:,:) / rcell_volume
512
513 call profiling_out("STRESS_FROM_XC_NLCC")
514 pop_sub(stress_from_xc_nlcc)
515 end subroutine stress_from_xc_nlcc
516
517 ! -------------------------------------------------------
536 subroutine stress_from_pseudo_nonloc(gr, st, hm, ions, stress_ps_nl)
537 type(grid_t), target, intent(in) :: gr
538 type(states_elec_t), intent(inout) :: st
539 type(hamiltonian_elec_t), intent(in) :: hm
540 type(ions_t), intent(in) :: ions
541 real(real64), intent(out) :: stress_ps_nl(3, 3)
542
543 integer :: ik, ist, idir, jdir
544 integer :: ib, minst, maxst
545 type(wfs_elec_t) :: psib, rvnl_psib(3), gpsib(3)
546 complex(real64), allocatable :: stress_tmp(:)
547
548 call profiling_in("STRESS_FROM_PSEUDO_NL")
550
551 assert(st%wfs_type == type_cmplx)
552
553 safe_allocate(stress_tmp(1:st%block_size))
554
555 stress_ps_nl = m_zero
557 do ik = st%d%kpt%start, st%d%kpt%end
558
559 if (st%kweights(ik) <= m_epsilon) cycle
560
561 do ib = st%group%block_start, st%group%block_end
562 minst = states_elec_block_min(st, ib)
563 maxst = states_elec_block_max(st, ib)
564
565 call hamiltonian_elec_copy_and_set_phase(hm, gr, st%d%kpt, st%group%psib(ib, ik), psib)
566
567 ! calculate the gradient
568 call zderivatives_batch_grad(gr%der, psib, gpsib, set_bc=.false.)
569
570
571 ! Get rV_NL |\psi> for all atoms
572 do idir = 1, gr%der%dim
573 call psib%copy_to(rvnl_psib(idir))
574 call batch_set_zero(rvnl_psib(idir))
575 end do
576 call hm%vnl%zr_vn_local(gr, st%d, gr%der%boundaries%spiral, psib, rvnl_psib)
577
578 do idir = 1, ions%space%periodic_dim
579 do jdir = idir, ions%space%periodic_dim
580 call zmesh_batch_dotp_vector(gr, gpsib(idir), rvnl_psib(jdir), stress_tmp)
581
582 do ist = minst, maxst
583 stress_ps_nl(idir, jdir) = stress_ps_nl(idir, jdir) &
584 + m_two * st%kweights(ik) * st%occ(ist, ik) * real(stress_tmp(ist-minst+1), real64)
585 end do
586
587 end do
588 end do
589
590 do idir = 1, gr%der%dim
591 call rvnl_psib(idir)%end()
592 call gpsib(idir)%end()
593 end do
594 call psib%end()
595 end do
596 end do
597
598 safe_deallocate_a(stress_tmp)
599
600 if (st%parallel_in_states .or. st%d%kpt%parallel) then
601 call comm_allreduce(st%st_kpt_mpi_grp, stress_ps_nl)
602 end if
603
604 ! Symmetrize the kinetic stress tensor
605 call upper_triangular_to_hermitian(ions%space%periodic_dim, stress_ps_nl)
606
607 ! Symmetrize the stress tensor if we use k-point symmetries
608 if (hm%kpoints%use_symmetries) then
609 call dsymmetrize_tensor_cart(gr%symm, stress_ps_nl, use_non_symmorphic=.true.)
610 end if
611
612 ! Add the nonlocal energy
613 do idir = 1, ions%space%periodic_dim
614 stress_ps_nl(idir, idir) = stress_ps_nl(idir, idir) + hm%energy%extern_non_local
615 end do
616
617 stress_ps_nl = stress_ps_nl/ions%latt%rcell_volume
618
619 call profiling_out("STRESS_FROM_PSEUDO_NL")
621
622 end subroutine stress_from_pseudo_nonloc
623
624
625 ! -------------------------------------------------------
642 subroutine stress_from_pseudo_local(gr, st, hm, ions, rho_total, vh, grad_vh, stress_ps_local)
643 type(grid_t), target, intent(in) :: gr
644 type(states_elec_t), intent(inout) :: st
645 type(hamiltonian_elec_t), intent(in) :: hm
646 type(ions_t), intent(in) :: ions
647 real(real64), contiguous, intent(inout) :: rho_total(:)
648 real(real64), intent(in) :: vh(:)
649 real(real64), intent(in) :: grad_vh(:,:)
650 real(real64), intent(out) :: stress_ps_local(3, 3)
651
652
653 real(real64) :: stress_SR(3, 3), stress_LR(3, 3)
654 real(real64) :: energy_ps_SR, charge, zi
655 real(real64), allocatable :: vloc(:), rvloc(:,:), rho_local_lr(:), rho_lr(:)
656 real(real64), allocatable :: grad_rho(:,:), rho_lr_x(:,:), vlr(:), grad_vlr(:,:)
657 integer :: idir, jdir, iatom
658 type(ps_t), pointer :: spec_ps
659
660 call profiling_in("STRESS_FROM_PSEUDO_LOC")
662
663 ! calculate stress from short-range local pseudopotentials
664 stress_sr = m_zero
665
666 safe_allocate(vloc(1:gr%np))
667 vloc = m_zero
668 safe_allocate(rvloc(1:gr%np, 1:gr%der%dim))
669 rvloc = m_zero
670 do iatom = 1, ions%natoms
671 call epot_local_pseudopotential_sr(gr, ions, iatom, vloc, rvloc)
672 end do
673 safe_deallocate_a(vloc)
674
675 safe_allocate(grad_rho(1:gr%np,1:gr%der%dim))
676 call dderivatives_grad(gr%der, rho_total, grad_rho)
677
678 energy_ps_sr = hm%energy%extern_local
679 do idir = 1, ions%space%periodic_dim
680 do jdir = idir, ions%space%periodic_dim
681 stress_sr(idir, jdir) = stress_sr(idir, jdir) &
682 +dmf_dotp(gr, rvloc(:, jdir), grad_rho(:, idir))
683 end do
684 stress_sr(idir,idir) = stress_sr(idir,idir) + energy_ps_sr
685 end do
686
687 call upper_triangular_to_hermitian(ions%space%periodic_dim, stress_sr)
688
689 stress_sr = stress_sr/ions%latt%rcell_volume
690
691 safe_deallocate_a(rvloc)
692 safe_deallocate_a(grad_rho)
693
694
695 ! calculate stress from long-range local pseudopotentials
696 stress_lr = m_zero
697
698 ! We treat the long-range part of the local potential as the Hartree term
699 ! We first sum the long range densities from atoms
700 safe_allocate(rho_lr(1:gr%np_part))
701 safe_allocate(rho_lr_x(1:gr%np, 1:gr%der%dim))
702 rho_lr = m_zero
703 rho_lr_x = m_zero
704 safe_allocate(rho_local_lr(1:gr%np))
705 do iatom = ions%atoms_dist%start, ions%atoms_dist%end
706 assert(ions%atom(iatom)%species%is_ps())
707 call species_get_long_range_density(ions%atom(iatom)%species, ions%namespace, ions%space, ions%latt, &
708 ions%pos(:, iatom), gr, rho_local_lr, nlr_x=rho_lr_x)
709
710 call lalg_axpy(gr%np, m_one, rho_local_lr, rho_lr)
711 end do
712 safe_deallocate_a(rho_local_lr)
713
714 if (ions%atoms_dist%parallel) then
715 call comm_allreduce(ions%atoms_dist%mpi_grp, rho_lr)
716 call comm_allreduce(ions%atoms_dist%mpi_grp, rho_lr_x)
717 end if
718
719 do idir = 1, ions%space%periodic_dim
720 do jdir = idir, ions%space%periodic_dim
721 stress_lr(idir, jdir) = stress_lr(idir, jdir) + dmf_dotp(gr, rho_lr_x(:,jdir), grad_vh(:, idir))
722 end do
723 end do
724 safe_deallocate_a(rho_lr_x)
725
726 safe_allocate(vlr(1:gr%np_part))
727 if (poisson_solver_is_iterative(hm%psolver)) then
728 ! vl has to be initialized before entering routine
729 ! and our best guess for the potential is zero
730 vlr(1:gr%np) = m_zero
731 end if
732 call dpoisson_solve(hm%psolver, ions%namespace, vlr, rho_lr, all_nodes = .true.)
733 safe_deallocate_a(rho_lr)
734
735 safe_allocate(grad_vlr(1:gr%np, 1:gr%der%dim))
736 call dderivatives_grad(gr%der, vlr, grad_vlr)
737 safe_deallocate_a(vlr)
738
739 do idir = 1, ions%space%periodic_dim
740 do jdir = idir, ions%space%periodic_dim
741 stress_lr(idir, jdir) = stress_lr(idir, jdir) - dmf_dotp(gr, grad_vh(:,idir), grad_vlr(:, jdir))/m_two/m_pi
742 end do
743 end do
744
745 call upper_triangular_to_hermitian(ions%space%periodic_dim, stress_lr)
746
747 safe_deallocate_a(grad_vlr)
748
749 ! Contribution from G=0 component of the long-range part
750 !
751 if (ions%space%periodic_dim == 3) then
752 charge = m_zero
753 do iatom = 1, ions%natoms
754 charge = charge + ions%atom(iatom)%species%get_zval()
755 end do
756
757 do iatom = 1, ions%natoms
758 select type(spec => ions%atom(iatom)%species)
759 type is(pseudopotential_t)
760 zi = spec%get_zval()
761 spec_ps => spec%ps
762
763 do idir = 1, ions%space%periodic_dim
764 stress_lr(idir, idir) = stress_lr(idir, idir) &
765 + m_two*m_pi*spec_ps%sigma_erf**2*charge*zi /ions%latt%rcell_volume
766 end do
767 end select
768 end do
769 end if
770
771 stress_lr = stress_lr/ions%latt%rcell_volume
772
773 stress_ps_local = stress_sr + stress_lr
774
775 call profiling_out("STRESS_FROM_PSEUDO_LOC")
777
778 end subroutine stress_from_pseudo_local
779
780 ! -------------------------------------------------------
781 subroutine epot_local_pseudopotential_sr(mesh, ions, iatom, vpsl, rvpsl)
782 class(mesh_t), intent(in) :: mesh
783 type(ions_t), intent(in) :: ions
784 integer, intent(in) :: iatom
785 real(real64), intent(inout) :: vpsl(:)
786 real(real64), intent(inout) :: rvpsl(:,:)
787
788 integer :: ip
789 real(real64) :: radius, vl_ip
790 type(submesh_t) :: sphere
791 type(ps_t), pointer :: ps
792
794
795 if (.not. ions%atom(iatom)%species%is_ps()) then
797 return
798 endif
799
800 call profiling_in("EPOT_LOCAL_PS_SR")
801
802 select type(spec=>ions%atom(iatom)%species)
803 type is(pseudopotential_t)
804
805 ps => spec%ps
806
807 radius = ps%vl%x_threshold*1.05_real64
808
809 call submesh_init(sphere, ions%space, mesh, ions%latt, ions%pos(:, iatom), radius)
810
811 ! Cannot be written (correctly) as a vector expression since for periodic systems,
812 ! there can be values ip, jp such that sphere%map(ip) == sphere%map(jp).
813 do ip = 1, sphere%np
814 vl_ip = spline_eval(ps%vl, sphere%r(ip))
815 vpsl(sphere%map(ip)) = vpsl(sphere%map(ip)) + vl_ip
816 rvpsl(sphere%map(ip), 1:ions%space%periodic_dim) = rvpsl(sphere%map(ip), 1:ions%space%periodic_dim) &
817 + sphere%rel_x(1:ions%space%periodic_dim, ip) * vl_ip
818 end do
819
820 call submesh_end(sphere)
821
822 nullify(ps)
823
824 end select
825
826 call profiling_out("EPOT_LOCAL_PS_SR")
828 end subroutine epot_local_pseudopotential_sr
829
830
831 ! -------------------------------------------------------
845 subroutine stress_from_hubbard(namespace, gr, st, hm, space, rcell_volume, stress_hubbard)
846 type(namespace_t), intent(in) :: namespace
847 type(grid_t), target, intent(in) :: gr
848 type(states_elec_t), intent(inout) :: st
849 type(hamiltonian_elec_t), intent(in) :: hm
850 type(space_t), intent(in) :: space
851 real(real64), intent(in) :: rcell_volume
852 real(real64), intent(out) :: stress_hubbard(3, 3)
853
854 integer :: ik, ist, idir, jdir
855 integer :: ib, minst, maxst
856 type(wfs_elec_t) :: psib, rvu_psib(3), gpsib(3)
857 complex(real64), allocatable :: stress_tmp(:)
858
859 if (hm%lda_u%level == dft_u_none) then
860 stress_hubbard = m_zero
861 return
862 end if
863
864 push_sub_with_profile(stress_from_hubbard)
865
866 assert(st%wfs_type == type_cmplx)
867
868 safe_allocate(stress_tmp(1:st%block_size))
869
870 stress_hubbard = m_zero
871
872 do ik = st%d%kpt%start, st%d%kpt%end
873
874 if (st%kweights(ik) <= m_epsilon) cycle
875
876 do ib = st%group%block_start, st%group%block_end
877 minst = states_elec_block_min(st, ib)
878 maxst = states_elec_block_max(st, ib)
879
880 call hamiltonian_elec_copy_and_set_phase(hm, gr, st%d%kpt, st%group%psib(ib, ik), psib)
881
882 ! calculate the gradient
883 call zderivatives_batch_grad(gr%der, psib, gpsib, set_bc=.false.)
884
885 ! Get rV_U |\psi> for all atoms
886 do idir = 1, gr%der%dim
887 call psib%copy_to(rvu_psib(idir))
888 call batch_set_zero(rvu_psib(idir))
889 end do
890
891 call zlda_u_rvu(hm%lda_u, gr, space, hm%d, namespace, psib, rvu_psib)
892
893 do idir = 1,3
894 do jdir = idir,3
895 call zmesh_batch_dotp_vector(gr, gpsib(idir), rvu_psib(jdir), stress_tmp)
896
897 do ist = minst, maxst
898 stress_hubbard(idir, jdir) = stress_hubbard(idir, jdir) &
899 + m_two * st%kweights(ik) * st%occ(ist, ik) * real(stress_tmp(ist-minst+1), real64)
900 end do
901
902 end do
903 end do
904
905 do idir = 1, gr%der%dim
906 call rvu_psib(idir)%end()
907 call gpsib(idir)%end()
908 end do
909 call psib%end()
910 end do
911 end do
912
913 safe_deallocate_a(stress_tmp)
914
915 if (st%parallel_in_states .or. st%d%kpt%parallel) then
916 call comm_allreduce(st%st_kpt_mpi_grp, stress_hubbard)
917 end if
918
919 ! Symmetrize the kinetic stress tensor
920 call upper_triangular_to_hermitian(gr%der%dim, stress_hubbard)
921
922 ! Symmetrize the stress tensor if we use k-point symmetries
923 if (hm%kpoints%use_symmetries) then
924 call dsymmetrize_tensor_cart(gr%symm, stress_hubbard)
925 end if
926
927 ! Add the Hubbard energy
928 do idir = 1,3
929 stress_hubbard(idir, idir) = stress_hubbard(idir, idir) + hm%energy%int_dft_u
930 end do
931
932 stress_hubbard = stress_hubbard/rcell_volume
933
934 pop_sub_with_profile(stress_from_hubbard)
935 end subroutine stress_from_hubbard
936
937
938 ! -------------------------------------------------------
939 subroutine output_stress(iunit, space_dim, stress_tensors, all_terms)
940 integer, intent(in) :: iunit
941 integer, intent(in) :: space_dim
942 type(stress_t), intent(in) :: stress_tensors
943 logical, optional, intent(in) :: all_terms
944
945 logical :: write_all_terms
946 character(len=16) :: stress_unit
947
948 if (present(all_terms)) then
949 write_all_terms = all_terms
950 else
951 write_all_terms = .true.
952 end if
953
954 write(stress_unit, '(4a,i1)') trim(units_abbrev(units_out%energy)), '/', &
955 trim(units_abbrev(units_out%length)), '^', space_dim
956
957 if (mpi_grp_is_root(mpi_world)) then
958
959 if (write_all_terms) then
960 write(iunit, '(3a)') 'Kinetic stress tensor [', trim(stress_unit), '] ='
961 call print_stress_tensor(iunit, space_dim, stress_tensors%kinetic)
962 if (space_dim == 3) then
963 write(iunit, '(a, es15.6, 3a)') 'Kinetic pressure sumrule violation: ', &
964 units_from_atomic(units_out%energy, stress_tensors%kinetic_sumrule), &
965 ' [', trim(units_abbrev(units_out%energy)), ']'
966 write(iunit,*)
967 end if
968
969
970 write(iunit, '(3a)') 'Hartree stress tensor [', trim(stress_unit), '] ='
971 call print_stress_tensor(iunit, space_dim, stress_tensors%Hartree)
972 if (space_dim == 3) then
973 write(iunit, '(a, es15.6, 3a)') 'Hartree pressure sumrule violation: ', &
974 units_from_atomic(units_out%energy, stress_tensors%hartree_sumrule), &
975 ' [', trim(units_abbrev(units_out%energy)), ']'
976 write(iunit,*)
977 end if
978
979 write(iunit, '(3a)') 'XC stress tensor [', trim(stress_unit), '] ='
980 call print_stress_tensor(iunit, space_dim, stress_tensors%xc)
981
982 write(iunit, '(3a)') 'Local pseudo. stress tensor [', trim(stress_unit), '] ='
983 call print_stress_tensor(iunit, space_dim, stress_tensors%ps_local)
984
985 write(iunit, '(3a)') 'Nonlocal pseudo. stress tensor [', trim(stress_unit), '] ='
986 call print_stress_tensor(iunit, space_dim, stress_tensors%ps_nl)
987
988 write(iunit, '(3a)') 'Ion-ion stress tensor [', trim(stress_unit), '] ='
989 call print_stress_tensor(iunit, space_dim, stress_tensors%ion_ion)
990
991 write(iunit, '(3a)') 'vdW stress tensor [', trim(stress_unit), '] ='
992 call print_stress_tensor(iunit, space_dim, stress_tensors%vdw)
993
994 write(iunit, '(3a)') 'Hubbard stress tensor [', trim(stress_unit), '] ='
995 call print_stress_tensor(iunit, space_dim, stress_tensors%hubbard)
996 end if
997
998 write(iunit, '(3a)') 'Total stress tensor [', trim(stress_unit), '] ='
999 call print_stress_tensor(iunit, space_dim, stress_tensors%total)
1000
1001 end if
1002 end subroutine output_stress
1003
1004
1005 subroutine output_pressure(iunit, space_dim, total_stress_tensor)
1006 integer, intent(in) :: iunit
1007 integer, intent(in) :: space_dim
1008 real(real64), intent(in) :: total_stress_tensor(3,3)
1009
1010 ! TODO(Alex). Issue 884. Move this to unit_system.F90
1011 real(real64), parameter :: au_to_GPa = 29421.02648438959_real64
1012
1013 integer :: idim
1014 real(real64) :: pressure = m_zero
1015 character(len=16) :: stress_unit
1016
1017 write(stress_unit, '(4a,i1)') trim(units_abbrev(units_out%energy)), '/', &
1018 trim(units_abbrev(units_out%length)), '^', space_dim
1019
1020 do idim = 1, space_dim
1021 pressure = pressure - total_stress_tensor(idim, idim) / real(space_dim, real64)
1022 end do
1023
1024 write(iunit,'(3a,es16.8)', advance="no") 'Pressure [', trim(stress_unit), '] = ', &
1025 units_from_atomic(units_out%energy/units_out%length**space_dim, pressure)
1026 if (space_dim == 3) then
1027 write(iunit,'(2x,a,f16.8)') 'Pressure [GPa] = ', pressure * au_to_gpa
1028 else
1029 write(iunit,*)
1030 end if
1031
1032 end subroutine output_pressure
1033
1034 subroutine print_stress_tensor(ounit, space_dim, tensor)
1035 integer, intent(in) :: ounit
1036 integer, intent(in) :: space_dim
1037 real(real64), intent(in) :: tensor(3,3)
1038
1039 real(real64) :: tensor_with_unit(3,3)
1040 integer :: idim, jdim
1041
1042 tensor_with_unit = units_from_atomic(units_out%energy/units_out%length**space_dim, tensor)
1043
1044 write(ounit,'(a9,2x)', advance="no")"T_{ij}"
1045 do jdim = 1, space_dim
1046 write(ounit,'(i18)', advance="no") jdim
1047 end do
1048 write(ounit,*)
1049 do idim = 1, space_dim
1050 write(ounit,'(i9,2x)', advance="no") idim
1051 do jdim = 1, space_dim
1052 write(ounit,'(es18.9)', advance="no") tensor_with_unit(idim, jdim)
1053 end do
1054 write(ounit,*)
1055 end do
1056 write(ounit,*)
1057
1058 end subroutine print_stress_tensor
1059
1060
1061end module stress_oct_m
1062
1063!! Local Variables:
1064!! mode: f90
1065!! coding: utf-8
1066!! End:
constant times a vector plus a vector
Definition: lalg_basic.F90:170
Copies a vector x, to a vector y.
Definition: lalg_basic.F90:185
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:116
subroutine, public batch_set_zero(this, np, async)
fill all mesh functions of the batch with zero
Definition: batch_ops.F90:240
Module implementing boundary conditions in Octopus.
Definition: boundaries.F90:122
This module implements a calculator for the density and defines related functions.
Definition: density.F90:120
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
Definition: epot.F90:166
integer, parameter, public fully_relativistic_zora
Definition: epot.F90:166
real(real64), parameter, public m_two
Definition: global.F90:189
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_four
Definition: global.F90:191
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:185
real(real64), parameter, public m_epsilon
Definition: global.F90:203
real(real64), parameter, public m_one
Definition: global.F90:188
This module implements the underlying real-space grid.
Definition: grid.F90:117
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.
integer, parameter, public independent_particles
integer, parameter, public kohn_sham_dft
This module defines classes and functions for interaction partners.
Definition: io.F90:114
subroutine, public ion_interaction_stress(this, space, latt, atom, natoms, pos, stress_ii)
Computes the contribution to the stress tensor the ion-ion energy.
integer, parameter, public dft_u_none
Definition: lda_u.F90:200
subroutine, public zlda_u_rvu(this, mesh, space, d, namespace, psib, gpsib)
This routine computes .
Definition: lda_u.F90:5352
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
subroutine, public dsymmetrize_matrix(nn, aa)
Definition: math.F90:1445
This module defines functions over batches of mesh functions.
Definition: mesh_batch.F90:116
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.
Definition: mesh.F90:118
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1125
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
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
logical pure function, public poisson_solver_is_iterative(this)
Definition: poisson.F90:1300
subroutine, public dpoisson_solve(this, namespace, pot, rho, all_nodes, kernel)
Calculates the Poisson equation. Given the density returns the corresponding potential.
Definition: poisson.F90:892
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
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)
Definition: splines.F90:441
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.
Definition: stress.F90:118
subroutine stress_from_kinetic(gr, space, hm, st, symm, rcell_volume, stress_kin)
Computes the contribution to the stress tensor from the kinetic energy.
Definition: stress.F90:390
subroutine stress_from_xc(energy, rcell_volume, periodic_dim, stress_xc)
Computes the contribution to the stress tensor from the xc energy.
Definition: stress.F90:528
subroutine print_stress_tensor(ounit, space_dim, tensor)
Definition: stress.F90:1128
subroutine, public output_pressure(iunit, space_dim, total_stress_tensor)
Definition: stress.F90:1099
subroutine epot_local_pseudopotential_sr(mesh, ions, iatom, vpsl, rvpsl)
Definition: stress.F90:875
subroutine, public stress_calculate(namespace, gr, hm, st, ions, ks, ext_partners)
This computes the total stress on the lattice.
Definition: stress.F90:184
subroutine stress_from_hubbard(namespace, gr, st, hm, space, rcell_volume, stress_hubbard)
Computes the contribution to the stress tensor from the Hubbard energy.
Definition: stress.F90:939
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.
Definition: stress.F90:557
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.
Definition: stress.F90:630
subroutine stress_from_hartree(gr, space, volume, vh, grad_vh, ehartree, stress_Hartree)
Computes the contribution to the stress tensor from the Hartree energy.
Definition: stress.F90:479
subroutine, public output_stress(iunit, space_dim, stress_tensors, all_terms)
Definition: stress.F90:1033
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.
Definition: stress.F90:736
subroutine, public submesh_end(this)
Definition: submesh.F90:737
subroutine, public submesh_init(this, space, mesh, latt, center, rc)
Definition: submesh.F90:282
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
Definition: types.F90:134
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:223
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
Definition: xc.F90:114
logical pure function, public xc_is_energy_functional(xcs)
Is one of the x or c functional is not an energy functional.
Definition: xc.F90:728
pure logical function, public in_family(family, xc_families)
Definition: xc.F90:620
A module that takes care of xc contribution from vdW interactions.
Definition: xc_vdw.F90:116
integer(int64), dimension(5), parameter, public d3_lib_options
VDWCORRECTION options that correspond to the DFT-D3 library.
Definition: xc_vdw.F90:169
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:168
Describes mesh distribution to nodes.
Definition: mesh.F90:186
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 ...
Definition: submesh.F90:177
batches of electronic states
Definition: wfs_elec.F90:138
int true(void)