Octopus
sternheimer.F90
Go to the documentation of this file.
1!! Copyright (C) 2004-2012 Xavier Andrade, Eugene S. Kadantsev (ekadants@mjs1.phy.queensu.ca), David Strubbe
2!! Copyright (C) 2021 Davis Welakuh
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
23 use batch_oct_m
25 use debug_oct_m
28 use global_oct_m
29 use grid_oct_m
31 use io_oct_m
32 use, intrinsic :: iso_fortran_env
38 use math_oct_m
39 use mesh_oct_m
42 use mix_oct_m
43 use mpi_oct_m
46 use parser_oct_m
53 use smear_oct_m
54 use space_oct_m
58 use unit_oct_m
60 use types_oct_m
62 use xc_oct_m
63 use xc_f03_lib_m
64
65 implicit none
66
67 private
68 public :: &
91 swap_sigma, &
94 dcalc_hvar, &
95 zcalc_hvar, &
96 dcalc_kvar, &
98
99 character(len=*), public, parameter :: EM_RESP_PHOTONS_DIR = "em_resp_photons/"
100
101 type sternheimer_t
102 private
103 type(linear_solver_t) :: solver
104 type(mix_t) :: mixer
105 type(scf_tol_t) :: scf_tol
106 real(real64) :: lrc_alpha
107 real(real64), allocatable :: fxc(:,:,:)
108 real(real64), allocatable :: kxc(:,:,:,:)
109 real(real64), pointer, contiguous :: drhs(:, :, :, :) => null()
110 complex(real64), pointer, contiguous :: zrhs(:, :, :, :) => null()
111 real(real64), pointer, contiguous :: dinhomog(:, :, :, :, :) => null()
112 complex(real64), pointer, contiguous :: zinhomog(:, :, :, :, :) => null()
113 logical :: add_fxc
114 logical :: add_hartree
115 logical :: ok
116 logical :: occ_response
117 logical :: last_occ_response
118 logical :: occ_response_by_sternheimer
119 logical :: preorthogonalization
120 logical, public :: has_photons
121 real(real64) :: domega
122 complex(real64) :: zomega
123 real(real64), allocatable, public :: dphoton_coord_q(:, :)
124 complex(real64), allocatable, public :: zphoton_coord_q(:, :)
125 real(real64) :: pt_eta
126 type(photon_mode_t) :: pt_modes
127 end type sternheimer_t
128
129
130contains
131
132 !-----------------------------------------------------------
133 subroutine sternheimer_init(this, namespace, space, gr, st, hm, xc, mc, wfs_are_cplx, set_ham_var, set_occ_response, &
134 set_last_occ_response, occ_response_by_sternheimer)
135 type(sternheimer_t), intent(out) :: this
136 type(namespace_t), intent(in) :: namespace
137 class(space_t), intent(in) :: space
138 type(grid_t), intent(inout) :: gr
139 type(states_elec_t), intent(in) :: st
140 type(hamiltonian_elec_t), intent(in) :: hm
141 type(xc_t), intent(in) :: xc
142 type(multicomm_t), intent(in) :: mc
143 logical, intent(in) :: wfs_are_cplx
144 integer, optional, intent(in) :: set_ham_var
145 logical, optional, intent(in) :: set_occ_response
146 logical, optional, intent(in) :: set_last_occ_response
147 logical, optional, intent(in) :: occ_response_by_sternheimer
148
149 integer :: ham_var, iunit
150 logical :: default_preorthog
151
152 push_sub(sternheimer_init)
153
154 if (st%smear%method == smear_fixed_occ) then
155 call messages_experimental("Sternheimer equation for arbitrary occupations", namespace=namespace)
156 end if
157 if (st%smear%method == smear_semiconductor .and. &
158 (abs(st%smear%ef_occ) > m_epsilon) .and. abs(st%smear%ef_occ - m_one) > m_epsilon) then
159 write(message(1),'(a,f12.6)') 'Partial occupation at the Fermi level: ', st%smear%ef_occ
160 message(2) = 'Semiconducting smearing cannot be used for Sternheimer in this situation.'
161 call messages_fatal(2, namespace=namespace)
162 end if
163
164 if (wfs_are_cplx) then
165 call mix_init(this%mixer, namespace, space, gr%der, gr%np, st%d%nspin, func_type_= type_cmplx)
166 else
167 call mix_init(this%mixer, namespace, space, gr%der, gr%np, st%d%nspin, func_type_= type_float)
168 end if
169
170 if (present(set_occ_response)) then
171 this%occ_response = set_occ_response
172 else
173 this%occ_response = .false.
174 end if
175
176 this%occ_response_by_sternheimer = optional_default(occ_response_by_sternheimer, .false.)
177
178 !%Variable Preorthogonalization
179 !%Type logical
180 !%Section Linear Response::Sternheimer
181 !%Description
182 !% Whether initial linear-response wavefunctions should be orthogonalized
183 !% or not against the occupied states, at the start of each SCF cycle.
184 !% Default is true only if <tt>SmearingFunction = semiconducting</tt>,
185 !% or if the <tt>Occupations</tt> block specifies all full or empty states,
186 !% and we are not solving for linear response in the occupied subspace too.
187 !%End
188 default_preorthog = (st%smear%method == smear_semiconductor .or. &
189 (st%smear%method == smear_fixed_occ .and. st%smear%integral_occs)) &
190 .and. .not. this%occ_response
191 call parse_variable(namespace, 'Preorthogonalization', default_preorthog, this%preorthogonalization)
193 !%Variable HamiltonianVariation
194 !%Type integer
195 !%Default hartree+fxc
196 !%Section Linear Response::Sternheimer
197 !%Description
198 !% The terms to be considered in the variation of the
199 !% Hamiltonian. The external potential (V_ext) is always considered. The default is to include
200 !% also the exchange-correlation and Hartree terms, which fully
201 !% takes into account local fields.
202 !% Just <tt>hartree</tt> gives you the random-phase approximation (RPA).
203 !% If you want to choose the exchange-correlation kernel, use the variable
204 !% <tt>XCKernel</tt>. For <tt>kdotp</tt> and magnetic <tt>em_resp</tt> modes,
205 !% or if <tt>TheoryLevel = independent_particles</tt>,
206 !% the value <tt>V_ext_only</tt> is used and this variable is ignored.
207 !%Option V_ext_only 0
208 !% Neither Hartree nor XC potentials included.
209 !%Option hartree 1
210 !% The variation of the Hartree potential only.
211 !%Option fxc 2
212 !% The exchange-correlation kernel (the variation of the
213 !% exchange-correlation potential) only.
214 !%End
216 if (present(set_ham_var)) then
217 ham_var = set_ham_var
218 else if (hm%theory_level /= independent_particles) then
219 call parse_variable(namespace, 'HamiltonianVariation', 3, ham_var)
220 else
221 ham_var = 0
222 end if
223
224 if (hm%theory_level /= independent_particles) then
225 this%add_fxc = ((ham_var / 2) == 1)
226 this%add_hartree = (mod(ham_var, 2) == 1)
227 else
228 this%add_fxc = .false.
229 this%add_hartree = .false.
230 end if
231
232 if (present(set_last_occ_response)) then
233 this%last_occ_response = set_last_occ_response
234 else
235 this%last_occ_response = .false.
236 end if
237
238 message(1) = "Variation of the Hamiltonian in Sternheimer equation: V_ext"
239 if (this%add_hartree) write(message(1), '(2a)') trim(message(1)), ' + hartree'
240 if (this%add_fxc) write(message(1), '(2a)') trim(message(1)), ' + fxc'
241
242 message(2) = "Solving Sternheimer equation for"
243 if (this%occ_response) then
244 write(message(2), '(2a)') trim(message(2)), ' full linear response.'
245 else
246 write(message(2), '(2a)') trim(message(2)), ' linear response in unoccupied subspace only.'
247 end if
248
249 message(3) = "Sternheimer preorthogonalization:"
250 if (this%preorthogonalization) then
251 write(message(3), '(2a)') trim(message(3)), ' yes'
252 else
253 write(message(3), '(2a)') trim(message(3)), ' no'
254 end if
255 call messages_info(3, namespace=namespace)
256
257 call linear_solver_init(this%solver, namespace, gr, states_are_real(st), mc, space)
258
259 ! will not converge for non-self-consistent calculation unless LRTolScheme = fixed
260 if (ham_var == 0) then
261 call scf_tol_init(this%scf_tol, namespace, st%qtot, tol_scheme = 0) ! fixed
262 else
263 call scf_tol_init(this%scf_tol, namespace, st%qtot)
264 end if
265
266 this%lrc_alpha = xc%kernel_lrc_alpha
267
268 if (this%add_fxc) call sternheimer_build_fxc(this, namespace, gr, st, xc)
269
270
271 ! This variable is documented in xc_oep_init.
272 call parse_variable(namespace, 'EnablePhotons', .false., this%has_photons)
273 call messages_print_var_value('EnablePhotons', this%has_photons, namespace=namespace)
274
275 if (this%has_photons) then
276 call messages_experimental('EnablePhotons = yes', namespace=namespace)
277 call photon_mode_init(this%pt_modes, namespace, space%dim)
278 call photon_mode_set_n_electrons(this%pt_modes, m_zero)
279 call photon_mode_compute_dipoles(this%pt_modes, gr)
280 call io_mkdir(em_resp_photons_dir, namespace)
281 iunit = io_open(em_resp_photons_dir // 'photon_modes', namespace, action='write')
282 call photon_mode_write_info(this%pt_modes, iunit=iunit)
283 safe_allocate(this%zphoton_coord_q(1:this%pt_modes%nmodes, 1:space%dim))
284 end if
285
286 !%Variable PhotonEta
287 !%Type float
288 !%Default 0.0000367
289 !%Section Linear Response::Sternheimer
290 !%Description
291 !% This variable provides the value for the broadening of the photonic spectra
292 !% when the coupling of electrons to photons is enabled in the frequency-dependent Sternheimer equation
293 !%End
294 call parse_variable(namespace, 'PhotonEta', 0.0000367_real64, this%pt_eta, units_inp%energy)
295 call messages_print_var_value('PhotonEta', this%pt_eta, units_inp%energy, namespace=namespace)
296
297 if (family_is_mgga(xc%family)) then
298 message(1) = "Using MGGA in Sternheimer calculation."
299 call messages_fatal(1)
300 end if
301
302 pop_sub(sternheimer_init)
303 end subroutine sternheimer_init
304
305 !-----------------------------------------------------------
306 subroutine sternheimer_end(this)
307 type(sternheimer_t), intent(inout) :: this
308
309 push_sub(sternheimer_end)
310
311 safe_deallocate_a(this%zphoton_coord_q)
312
313 call linear_solver_end(this%solver)
314 call scf_tol_end(this%scf_tol)
315 call mix_end(this%mixer)
316
317 safe_deallocate_a(this%fxc)
318
319 pop_sub(sternheimer_end)
320 end subroutine sternheimer_end
321
322
323 !-----------------------------------------------------------
324 subroutine sternheimer_build_fxc(this, namespace, mesh, st, xc)
325 type(sternheimer_t), intent(inout) :: this
326 type(namespace_t), intent(in) :: namespace
327 class(mesh_t), intent(in) :: mesh
328 type(states_elec_t), intent(in) :: st
329 type(xc_t), intent(in) :: xc
330
331 real(real64), allocatable :: rho(:, :)
332
333 push_sub(sternheimer_build_fxc)
334
335 safe_allocate(this%fxc(1:mesh%np, 1:st%d%nspin, 1:st%d%nspin))
336 this%fxc = m_zero
337
338 safe_allocate(rho(1:mesh%np, 1:st%d%nspin))
339 call states_elec_total_density(st, mesh, rho)
340 call xc_get_fxc(xc, mesh, namespace, rho, st%d%ispin, this%fxc)
341 safe_deallocate_a(rho)
342
343 pop_sub(sternheimer_build_fxc)
344
345 end subroutine sternheimer_build_fxc
346
347
348 !-----------------------------------------------------------
349 subroutine sternheimer_build_kxc(this, namespace, mesh, st, xc)
350 type(sternheimer_t), intent(inout) :: this
351 type(namespace_t), intent(in) :: namespace
352 class(mesh_t), intent(in) :: mesh
353 type(states_elec_t), intent(in) :: st
354 type(xc_t), intent(in) :: xc
355
356 real(real64), allocatable :: rho(:, :)
357
358 push_sub(sternheimer_build_kxc)
359
360 if (this%add_fxc) then
361 safe_allocate(this%kxc(1:mesh%np, 1:st%d%nspin, 1:st%d%nspin, 1:st%d%nspin))
362 this%kxc = m_zero
363
364 safe_allocate(rho(1:mesh%np, 1:st%d%nspin))
365 call states_elec_total_density(st, mesh, rho)
366 call xc_get_kxc(xc, mesh, namespace, rho, st%d%ispin, this%kxc)
367 safe_deallocate_a(rho)
368 end if
369
370 pop_sub(sternheimer_build_kxc)
371
372 end subroutine sternheimer_build_kxc
373
374 !---------------------------------------------
375 subroutine sternheimer_unset_kxc(this)
376 type(sternheimer_t), intent(inout) :: this
377
378 push_sub(sternheimer_unset_kxc)
379
380 safe_deallocate_a(this%kxc)
381
382 pop_sub(sternheimer_unset_kxc)
383 end subroutine sternheimer_unset_kxc
384
385 !-----------------------------------------------------------
386 logical function sternheimer_add_fxc(this) result(rr)
387 type(sternheimer_t), intent(in) :: this
388 rr = this%add_fxc
389 end function sternheimer_add_fxc
390
391
392 !-----------------------------------------------------------
393 logical function sternheimer_add_hartree(this) result(rr)
394 type(sternheimer_t), intent(in) :: this
395 rr = this%add_hartree
396 end function sternheimer_add_hartree
397
398
399 !-----------------------------------------------------------
400 logical function sternheimer_has_converged(this) result(rr)
401 type(sternheimer_t), intent(in) :: this
402 rr = this%ok
403 end function sternheimer_has_converged
404
405 !-----------------------------------------------------------
406 logical pure function sternheimer_have_rhs(this) result(have)
407 type(sternheimer_t), intent(in) :: this
408 have = associated(this%drhs) .or. associated(this%zrhs)
409 end function sternheimer_have_rhs
410
411 !-----------------------------------------------------------
412 subroutine sternheimer_unset_rhs(this)
413 type(sternheimer_t), intent(inout) :: this
414
415 push_sub(sternheimer_unset_rhs)
416
417 nullify(this%drhs)
418 nullify(this%zrhs)
419
420 pop_sub(sternheimer_unset_rhs)
421 end subroutine sternheimer_unset_rhs
422
423 !-----------------------------------------------------------
424 logical pure function sternheimer_have_inhomog(this) result(have)
425 type(sternheimer_t), intent(in) :: this
426 have = associated(this%dinhomog) .or. associated(this%zinhomog)
427 end function sternheimer_have_inhomog
428
429 !-----------------------------------------------------------
430 subroutine sternheimer_unset_inhomog(this)
431 type(sternheimer_t), intent(inout) :: this
432
434
435 nullify(this%dinhomog)
436 nullify(this%zinhomog)
437
439 end subroutine sternheimer_unset_inhomog
440
441 !-----------------------------------------------------------
442 integer pure function swap_sigma(sigma)
443 integer, intent(in) :: sigma
444
445 if (sigma == 1) then
446 swap_sigma = 2
447 else
448 swap_sigma = 1
449 end if
450
451 end function swap_sigma
452
453! ---------------------------------------------------------
454 character(len=100) function wfs_tag_sigma(namespace, base_name, isigma) result(str)
455 type(namespace_t), intent(in) :: namespace
456 character(len=*), intent(in) :: base_name
457 integer, intent(in) :: isigma
458
459 character :: sigma_char
460
461 push_sub(wfs_tag_sigma)
462
463 select case (isigma)
464 case (1)
465 sigma_char = '+'
466 case (2)
467 sigma_char = '-'
468 case default
469 write(message(1),'(a,i2)') "Illegal integer isigma passed to wfs_tag_sigma: ", isigma
470 call messages_fatal(1, namespace=namespace)
471 end select
472
473 str = trim(base_name) // sigma_char
474
475 pop_sub(wfs_tag_sigma)
476
477 end function wfs_tag_sigma
478
479 ! --------------------------------------------------------
480
481 subroutine sternheimer_obsolete_variables(namespace, old_prefix, new_prefix)
482 type(namespace_t), intent(in) :: namespace
483 character(len=*), intent(in) :: old_prefix
484 character(len=*), intent(in) :: new_prefix
485
487
488 call messages_obsolete_variable(namespace, trim(old_prefix)//'Preorthogonalization', trim(new_prefix)//'Preorthogonalization')
489 call messages_obsolete_variable(namespace, trim(old_prefix)//'HamiltonianVariation', trim(new_prefix)//'HamiltonianVariation')
490
491 call linear_solver_obsolete_variables(namespace, old_prefix, new_prefix)
492 call scf_tol_obsolete_variables(namespace, old_prefix, new_prefix)
495 end subroutine sternheimer_obsolete_variables
496
497 !--------------------------------------------------------------
498 subroutine calc_hvar_photons(this, mesh, nspin, lr_rho, nsigma, hvar, idir)
499 type(sternheimer_t), intent(inout) :: this
500 class(mesh_t), intent(in) :: mesh
501 integer, intent(in) :: nspin
502 integer, intent(in) :: nsigma
503 complex(real64), intent(in) :: lr_rho(:,:)
504 complex(real64), intent(inout) :: hvar(:,:,:)
505 integer, optional, intent(in) :: idir
506
507 real(real64), allocatable :: lambda_dot_r(:)
508 complex(real64), allocatable :: s_lr_rho(:), vp_dip_self_ener(:), vp_bilinear_el_pt(:)
509 integer :: nm, is, ii
510 complex(real64) :: first_moments
511
512 push_sub(calc_hvar_photons)
513 call profiling_in('CALC_HVAR_PHOTONS')
514
515 nm = this%pt_modes%nmodes
516
517 ! photonic terms
518 safe_allocate(s_lr_rho(1:mesh%np))
519 safe_allocate(lambda_dot_r(1:mesh%np))
520 safe_allocate(vp_dip_self_ener(1:mesh%np))
521 safe_allocate(vp_bilinear_el_pt(1:mesh%np))
522
523 ! spin summed density
524 s_lr_rho = m_zero
525 do is = 1, nspin
526 s_lr_rho = s_lr_rho + lr_rho(:, is)
527 end do
528
529 ! Compute electron-photon potentials
530 vp_dip_self_ener = m_zero
531 vp_bilinear_el_pt = m_zero
532 do ii = 1, nm
533 lambda_dot_r(1:mesh%np) = this%pt_modes%lambda(ii)*this%pt_modes%pol_dipole(1:mesh%np, ii)
534 first_moments = zmf_integrate(mesh, lambda_dot_r(1:mesh%np)*s_lr_rho(1:mesh%np))
536 ! Compute photon displacement coordinate q_{\alpha}s
537 this%zphoton_coord_q(ii, idir) = (m_one/(m_two*(this%pt_modes%omega(ii))**2)) * &
538 ((m_one/(this%zomega - cmplx(this%pt_modes%omega(ii),-this%pt_eta, real64))) - &
539 (m_one/(this%zomega + cmplx(this%pt_modes%omega(ii), this%pt_eta, real64)))) * &
540 (-(this%pt_modes%omega(ii))**2)*first_moments
541
542 ! Compute potential for bilinear el-pt interaction
543 vp_bilinear_el_pt = vp_bilinear_el_pt - &
544 this%pt_modes%omega(ii)*lambda_dot_r(1:mesh%np)*this%zphoton_coord_q(ii, idir)
545
546 ! Compute potential with dipole-self energy term
547 vp_dip_self_ener = vp_dip_self_ener + first_moments*lambda_dot_r(1:mesh%np)
548 end do
549
550 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)
551
552 safe_deallocate_a(s_lr_rho)
553 safe_deallocate_a(lambda_dot_r)
554 safe_deallocate_a(vp_dip_self_ener)
555 safe_deallocate_a(vp_bilinear_el_pt)
556
557 if (nsigma == 2) hvar(1:mesh%np, 1:nspin, 2) = conjg(hvar(1:mesh%np, 1:nspin, 1))
558
559 call profiling_out('CALC_HVAR_PHOTONS')
560 pop_sub(calc_hvar_photons)
561 end subroutine calc_hvar_photons
562
563
564#include "complex.F90"
565#include "sternheimer_inc.F90"
566
567#include "undef.F90"
568
569#include "real.F90"
570#include "sternheimer_inc.F90"
571
572end module sternheimer_oct_m
573
574!! Local Variables:
575!! mode: f90
576!! coding: utf-8
577!! End:
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:180
This module implements batches of mesh functions.
Definition: batch.F90:133
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:116
This module implements a calculator for the density and defines related functions.
Definition: density.F90:120
subroutine, public states_elec_total_density(st, mesh, total_rho)
This routine calculates the total electronic density.
Definition: density.F90:850
real(real64), parameter, public m_zero
Definition: global.F90:188
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
Definition: io.F90:114
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:311
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:352
A module to handle KS potential, without the external potential.
integer, parameter, public independent_particles
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.
Definition: math.F90:115
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
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
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1085
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
subroutine, public mix_init(smix, namespace, space, der, d1, d2, def_, func_type_, prefix_)
Initialise mix_t instance.
Definition: mix.F90:265
subroutine, public mix_end(smix)
Definition: mix.F90:522
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
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)
Definition: scf_tol.F90:160
subroutine, public scf_tol_end(this)
Definition: scf_tol.F90:365
integer, parameter, public smear_semiconductor
Definition: smear.F90:171
integer, parameter, public smear_fixed_occ
Definition: smear.F90:171
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 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, 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 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)
logical function, public 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 zcalc_hvar(namespace, add_hartree, mesh, hm, lrc_alpha, lr_rho, nsigma, hvar, fxc)
subroutine, public sternheimer_unset_kxc(this)
subroutine, public sternheimer_init(this, namespace, space, gr, st, hm, xc, mc, wfs_are_cplx, set_ham_var, set_occ_response, set_last_occ_response, occ_response_by_sternheimer)
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.
subroutine sternheimer_build_fxc(this, namespace, mesh, st, xc)
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 dcalc_hvar(namespace, add_hartree, mesh, hm, lrc_alpha, lr_rho, nsigma, hvar, fxc)
subroutine, public dsternheimer_set_inhomog(this, inhomog)
logical function, public sternheimer_add_hartree(this)
subroutine, public sternheimer_unset_rhs(this)
type(type_t), public type_float
Definition: types.F90:133
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
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
Describes mesh distribution to nodes.
Definition: mesh.F90:186
The states_elec_t class contains all electronic wave functions.