34 use,
intrinsic :: iso_fortran_env
104 character(len=*),
public,
parameter :: EM_RESP_PHOTONS_DIR =
"em_resp_photons/"
108 type(linear_solver_t) :: solver
110 type(mixfield_t),
pointer :: mixfield
111 type(scf_tol_t) :: scf_tol
112 real(real64),
allocatable :: fxc(:,:,:)
113 real(real64),
allocatable :: kxc(:,:,:,:)
114 real(real64),
pointer,
contiguous :: drhs(:, :, :, :) => null()
115 complex(real64),
pointer,
contiguous :: zrhs(:, :, :, :) => null()
116 real(real64),
pointer,
contiguous :: dinhomog(:, :, :, :, :) => null()
117 complex(real64),
pointer,
contiguous :: zinhomog(:, :, :, :, :) => null()
118 logical :: add_fxc_kernel
120 logical :: occ_response
121 logical :: last_occ_response
122 logical :: occ_response_by_sternheimer
123 logical :: preorthogonalization
124 logical,
public :: has_photons
125 real(real64) :: domega
126 complex(real64) :: zomega
127 real(real64),
allocatable,
public :: dphoton_coord_q(:, :)
128 complex(real64),
allocatable,
public :: zphoton_coord_q(:, :)
129 real(real64) :: pt_eta
130 type(photon_mode_t) :: pt_modes
132 real(real64) :: coeff_hartree
142 subroutine sternheimer_init(this, namespace, space, gr, st, hm, ks, mc, wfs_are_cplx, set_ham_var, set_occ_response, &
143 set_last_occ_response, occ_response_by_sternheimer)
144 type(sternheimer_t),
intent(out) :: this
145 type(namespace_t),
intent(in) :: namespace
146 class(space_t),
intent(in) :: space
147 type(grid_t),
intent(inout) :: gr
148 type(states_elec_t),
intent(in) :: st
149 type(hamiltonian_elec_t),
intent(in) :: hm
150 type(v_ks_t),
intent(in) :: ks
151 type(multicomm_t),
intent(in) :: mc
152 logical,
intent(in) :: wfs_are_cplx
153 integer,
optional,
intent(in) :: set_ham_var
154 logical,
optional,
intent(in) :: set_occ_response
155 logical,
optional,
intent(in) :: set_last_occ_response
156 logical,
optional,
intent(in) :: occ_response_by_sternheimer
158 logical :: add_hartree
159 integer :: ham_var, iunit
160 logical :: default_preorthog
184 write(
message(1),
'(a,f12.6)')
'Partial occupation at the Fermi level: ', st%smear%ef_occ
185 message(2) =
'Semiconducting smearing cannot be used for Sternheimer in this situation.'
189 if (wfs_are_cplx)
then
190 call mix_init(this%mixer, namespace, space, gr%der, gr%np, st%d%nspin, func_type_=
type_cmplx)
192 call mix_init(this%mixer, namespace, space, gr%der, gr%np, st%d%nspin, func_type_=
type_float)
197 this%occ_response_by_sternheimer =
optional_default(occ_response_by_sternheimer, .false.)
212 .and. .not. this%occ_response
213 call parse_variable(namespace,
'Preorthogonalization', default_preorthog, this%preorthogonalization)
237 if (
present(set_ham_var))
then
238 ham_var = set_ham_var
240 call parse_variable(namespace,
'HamiltonianVariation', 3, ham_var)
246 this%add_fxc_kernel = ((ham_var / 2) == 1)
247 add_hartree = (mod(ham_var, 2) == 1)
249 this%add_fxc_kernel = .false.
250 add_hartree = .false.
254 message(1) =
"Variation of the Hamiltonian in Sternheimer equation: V_ext"
255 if (add_hartree)
write(
message(1),
'(2a)') trim(
message(1)),
' + hartree'
256 if (this%add_fxc_kernel)
write(
message(1),
'(2a)') trim(
message(1)),
' + fxc'
258 message(2) =
"Solving Sternheimer equation for"
259 if (this%occ_response)
then
260 write(
message(2),
'(2a)') trim(
message(2)),
' full linear response.'
262 write(
message(2),
'(2a)') trim(
message(2)),
' linear response in unoccupied subspace only.'
265 message(3) =
"Sternheimer preorthogonalization:"
266 if (this%preorthogonalization)
then
276 if (ham_var == 0)
then
277 call scf_tol_init(this%scf_tol, namespace, st%qtot, tol_scheme = 0)
283 this%coeff_hartree =
m_zero
284 if (this%add_fxc_kernel)
then
285 this%coeff_hartree = this%coeff_hartree - ks%xc%lrc%alpha / (
m_four *
m_pi)
287 if (add_hartree) this%coeff_hartree = this%coeff_hartree +
m_one
288 if (this%add_fxc_kernel)
then
299 call parse_variable(namespace,
'EnablePhotons', .false., this%has_photons)
302 if (this%has_photons)
then
307 call io_mkdir(em_resp_photons_dir, namespace)
308 iunit =
io_open(em_resp_photons_dir //
'photon_modes', namespace, action=
'write')
310 safe_allocate(this%zphoton_coord_q(1:this%pt_modes%nmodes, 1:space%dim))
334 safe_deallocate_a(this%zphoton_coord_q)
340 nullify(this%mixfield)
342 safe_deallocate_a(this%fxc)
353 type(
grid_t),
intent(in) :: gr
355 type(
xc_t),
intent(in) :: xc
357 real(real64),
allocatable :: rho(:, :)
361 safe_allocate(this%fxc(1:gr%np, 1:st%d%nspin, 1:st%d%nspin))
362 safe_allocate(rho(1:gr%np_part, 1:st%d%nspin))
364 call xc_get_fxc(xc, gr, namespace, rho, st%d%ispin, this%fxc)
365 safe_deallocate_a(rho)
375 class(
mesh_t),
intent(in) :: mesh
377 type(
xc_t),
intent(in) :: xc
379 real(real64),
allocatable :: rho(:, :)
383 if (this%add_fxc())
then
384 safe_allocate(this%kxc(1:mesh%np, 1:st%d%nspin, 1:st%d%nspin, 1:st%d%nspin))
387 safe_allocate(rho(1:mesh%np, 1:st%d%nspin))
389 call xc_get_kxc(xc, mesh, namespace, rho, st%d%ispin, this%kxc)
390 safe_deallocate_a(rho)
403 safe_deallocate_a(this%kxc)
411 rr = this%add_fxc_kernel
431 have =
associated(this%drhs) .or.
associated(this%zrhs)
447 logical pure function sternheimer_have_inhomog(this) result(have)
449 have =
associated(this%dinhomog) .or.
associated(this%zinhomog)
458 nullify(this%dinhomog)
459 nullify(this%zinhomog)
465 integer pure function swap_sigma(sigma)
466 integer,
intent(in) :: sigma
477 character(len=100) function wfs_tag_sigma(namespace, base_name, isigma)
result(str)
478 type(namespace_t),
intent(in) :: namespace
479 character(len=*),
intent(in) :: base_name
480 integer,
intent(in) :: isigma
482 character :: sigma_char
492 write(message(1),
'(a,i2)')
"Illegal integer isigma passed to wfs_tag_sigma: ", isigma
493 call messages_fatal(1, namespace=namespace)
496 str = trim(base_name) // sigma_char
505 type(namespace_t),
intent(in) :: namespace
506 character(len=*),
intent(in) :: old_prefix
507 character(len=*),
intent(in) :: new_prefix
511 call messages_obsolete_variable(namespace, trim(old_prefix)//
'Preorthogonalization', trim(new_prefix)//
'Preorthogonalization')
512 call messages_obsolete_variable(namespace, trim(old_prefix)//
'HamiltonianVariation', trim(new_prefix)//
'HamiltonianVariation')
514 call linear_solver_obsolete_variables(namespace, old_prefix, new_prefix)
515 call scf_tol_obsolete_variables(namespace, old_prefix, new_prefix)
523 class(mesh_t),
intent(in) :: mesh
524 integer,
intent(in) :: nspin
525 integer,
intent(in) :: nsigma
526 complex(real64),
intent(in) :: lr_rho(:,:)
527 complex(real64),
intent(inout) :: hvar(:,:,:)
528 integer,
optional,
intent(in) :: idir
530 real(real64),
allocatable :: lambda_dot_r(:)
531 complex(real64),
allocatable :: s_lr_rho(:), vp_dip_self_ener(:), vp_bilinear_el_pt(:)
532 integer :: nm, is, ii
533 complex(real64) :: first_moments
536 call profiling_in(
'CALC_HVAR_PHOTONS')
538 nm = this%pt_modes%nmodes
541 safe_allocate(s_lr_rho(1:mesh%np))
542 safe_allocate(lambda_dot_r(1:mesh%np))
543 safe_allocate(vp_dip_self_ener(1:mesh%np))
544 safe_allocate(vp_bilinear_el_pt(1:mesh%np))
549 s_lr_rho = s_lr_rho + lr_rho(:, is)
553 vp_dip_self_ener = m_zero
554 vp_bilinear_el_pt = m_zero
556 lambda_dot_r(1:mesh%np) = this%pt_modes%lambda(ii)*this%pt_modes%pol_dipole(1:mesh%np, ii)
557 first_moments = zmf_integrate(mesh, lambda_dot_r(1:mesh%np)*s_lr_rho(1:mesh%np))
560 this%zphoton_coord_q(ii, idir) = (m_one/(m_two*(this%pt_modes%omega(ii))**2)) * &
561 ((m_one/(this%zomega - cmplx(this%pt_modes%omega(ii),-this%pt_eta, real64))) - &
562 (m_one/(this%zomega + cmplx(this%pt_modes%omega(ii), this%pt_eta, real64)))) * &
563 (-(this%pt_modes%omega(ii))**2)*first_moments
566 vp_bilinear_el_pt = vp_bilinear_el_pt - &
567 this%pt_modes%omega(ii)*lambda_dot_r(1:mesh%np)*this%zphoton_coord_q(ii, idir)
570 vp_dip_self_ener = vp_dip_self_ener + first_moments*lambda_dot_r(1:mesh%np)
573 hvar(1:mesh%np, 1, 1) = hvar(1:mesh%np, 1, 1) + vp_dip_self_ener(1:mesh%np) + vp_bilinear_el_pt(1:mesh%np)
575 safe_deallocate_a(s_lr_rho)
576 safe_deallocate_a(lambda_dot_r)
577 safe_deallocate_a(vp_dip_self_ener)
578 safe_deallocate_a(vp_bilinear_el_pt)
580 if (nsigma == 2) hvar(1:mesh%np, 1:nspin, 2) = conjg(hvar(1:mesh%np, 1:nspin, 1))
582 call profiling_out(
'CALC_HVAR_PHOTONS')
587#include "complex.F90"
588#include "sternheimer_inc.F90"
593#include "sternheimer_inc.F90"
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
This module implements batches of mesh functions.
This module implements common operations on batches of mesh functions.
This module implements a calculator for the density and defines related functions.
subroutine, public states_elec_total_density(st, mesh, total_rho)
This routine calculates the total electronic density.
real(real64), parameter, public m_zero
real(real64), parameter, public m_four
real(real64), parameter, public m_pi
some mathematical constants
integer, parameter, public independent_particles
Theory level.
real(real64), parameter, public m_epsilon
real(real64), parameter, public m_one
This module implements the underlying real-space grid.
subroutine, public io_mkdir(fname, namespace, parents)
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
A module to handle KS potential, without the external potential.
integer, parameter, public dft_u_none
subroutine, public linear_solver_end(this)
subroutine, public linear_solver_init(this, namespace, gr, states_are_real, mc, space)
This module is intended to contain "only mathematical" functions and procedures.
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
subroutine, public messages_not_implemented(feature, namespace)
character(len=512), private msg
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_experimental(name, namespace)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
subroutine, public mix_get_field(this, mixfield)
subroutine, public mix_init(smix, namespace, space, der, d1, d2, def_, func_type_, prefix_)
Initialise mix_t instance.
subroutine, public mix_end(smix)
This module handles the communicators for the various parallelization strategies.
subroutine, public photon_mode_compute_dipoles(this, mesh)
Computes the polarization dipole.
subroutine, public photon_mode_write_info(this, iunit, namespace)
subroutine, public photon_mode_set_n_electrons(this, qtot)
subroutine, public photon_mode_init(this, namespace, dim, photon_free)
subroutine, public scf_tol_init(this, namespace, qtot, def_maximumiter, tol_scheme)
subroutine, public scf_tol_end(this)
integer, parameter, public smear_semiconductor
integer, parameter, public smear_fixed_occ
pure logical function, public states_are_real(st)
This module handles reading and writing restart information for the states_elec_t.
subroutine, public zsternheimer_set_inhomog(this, inhomog)
subroutine, public zsternheimer_calc_hvar(this, namespace, mesh, hm, lr, nsigma, hvar, idir)
subroutine, public zcalc_hvar(namespace, coeff_hartree, mesh, hm, lr_rho, nsigma, hvar, fxc)
Computes the first-order variation of the Kohn-Sham Hamiltonian.
subroutine, public dsternheimer_solve(this, namespace, space, gr, kpoints, st, hm, mc, lr, nsigma, omega, perturbation, restart, rho_tag, wfs_tag, idir, have_restart_rho, have_exact_freq)
This routine calculates the first-order variations of the wavefunctions for an applied perturbation.
subroutine, public dsternheimer_calc_hvar(this, namespace, mesh, hm, lr, nsigma, hvar, idir)
integer pure function, public swap_sigma(sigma)
subroutine sternheimer_build_fxc(this, namespace, gr, st, xc)
Builds the exchange-correlation kernel for computing the density response.
subroutine, public dcalc_kvar(this, mesh, st, lr_rho1, lr_rho2, nsigma, kvar)
logical function, public sternheimer_has_converged(this)
subroutine, public zsternheimer_set_rhs(this, rhs)
character(len=100) function, public wfs_tag_sigma(namespace, base_name, isigma)
subroutine, public dcalc_hvar(namespace, coeff_hartree, mesh, hm, lr_rho, nsigma, hvar, fxc)
Computes the first-order variation of the Kohn-Sham Hamiltonian.
subroutine calc_hvar_photons(this, mesh, nspin, lr_rho, nsigma, hvar, idir)
logical pure function, public sternheimer_have_rhs(this)
subroutine, public dsternheimer_solve_order2(sh1, sh2, sh_2ndorder, namespace, space, gr, kpoints, st, hm, mc, lr1, lr2, nsigma, omega1, omega2, pert1, pert2, lr_2ndorder, pert_2ndorder, restart, rho_tag, wfs_tag, have_restart_rho, have_exact_freq, give_pert1psi2, give_dl_eig1)
subroutine, public dsternheimer_set_rhs(this, rhs)
pure logical function sternheimer_add_fxc(this)
subroutine, public zsternheimer_solve_order2(sh1, sh2, sh_2ndorder, namespace, space, gr, kpoints, st, hm, mc, lr1, lr2, nsigma, omega1, omega2, pert1, pert2, lr_2ndorder, pert_2ndorder, restart, rho_tag, wfs_tag, have_restart_rho, have_exact_freq, give_pert1psi2, give_dl_eig1)
subroutine, public sternheimer_unset_kxc(this)
subroutine, public zsternheimer_solve(this, namespace, space, gr, kpoints, st, hm, mc, lr, nsigma, omega, perturbation, restart, rho_tag, wfs_tag, idir, have_restart_rho, have_exact_freq)
This routine calculates the first-order variations of the wavefunctions for an applied perturbation.
pure logical function sternheimer_add_hartree(this)
subroutine, public sternheimer_unset_inhomog(this)
logical pure function, public sternheimer_have_inhomog(this)
subroutine, public zcalc_kvar(this, mesh, st, lr_rho1, lr_rho2, nsigma, kvar)
subroutine, public sternheimer_end(this)
subroutine, public sternheimer_obsolete_variables(namespace, old_prefix, new_prefix)
subroutine, public sternheimer_build_kxc(this, namespace, mesh, st, xc)
subroutine, public dsternheimer_set_inhomog(this, inhomog)
subroutine, public sternheimer_init(this, namespace, space, gr, st, hm, ks, mc, wfs_are_cplx, set_ham_var, set_occ_response, set_last_occ_response, occ_response_by_sternheimer)
subroutine, public sternheimer_unset_rhs(this)
type(type_t), public type_float
type(type_t), public type_cmplx
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
This module defines the unit system, used for input and output.
type(unit_system_t), public units_inp
the units systems for reading and writing
This module provices a simple timer class which can be used to trigger the writing of a restart file ...
subroutine, public xc_get_fxc(xcs, gr, namespace, rho, ispin, fxc, fxc_grad, fxc_grad_spin)
Returns the exchange-correlation kernel.
subroutine, public xc_get_kxc(xcs, mesh, namespace, rho, ispin, kxc)
pure logical function, public family_is_mgga(family, only_collinear)
Is the xc function part of the mGGA family.
logical pure function, public family_is_hybrid(xcs)
Returns true if the functional is an hybrid functional.
subroutine, public xc_write_fxc_info(xcs, iunit, namespace)
integer, parameter, public sic_none
no self-interaction correction
subroutine, public xc_sic_write_info(sic, iunit, namespace)
Description of the grid, containing information on derivatives, stencil, and symmetries.
Describes mesh distribution to nodes.
The states_elec_t class contains all electronic wave functions.