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
47 use lda_u_oct_m
50 use math_oct_m
51 use mesh_oct_m
55 use mpi_oct_m
60 use ps_oct_m
62 use space_oct_m
72 use types_oct_m
73 use unit_oct_m
75 use v_ks_oct_m
77 use xc_oct_m
78 use xc_f03_lib_m
80 implicit none
81
82 private
83 public :: &
87
88contains
89
90 ! ---------------------------------------------------------
92 subroutine stress_calculate(namespace, gr, hm, st, ions, ks, ext_partners)
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
100
101 real(real64), allocatable :: rho_total(:)
102 real(real64) :: stress(3,3) ! stress tensor in Cartesian coordinate
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)
106 integer :: ip
107 real(real64), allocatable :: vh(:)
108 real(real64), allocatable :: grad_vh(:,:)
109 real(real64) :: ehartree
110 real(real64), contiguous, pointer :: rho(:)
111
112 call profiling_in("STRESS_CALCULATE")
113 push_sub(stress_calculate)
114
115 if (st%wfs_type /= type_cmplx) then
116 write(message(1),'(a)') 'The stress tensors for real wavefunctions has not been implemented!'
117
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'
120 call messages_fatal(2, namespace=namespace)
121 end if
122
123 call messages_fatal(1, namespace=namespace)
124 end if
125
126 if (ions%space%periodic_dim == 1) then
127 call messages_not_implemented("Stress tensor for 1D periodic systems", namespace=namespace)
128 end if
129
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'
132 call messages_fatal(1, namespace=namespace)
133 end if
134
135 if (hm%pcm%run_pcm) then
136 call messages_not_implemented('Stress tensor with PCM')
137 end if
138
139 if (allocated(hm%v_static)) then
140 call messages_not_implemented('Stress tensor with static electric fields')
141 end if
142
143 if (ks%has_photons) then
144 call messages_not_implemented('Stress tensor with photon modes')
145 end if
146
147 if (.not. hm%vnl%apply_projector_matrices) then
148 call messages_not_implemented('Stress tensor with relativistic Kleinman-Bylander pseudopotential')
149 end if
150
151 if (hm%ep%reltype == scalar_relativistic_zora .or. hm%ep%reltype == fully_relativistic_zora) then
152 call messages_not_implemented('Stress tensor with ZORA')
153 end if
154
155 ! Checks for the xc part of KS-DFT and GKS-DFT
156 if (ks%theory_level == kohn_sham_dft .or. ks%theory_level == generalized_kohn_sham_dft) then
157 if (.not. xc_is_energy_functional(hm%xc)) then
158 call messages_not_implemented("Stress tensor with xc functionals that are not energy functionals")
159 end if
160
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.'
164 call messages_fatal(2, namespace=namespace)
165 end if
166
167 if (in_family(hm%xc%family, [xc_family_gga]) .and. st%d%ispin == spinors) then
168 call messages_not_implemented("Stress tensor for GGAs with spinors", namespace=namespace)
169 end if
170 end if
171
172 if (hm%magnetic_constrain%level /= constrain_none) then
173 call messages_not_implemented("Stress tensor with MagneticConstrain /= constrain_none")
174 end if
175
176 stress(:,:) = m_zero
177
178 safe_allocate(rho_total(1:gr%np_part))
179 do ip = 1, gr%np
180 rho_total(ip) = sum(st%rho(ip, 1:st%d%nspin))
181 end do
182
183 ! As we rely on some of the full energy components, we need to recompute it first
184 ! TODO: We should restrict the components of the energy needed to be computed
185 call energy_calc_total(namespace, ions%space, hm, gr, st, ext_partners, iunit = -1, full = .true.)
186
187 ! In order to get the electrostatic part (Hartree and local pseudopotential part),
188 ! we need to get the Hartree potential and its gradient
189 safe_allocate(vh(1:gr%np_part))
190 safe_allocate(grad_vh(1:gr%np, 1:gr%der%dim))
191 if (ks%theory_level /= independent_particles) then
192 call lalg_copy(gr%np, hm%ks_pot%vhartree, vh)
193 else
194 if (hm%d%spin_channels > 1) then
195 safe_allocate(rho(1:gr%np_part))
196 call lalg_copy(gr%np, st%rho(:,1), rho)
197 call lalg_axpy(gr%np, m_one, st%rho(:,2), rho)
198 else
199 rho => st%rho(:,1)
200 end if
201 ! In the case of independent particles, we use the electron density without NLCC
202 call dpoisson_solve(hm%psolver, ions%namespace, vh, rho, all_nodes = .true.)
203 if (hm%d%spin_channels > 1) then
204 safe_deallocate_p(rho)
205 else
206 nullify(rho)
207 end if
208 end if
209 ehartree = hm%energy%hartree
210 ! We also compute the gradient here
211 call dderivatives_grad(gr%der, vh, grad_vh)
212
213 ! We now compute the various contributions to the stress tensor
214
215 ! Stress from kinetic energy of electrons
216 call stress_from_kinetic(gr, ions%space, hm, st, gr%symm, ions%latt%rcell_volume, stress_kin)
217 stress = stress + stress_kin
218
219 if (ks%theory_level == independent_particles) then
220 stress_hartree = m_zero
221 stress_xc = m_zero
222 else
223 call stress_from_hartree(gr, ions%space, ions%latt%rcell_volume, vh, grad_vh, ehartree, stress_hartree)
224 stress = stress + stress_hartree
225
226 call stress_from_xc(hm%energy, ions%latt%rcell_volume, ions%space%periodic_dim, stress_xc)
227
228 ! Nonlinear core correction contribution
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
232 end if
233 ! Adds the beyond LDA contribution to the stress tensor
234 stress_xc = stress_xc + ks%stress_xc_gga / ions%latt%rcell_volume
235 stress = stress + stress_xc
236 end if
237
238 call stress_from_pseudo_local(gr, st, hm, ions, rho_total, vh, grad_vh, stress_ps_local)
239 stress_ps = stress_ps_local
240 stress = stress + stress_ps_local
241
242 safe_deallocate_a(vh)
243 safe_deallocate_a(grad_vh)
244
245 call stress_from_pseudo_nonloc(gr, st, hm, ions, stress_ps_nl)
246 stress_ps = stress_ps + stress_ps_nl
247 stress = stress + stress_ps_nl
248
249 call stress_from_hubbard(namespace, gr, st, hm, ions%space, ions%latt%rcell_volume, stress_hubbard)
250 stress = stress + stress_hubbard
251
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
254 ! Stress from kinetic energy of ion
255 ! Stress from ion-field interaction
256
257 ! Sign changed to fit conventional definition
258 stress = -stress
259
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
267
268 ! Stress contribution from vdW D3
269 if (ks%vdw%vdw_correction /= option__vdwcorrection__none) then
270 st%stress_tensors%vdw = hm%ep%vdw_stress
271 else
272 st%stress_tensors%vdw = m_zero
273 end if
274 stress = stress + st%stress_tensors%vdw
275
276 ! Symmetrize the stress tensor if we use k-point symmetries
277 if (hm%kpoints%use_symmetries) then
278 call dsymmetrize_tensor_cart(gr%symm, stress, use_non_symmorphic=.true.)
279 ! We guarantee that the matrix is truely symmetric. There could be small numerical assymetries after symmetrization
280 call dsymmetrize_matrix(ions%space%periodic_dim, stress)
281 end if
282
283 st%stress_tensors%total = stress
284
285 ! Some sumrule for validation
286 ! Sumrule is -3P_{kin}\Omega = 2 E_{kin}
287 st%stress_tensors%kinetic_sumrule = m_zero
288 ! Sumrule is -3P_{Hartree}\Omega = E_{Hartree}
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
293
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
296 end if
297
298 safe_deallocate_a(rho_total)
299
300 pop_sub(stress_calculate)
301 call profiling_out("STRESS_CALCULATE")
302 end subroutine stress_calculate
303
304 ! -------------------------------------------------------
319 subroutine stress_from_kinetic(gr, space, hm, st, symm, rcell_volume, stress_kin)
320 type(grid_t), intent(in) :: gr
321 class(space_t), intent(in) :: space
322 type(hamiltonian_elec_t), intent(in) :: hm
323 type(states_elec_t), intent(inout) :: st
324 type(symmetries_t), intent(in) :: symm
325 real(real64), intent(in) :: rcell_volume
326 real(real64), intent(out) :: stress_kin(3, 3)
327
328 integer :: ik, ist, idir, jdir, ib, minst, maxst
329 complex(real64), allocatable :: stress_l_block(:)
330 type(wfs_elec_t) :: psib, gpsib(space%dim)
331
332 call profiling_in("STRESS_FROM_KINETIC")
333 push_sub(stress_from_kinetic)
334
335 stress_kin(:,:) = m_zero
336
337 safe_allocate(stress_l_block(1:st%block_size))
338
339 do ik = st%d%kpt%start, st%d%kpt%end
340 if (st%kweights(ik) <= m_epsilon) cycle
341
342 do ib = st%group%block_start, st%group%block_end
343 minst = states_elec_block_min(st, ib)
344 maxst = states_elec_block_max(st, ib)
345
346 call hamiltonian_elec_copy_and_set_phase(hm, gr, st%d%kpt, st%group%psib(ib, ik), psib)
347
348 ! calculate the gradient
349 call zderivatives_batch_grad(gr%der, psib, gpsib, set_bc=.false.)
350
351 ! Accumulate the result
352 do idir = 1, space%periodic_dim
353 do jdir = idir, space%periodic_dim
354 call zmesh_batch_dotp_vector(gr, gpsib(idir), gpsib(jdir), stress_l_block)
355
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)
360 end do
361 end do
362 end do
363
364 do idir = 1, space%dim
365 call gpsib(idir)%end()
366 end do
367 call psib%end()
368
369 end do
370 end do
371
372 if (st%parallel_in_states .or. st%d%kpt%parallel) then
373 call comm_allreduce(st%st_kpt_mpi_grp, stress_kin)
374 end if
375
376
377 ! Symmetrize the kinetic stress tensor
378 call upper_triangular_to_hermitian(space%periodic_dim, stress_kin)
379
380 ! Symmetrize the stress tensor if we use k-point symmetries
381 if (hm%kpoints%use_symmetries) then
382 call dsymmetrize_tensor_cart(symm, stress_kin, use_non_symmorphic=.true.)
383 end if
384
385 stress_kin = stress_kin / rcell_volume
386
387 call profiling_out("STRESS_FROM_KINETIC")
388 pop_sub(stress_from_kinetic)
389 end subroutine stress_from_kinetic
390
391 ! -------------------------------------------------------
408 subroutine stress_from_hartree(gr, space, volume, vh, grad_vh, ehartree, stress_Hartree)
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)
416
417 integer :: idir, jdir
418
419 call profiling_in("STRESS_FROM_HARTREE")
420 push_sub(stress_from_hartree)
421
422 stress_hartree(:,:) = m_zero
423
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
427 end do
428 stress_hartree(idir, idir) = stress_hartree(idir, idir) + ehartree
429 end do
430
431 call upper_triangular_to_hermitian(space%periodic_dim, stress_hartree)
432
433 stress_hartree = stress_hartree/volume
434
435 call profiling_out("STRESS_FROM_HARTREE")
436 pop_sub(stress_from_hartree)
437 end subroutine stress_from_hartree
438
439
440 ! -------------------------------------------------------
454 !
455 ! Note: We assume hm%energy%echange, correlation, and intnvxc
456 ! have already been calculated somewhere else.
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)
462
463 integer :: idir
464
465 call profiling_in("STRESS_FROM_XC")
466 push_sub(stress_from_xc)
467
468 stress_xc = m_zero
469 do idir = 1, periodic_dim
470 stress_xc(idir, idir) = - energy%exchange - energy%correlation + energy%intnvxc
471 end do
472 stress_xc(:,:) = stress_xc(:,:) / rcell_volume
473
474 call profiling_out("STRESS_FROM_XC")
475 pop_sub(stress_from_xc)
476 end subroutine stress_from_xc
477
478
479 ! -------------------------------------------------------
486 subroutine stress_from_xc_nlcc(rcell_volume, gr, st, ions, vxc, stress_xc_nlcc)
487 real(real64), intent(in) :: rcell_volume
488 type(grid_t), intent(in) :: gr
489 type(states_elec_t), intent(in) :: st
490 type(ions_t), intent(in) :: ions
491 real(real64), intent(in) :: vxc(:,:)
492 real(real64), intent(out) :: stress_xc_nlcc(3, 3)
493
494 integer :: idir, jdir, iat
495 real(real64), allocatable :: gnlcc(:,:), gnlcc_x(:,:,:), vxc_tot(:)
496
497 call profiling_in("STRESS_FROM_XC_NLCC")
498 push_sub(stress_from_xc_nlcc)
499
500 assert(allocated(st%rho_core))
502 stress_xc_nlcc = m_zero
503
504 ! We first accumulate the contribution from all the pseudo-ions
505 safe_allocate(gnlcc(gr%np, gr%der%dim))
506 safe_allocate(gnlcc_x(gr%np, gr%der%dim, gr%der%dim))
507 gnlcc_x = m_zero
508 do iat = ions%atoms_dist%start, ions%atoms_dist%end
509 assert(ions%atom(iat)%species%is_ps())
510 call species_get_nlcc_grad(ions%atom(iat)%species, ions%space, ions%latt, &
511 ions%pos(:,iat), gr, gnlcc, gnlcc_x)
512 end do
513 safe_deallocate_a(gnlcc)
514
515 if (ions%atoms_dist%parallel) then
516 call comm_allreduce(ions%atoms_dist%mpi_grp, gnlcc_x)
517 end if
518
519 ! Sum over spin of the xc potential
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)
523
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))
527 end do
528 end do
529 safe_deallocate_a(vxc_tot)
530 safe_deallocate_a(gnlcc_x)
531
532 call upper_triangular_to_hermitian(ions%space%periodic_dim, stress_xc_nlcc)
533
534 stress_xc_nlcc(:,:) = stress_xc_nlcc(:,:) / rcell_volume
535
536 call profiling_out("STRESS_FROM_XC_NLCC")
537 pop_sub(stress_from_xc_nlcc)
538 end subroutine stress_from_xc_nlcc
539
540 ! -------------------------------------------------------
559 subroutine stress_from_pseudo_nonloc(gr, st, hm, ions, stress_ps_nl)
560 type(grid_t), target, intent(in) :: gr
561 type(states_elec_t), intent(inout) :: st
562 type(hamiltonian_elec_t), intent(in) :: hm
563 type(ions_t), intent(in) :: ions
564 real(real64), intent(out) :: stress_ps_nl(3, 3)
565
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(:)
570
571 call profiling_in("STRESS_FROM_PSEUDO_NL")
573
574 assert(st%wfs_type == type_cmplx)
575
576 safe_allocate(stress_tmp(1:st%block_size))
577
578 stress_ps_nl = m_zero
580 do ik = st%d%kpt%start, st%d%kpt%end
581
582 if (st%kweights(ik) <= m_epsilon) cycle
583
584 do ib = st%group%block_start, st%group%block_end
585 minst = states_elec_block_min(st, ib)
586 maxst = states_elec_block_max(st, ib)
587
588 call hamiltonian_elec_copy_and_set_phase(hm, gr, st%d%kpt, st%group%psib(ib, ik), psib)
589
590 ! calculate the gradient
591 call zderivatives_batch_grad(gr%der, psib, gpsib, set_bc=.false.)
592
593
594 ! Get rV_NL |\psi> for all atoms
595 do idir = 1, gr%der%dim
596 call psib%copy_to(rvnl_psib(idir))
597 call batch_set_zero(rvnl_psib(idir))
598 end do
599 call hm%vnl%zr_vn_local(gr, st%d, gr%der%boundaries%spiral, psib, rvnl_psib)
600
601 do idir = 1, ions%space%periodic_dim
602 do jdir = idir, ions%space%periodic_dim
603 call zmesh_batch_dotp_vector(gr, gpsib(idir), rvnl_psib(jdir), stress_tmp)
604
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)
608 end do
609
610 end do
611 end do
612
613 do idir = 1, gr%der%dim
614 call rvnl_psib(idir)%end()
615 call gpsib(idir)%end()
616 end do
617 call psib%end()
618 end do
619 end do
620
621 safe_deallocate_a(stress_tmp)
622
623 if (st%parallel_in_states .or. st%d%kpt%parallel) then
624 call comm_allreduce(st%st_kpt_mpi_grp, stress_ps_nl)
625 end if
626
627 ! Symmetrize the kinetic stress tensor
628 call upper_triangular_to_hermitian(ions%space%periodic_dim, stress_ps_nl)
629
630 ! Symmetrize the stress tensor if we use k-point symmetries
631 if (hm%kpoints%use_symmetries) then
632 call dsymmetrize_tensor_cart(gr%symm, stress_ps_nl, use_non_symmorphic=.true.)
633 end if
634
635 ! Add the nonlocal energy
636 do idir = 1, ions%space%periodic_dim
637 stress_ps_nl(idir, idir) = stress_ps_nl(idir, idir) + hm%energy%extern_non_local
638 end do
639
640 stress_ps_nl = stress_ps_nl/ions%latt%rcell_volume
641
642 call profiling_out("STRESS_FROM_PSEUDO_NL")
644
645 end subroutine stress_from_pseudo_nonloc
646
647
648 ! -------------------------------------------------------
665 subroutine stress_from_pseudo_local(gr, st, hm, ions, rho_total, vh, grad_vh, stress_ps_local)
666 type(grid_t), target, intent(in) :: gr
667 type(states_elec_t), intent(inout) :: st
668 type(hamiltonian_elec_t), intent(in) :: hm
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)
674
675
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
682
683 call profiling_in("STRESS_FROM_PSEUDO_LOC")
685
686 ! calculate stress from short-range local pseudopotentials
687 stress_sr = m_zero
688
689 safe_allocate(vloc(1:gr%np))
690 vloc = m_zero
691 safe_allocate(rvloc(1:gr%np, 1:gr%der%dim))
692 rvloc = m_zero
693 do iatom = 1, ions%natoms
694 call epot_local_pseudopotential_sr(gr, ions, iatom, vloc, rvloc)
695 end do
696 safe_deallocate_a(vloc)
697
698 safe_allocate(grad_rho(1:gr%np,1:gr%der%dim))
699 call dderivatives_grad(gr%der, rho_total, grad_rho)
700
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))
706 end do
707 stress_sr(idir,idir) = stress_sr(idir,idir) + energy_ps_sr
708 end do
709
710 call upper_triangular_to_hermitian(ions%space%periodic_dim, stress_sr)
711
712 stress_sr = stress_sr/ions%latt%rcell_volume
713
714 safe_deallocate_a(rvloc)
715 safe_deallocate_a(grad_rho)
716
717
718 ! calculate stress from long-range local pseudopotentials
719 stress_lr = m_zero
720
721 ! We treat the long-range part of the local potential as the Hartree term
722 ! We first sum the long range densities from atoms
723 safe_allocate(rho_lr(1:gr%np_part))
724 safe_allocate(rho_lr_x(1:gr%np, 1:gr%der%dim))
725 rho_lr = m_zero
726 rho_lr_x = m_zero
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())
730 call species_get_long_range_density(ions%atom(iatom)%species, ions%namespace, ions%space, ions%latt, &
731 ions%pos(:, iatom), gr, rho_local_lr, nlr_x=rho_lr_x)
732
733 call lalg_axpy(gr%np, m_one, rho_local_lr, rho_lr)
734 end do
735 safe_deallocate_a(rho_local_lr)
736
737 if (ions%atoms_dist%parallel) then
738 call comm_allreduce(ions%atoms_dist%mpi_grp, rho_lr)
739 call comm_allreduce(ions%atoms_dist%mpi_grp, rho_lr_x)
740 end if
741
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))
745 end do
746 end do
747 safe_deallocate_a(rho_lr_x)
748
749 safe_allocate(vlr(1:gr%np_part))
750 call dpoisson_solve(hm%psolver, ions%namespace, vlr, rho_lr, all_nodes = .true.)
751 safe_deallocate_a(rho_lr)
752
753 safe_allocate(grad_vlr(1:gr%np, 1:gr%der%dim))
754 call dderivatives_grad(gr%der, vlr, grad_vlr)
755 safe_deallocate_a(vlr)
756
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
760 end do
761 end do
762
763 call upper_triangular_to_hermitian(ions%space%periodic_dim, stress_lr)
764
765 safe_deallocate_a(grad_vlr)
766
767 ! Contribution from G=0 component of the long-range part
768 !
769 if (ions%space%periodic_dim == 3) then
770 charge = m_zero
771 do iatom = 1, ions%natoms
772 charge = charge + ions%atom(iatom)%species%get_zval()
773 end do
774
775 do iatom = 1, ions%natoms
776 select type(spec => ions%atom(iatom)%species)
777 type is(pseudopotential_t)
778 zi = spec%get_zval()
779 spec_ps => spec%ps
780
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
784 end do
785 end select
786 end do
787 end if
788
789 stress_lr = stress_lr/ions%latt%rcell_volume
790
791 stress_ps_local = stress_sr + stress_lr
792
793 call profiling_out("STRESS_FROM_PSEUDO_LOC")
795
796 end subroutine stress_from_pseudo_local
797
798 ! -------------------------------------------------------
799 subroutine epot_local_pseudopotential_sr(mesh, ions, iatom, vpsl, rvpsl)
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(:,:)
805
806 integer :: ip
807 real(real64) :: radius, vl_ip
808 type(submesh_t) :: sphere
809 type(ps_t), pointer :: ps
810
812
813 if (.not. ions%atom(iatom)%species%is_ps()) then
815 return
816 endif
817
818 call profiling_in("EPOT_LOCAL_PS_SR")
819
820 select type(spec=>ions%atom(iatom)%species)
821 type is(pseudopotential_t)
822
823 ps => spec%ps
824
825 radius = ps%vl%x_threshold*1.05_real64
826
827 call submesh_init(sphere, ions%space, mesh, ions%latt, ions%pos(:, iatom), radius)
828
829 ! Cannot be written (correctly) as a vector expression since for periodic systems,
830 ! there can be values ip, jp such that sphere%map(ip) == sphere%map(jp).
831 do ip = 1, sphere%np
832 vl_ip = spline_eval(ps%vl, sphere%r(ip))
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
836 end do
837
838 call submesh_end(sphere)
839
840 nullify(ps)
841
842 end select
843
844 call profiling_out("EPOT_LOCAL_PS_SR")
846 end subroutine epot_local_pseudopotential_sr
847
848
849 ! -------------------------------------------------------
863 subroutine stress_from_hubbard(namespace, gr, st, hm, space, rcell_volume, stress_hubbard)
864 type(namespace_t), intent(in) :: namespace
865 type(grid_t), target, intent(in) :: gr
866 type(states_elec_t), intent(inout) :: st
867 type(hamiltonian_elec_t), intent(in) :: hm
868 type(space_t), intent(in) :: space
869 real(real64), intent(in) :: rcell_volume
870 real(real64), intent(out) :: stress_hubbard(3, 3)
871
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(:)
876
877 if (hm%lda_u%level == dft_u_none) then
878 stress_hubbard = m_zero
879 return
880 end if
881
882 push_sub_with_profile(stress_from_hubbard)
883
884 assert(st%wfs_type == type_cmplx)
885
886 safe_allocate(stress_tmp(1:st%block_size))
887
888 stress_hubbard = m_zero
889
890 do ik = st%d%kpt%start, st%d%kpt%end
891
892 if (st%kweights(ik) <= m_epsilon) cycle
893
894 do ib = st%group%block_start, st%group%block_end
895 minst = states_elec_block_min(st, ib)
896 maxst = states_elec_block_max(st, ib)
897
898 call hamiltonian_elec_copy_and_set_phase(hm, gr, st%d%kpt, st%group%psib(ib, ik), psib)
899
900 ! calculate the gradient
901 call zderivatives_batch_grad(gr%der, psib, gpsib, set_bc=.false.)
902
903 ! Get rV_U |\psi> for all atoms
904 do idir = 1, gr%der%dim
905 call psib%copy_to(rvu_psib(idir))
906 call batch_set_zero(rvu_psib(idir))
907 end do
908
909 call zlda_u_rvu(hm%lda_u, gr, space, hm%d, namespace, psib, rvu_psib)
910
911 do idir = 1,3
912 do jdir = idir,3
913 call zmesh_batch_dotp_vector(gr, gpsib(idir), rvu_psib(jdir), stress_tmp)
914
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)
918 end do
919
920 end do
921 end do
922
923 do idir = 1, gr%der%dim
924 call rvu_psib(idir)%end()
925 call gpsib(idir)%end()
926 end do
927 call psib%end()
928 end do
929 end do
930
931 safe_deallocate_a(stress_tmp)
932
933 if (st%parallel_in_states .or. st%d%kpt%parallel) then
934 call comm_allreduce(st%st_kpt_mpi_grp, stress_hubbard)
935 end if
936
937 ! Symmetrize the kinetic stress tensor
938 call upper_triangular_to_hermitian(gr%der%dim, stress_hubbard)
939
940 ! Symmetrize the stress tensor if we use k-point symmetries
941 if (hm%kpoints%use_symmetries) then
942 call dsymmetrize_tensor_cart(gr%symm, stress_hubbard)
943 end if
944
945 ! Add the Hubbard energy
946 do idir = 1,3
947 stress_hubbard(idir, idir) = stress_hubbard(idir, idir) + hm%energy%int_dft_u
948 end do
949
950 stress_hubbard = stress_hubbard/rcell_volume
951
952 pop_sub_with_profile(stress_from_hubbard)
953 end subroutine stress_from_hubbard
954
955
956 ! -------------------------------------------------------
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
962
963 logical :: write_all_terms
964 character(len=16) :: stress_unit
965
966 write_all_terms = optional_default(all_terms, .true.)
967
968 write(stress_unit, '(4a,i1)') trim(units_abbrev(units_out%energy)), '/', &
969 trim(units_abbrev(units_out%length)), '^', space_dim
970
971 if (mpi_grp_is_root(mpi_world)) then
972
973 if (write_all_terms) then
974 write(iunit, '(3a)') 'Kinetic stress tensor [', trim(stress_unit), '] ='
975 call print_stress_tensor(iunit, space_dim, stress_tensors%kinetic)
976 if (space_dim == 3) then
977 write(iunit, '(a, es15.6, 3a)') 'Kinetic pressure sumrule violation: ', &
978 units_from_atomic(units_out%energy, stress_tensors%kinetic_sumrule), &
979 ' [', trim(units_abbrev(units_out%energy)), ']'
980 write(iunit,*)
981 end if
982
983
984 write(iunit, '(3a)') 'Hartree stress tensor [', trim(stress_unit), '] ='
985 call print_stress_tensor(iunit, space_dim, stress_tensors%Hartree)
986 if (space_dim == 3) then
987 write(iunit, '(a, es15.6, 3a)') 'Hartree pressure sumrule violation: ', &
988 units_from_atomic(units_out%energy, stress_tensors%hartree_sumrule), &
989 ' [', trim(units_abbrev(units_out%energy)), ']'
990 write(iunit,*)
991 end if
992
993 write(iunit, '(3a)') 'XC stress tensor [', trim(stress_unit), '] ='
994 call print_stress_tensor(iunit, space_dim, stress_tensors%xc)
995
996 write(iunit, '(3a)') 'Local pseudo. stress tensor [', trim(stress_unit), '] ='
997 call print_stress_tensor(iunit, space_dim, stress_tensors%ps_local)
998
999 write(iunit, '(3a)') 'Nonlocal pseudo. stress tensor [', trim(stress_unit), '] ='
1000 call print_stress_tensor(iunit, space_dim, stress_tensors%ps_nl)
1001
1002 write(iunit, '(3a)') 'Ion-ion stress tensor [', trim(stress_unit), '] ='
1003 call print_stress_tensor(iunit, space_dim, stress_tensors%ion_ion)
1004
1005 write(iunit, '(3a)') 'vdW stress tensor [', trim(stress_unit), '] ='
1006 call print_stress_tensor(iunit, space_dim, stress_tensors%vdw)
1007
1008 write(iunit, '(3a)') 'Hubbard stress tensor [', trim(stress_unit), '] ='
1009 call print_stress_tensor(iunit, space_dim, stress_tensors%hubbard)
1010 end if
1011
1012 write(iunit, '(3a)') 'Total stress tensor [', trim(stress_unit), '] ='
1013 call print_stress_tensor(iunit, space_dim, stress_tensors%total)
1014
1015 end if
1016 end subroutine output_stress
1017
1018
1019 subroutine output_pressure(iunit, space_dim, total_stress_tensor)
1020 integer, intent(in) :: iunit
1021 integer, intent(in) :: space_dim
1022 real(real64), intent(in) :: total_stress_tensor(3,3)
1023
1024 ! TODO(Alex). Issue 884. Move this to unit_system.F90
1025 real(real64), parameter :: au_to_GPa = 29421.02648438959_real64
1026
1027 integer :: idim
1028 real(real64) :: pressure = m_zero
1029 character(len=16) :: stress_unit
1030
1031 write(stress_unit, '(4a,i1)') trim(units_abbrev(units_out%energy)), '/', &
1032 trim(units_abbrev(units_out%length)), '^', space_dim
1033
1034 do idim = 1, space_dim
1035 pressure = pressure - total_stress_tensor(idim, idim) / real(space_dim, real64)
1036 end do
1037
1038 write(iunit,'(3a,es16.8)', advance="no") 'Pressure [', trim(stress_unit), '] = ', &
1039 units_from_atomic(units_out%energy/units_out%length**space_dim, pressure)
1040 if (space_dim == 3) then
1041 write(iunit,'(2x,a,f16.8)') 'Pressure [GPa] = ', pressure * au_to_gpa
1042 else
1043 write(iunit,*)
1044 end if
1045
1046 end subroutine output_pressure
1047
1048 subroutine print_stress_tensor(ounit, space_dim, tensor)
1049 integer, intent(in) :: ounit
1050 integer, intent(in) :: space_dim
1051 real(real64), intent(in) :: tensor(3,3)
1052
1053 real(real64) :: tensor_with_unit(3,3)
1054 integer :: idim, jdim
1055
1056 tensor_with_unit = units_from_atomic(units_out%energy/units_out%length**space_dim, tensor)
1057
1058 write(ounit,'(a9,2x)', advance="no")"T_{ij}"
1059 do jdim = 1, space_dim
1060 write(ounit,'(i18)', advance="no") jdim
1061 end do
1062 write(ounit,*)
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)
1067 end do
1068 write(ounit,*)
1069 end do
1070 write(ounit,*)
1071
1072 end subroutine print_stress_tensor
1073
1074
1075end module stress_oct_m
1076
1077!! Local Variables:
1078!! mode: f90
1079!! coding: utf-8
1080!! End:
constant times a vector plus a vector
Definition: lalg_basic.F90:171
Copies a vector x, to a vector y.
Definition: lalg_basic.F90:186
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:242
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:190
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_four
Definition: global.F90:192
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:186
real(real64), parameter, public m_epsilon
Definition: global.F90:204
real(real64), parameter, public m_one
Definition: global.F90:189
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.
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.
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
Definition: lda_u.F90:201
subroutine, public zlda_u_rvu(this, mesh, space, d, namespace, psib, gpsib)
This routine computes .
Definition: lda_u.F90:5345
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.
Definition: math.F90:115
subroutine, public dsymmetrize_matrix(nn, aa)
Definition: math.F90:1446
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:1113
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:414
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
subroutine, public dpoisson_solve(this, namespace, pot, rho, all_nodes, kernel, reset)
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:413
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:551
subroutine print_stress_tensor(ounit, space_dim, tensor)
Definition: stress.F90:1142
subroutine, public output_pressure(iunit, space_dim, total_stress_tensor)
Definition: stress.F90:1113
subroutine epot_local_pseudopotential_sr(mesh, ions, iatom, vpsl, rvpsl)
Definition: stress.F90:893
subroutine, public stress_calculate(namespace, gr, hm, st, ions, ks, ext_partners)
This computes the total stress on the lattice.
Definition: stress.F90:186
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:957
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:580
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:653
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:502
subroutine, public output_stress(iunit, space_dim, stress_tensors, all_terms)
Definition: stress.F90:1051
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:759
subroutine, public submesh_end(this)
Definition: submesh.F90:735
subroutine, public submesh_init(this, space, mesh, latt, center, rc)
Definition: submesh.F90:280
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
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
A type storing the information and data about a pseudopotential.
Definition: ps.F90:184
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:175
batches of electronic states
Definition: wfs_elec.F90:138
int true(void)