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