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
38 use lda_u_oct_m
41 use math_oct_m
42 use mesh_oct_m
45 use mix_oct_m
46 use mpi_oct_m
49 use parser_oct_m
56 use smear_oct_m
57 use space_oct_m
61 use unit_oct_m
63 use types_oct_m
64 use v_ks_oct_m
66 use xc_oct_m
67 use xc_sic_oct_m
68 use xc_f03_lib_m
71
72 implicit none
73
74 private
75 public :: &
96 swap_sigma, &
99 dcalc_hvar, &
100 zcalc_hvar, &
101 dcalc_kvar, &
103
104 character(len=*), public, parameter :: EM_RESP_PHOTONS_DIR = "em_resp_photons/"
105
106 type sternheimer_t
107 private
108 type(linear_solver_t) :: solver
109 type(mix_t) :: mixer
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
119 logical :: ok
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
131
132 real(real64) :: coeff_hartree
133 contains
134 procedure :: add_fxc => sternheimer_add_fxc
135 procedure :: add_hartree => sternheimer_add_hartree
136 end type sternheimer_t
137
138
139contains
140
141 !-----------------------------------------------------------
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
157
158 logical :: add_hartree
159 integer :: ham_var, iunit
160 logical :: default_preorthog
161
162 push_sub(sternheimer_init)
163
164 if (family_is_mgga(ks%xc%family)) then
165 call messages_not_implemented("MGGA functionals with Sternheimer calculation")
166 end if
167 if (family_is_hybrid(ks%xc)) then
168 call messages_not_implemented("Hybrid functionals with Sternheimer calculation")
169 end if
170
171 if (ks%sic%level /= sic_none) then
172 call messages_not_implemented("SIC with Sternheimer calculation")
173 end if
174
175 if (hm%lda_u_level /= dft_u_none) then
176 call messages_not_implemented("DFT+U with Sterheimer calculation")
177 end if
178
179 if (st%smear%method == smear_fixed_occ) then
180 call messages_experimental("Sternheimer equation for arbitrary occupations", namespace=namespace)
181 end if
182 if (st%smear%method == smear_semiconductor .and. &
183 (abs(st%smear%ef_occ) > m_epsilon) .and. abs(st%smear%ef_occ - m_one) > m_epsilon) then
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.'
186 call messages_fatal(2, namespace=namespace)
187 end if
188
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)
191 else
192 call mix_init(this%mixer, namespace, space, gr%der, gr%np, st%d%nspin, func_type_= type_float)
193 end if
194 call mix_get_field(this%mixer, this%mixfield)
195
196 this%occ_response = optional_default(set_occ_response, .false.)
197 this%occ_response_by_sternheimer = optional_default(occ_response_by_sternheimer, .false.)
198 this%last_occ_response = optional_default(set_last_occ_response, .false.)
200 !%Variable Preorthogonalization
201 !%Type logical
202 !%Section Linear Response::Sternheimer
203 !%Description
204 !% Whether initial linear-response wavefunctions should be orthogonalized
205 !% or not against the occupied states, at the start of each SCF cycle.
206 !% Default is true only if <tt>SmearingFunction = semiconducting</tt>,
207 !% or if the <tt>Occupations</tt> block specifies all full or empty states,
208 !% and we are not solving for linear response in the occupied subspace too.
209 !%End
210 default_preorthog = (st%smear%method == smear_semiconductor .or. &
211 (st%smear%method == smear_fixed_occ .and. st%smear%integral_occs)) &
212 .and. .not. this%occ_response
213 call parse_variable(namespace, 'Preorthogonalization', default_preorthog, this%preorthogonalization)
215 !%Variable HamiltonianVariation
216 !%Type integer
217 !%Default hartree+fxc
218 !%Section Linear Response::Sternheimer
219 !%Description
220 !% The terms to be considered in the variation of the
221 !% Hamiltonian. The external potential (V_ext) is always considered. The default is to include
222 !% also the exchange-correlation and Hartree terms, which fully
223 !% takes into account local fields.
224 !% Just <tt>hartree</tt> gives you the random-phase approximation (RPA).
225 !% If you want to choose the exchange-correlation kernel, use the variable
226 !% <tt>XCKernel</tt>. For <tt>kdotp</tt> and magnetic <tt>em_resp</tt> modes,
227 !% or if <tt>TheoryLevel = independent_particles</tt>,
228 !% the value <tt>V_ext_only</tt> is used and this variable is ignored.
229 !%Option V_ext_only 0
230 !% Neither Hartree nor XC potentials included.
231 !%Option hartree 1
232 !% The variation of the Hartree potential only.
233 !%Option fxc 2
234 !% The exchange-correlation kernel (the variation of the
235 !% exchange-correlation potential) only.
236 !%End
237 if (present(set_ham_var)) then
238 ham_var = set_ham_var
239 else if (hm%theory_level /= independent_particles) then
240 call parse_variable(namespace, 'HamiltonianVariation', 3, ham_var)
241 else
242 ham_var = 0
243 end if
244
245 if (hm%theory_level /= independent_particles) then
246 this%add_fxc_kernel = ((ham_var / 2) == 1)
247 add_hartree = (mod(ham_var, 2) == 1)
248 else
249 this%add_fxc_kernel = .false.
250 add_hartree = .false.
251 end if
252
253
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'
257
258 message(2) = "Solving Sternheimer equation for"
259 if (this%occ_response) then
260 write(message(2), '(2a)') trim(message(2)), ' full linear response.'
261 else
262 write(message(2), '(2a)') trim(message(2)), ' linear response in unoccupied subspace only.'
263 end if
264
265 message(3) = "Sternheimer preorthogonalization:"
266 if (this%preorthogonalization) then
267 write(message(3), '(2a)') trim(message(3)), ' yes'
268 else
269 write(message(3), '(2a)') trim(message(3)), ' no'
270 end if
271 call messages_info(3, namespace=namespace)
272
273 call linear_solver_init(this%solver, namespace, gr, states_are_real(st), mc, space)
274
275 ! will not converge for non-self-consistent calculation unless LRTolScheme = fixed
276 if (ham_var == 0) then
277 call scf_tol_init(this%scf_tol, namespace, st%qtot, tol_scheme = 0) ! fixed
278 else
279 call scf_tol_init(this%scf_tol, namespace, st%qtot)
280 end if
281
282 ! Build the fxc from the GS density
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)
286 end if
287 if (add_hartree) this%coeff_hartree = this%coeff_hartree + m_one
288 if (this%add_fxc_kernel) then
289 call sternheimer_build_fxc(this, namespace, gr, st, ks%xc)
290
291 call messages_print_with_emphasis(msg="XC Kernel level", namespace=namespace)
292 call xc_write_fxc_info(ks%xc, namespace=namespace)
293 call xc_sic_write_info(ks%sic, namespace=namespace)
294 call messages_print_with_emphasis(namespace=namespace)
295 end if
296
297
298 ! This variable is documented in xc_oep_init.
299 call parse_variable(namespace, 'EnablePhotons', .false., this%has_photons)
300 call messages_print_var_value('EnablePhotons', this%has_photons, namespace=namespace)
301
302 if (this%has_photons) then
303 call messages_experimental('EnablePhotons = yes', namespace=namespace)
304 call photon_mode_init(this%pt_modes, namespace, space%dim)
305 call photon_mode_set_n_electrons(this%pt_modes, m_zero)
306 call photon_mode_compute_dipoles(this%pt_modes, gr)
307 call io_mkdir(em_resp_photons_dir, namespace)
308 iunit = io_open(em_resp_photons_dir // 'photon_modes', namespace, action='write')
309 call photon_mode_write_info(this%pt_modes, iunit=iunit)
310 safe_allocate(this%zphoton_coord_q(1:this%pt_modes%nmodes, 1:space%dim))
311 end if
312
313 !%Variable PhotonEta
314 !%Type float
315 !%Default 0.0000367
316 !%Section Linear Response::Sternheimer
317 !%Description
318 !% This variable provides the value for the broadening of the photonic spectra
319 !% when the coupling of electrons to photons is enabled in the frequency-dependent Sternheimer equation
320 !%End
321 call parse_variable(namespace, 'PhotonEta', 0.0000367_real64, this%pt_eta, units_inp%energy)
322 call messages_print_var_value('PhotonEta', this%pt_eta, units_inp%energy, namespace=namespace)
323
324
325 pop_sub(sternheimer_init)
326 end subroutine sternheimer_init
327
328 !-----------------------------------------------------------
329 subroutine sternheimer_end(this)
330 type(sternheimer_t), intent(inout) :: this
331
332 push_sub(sternheimer_end)
333
334 safe_deallocate_a(this%zphoton_coord_q)
335
336 call linear_solver_end(this%solver)
337 call scf_tol_end(this%scf_tol)
338 call mix_end(this%mixer)
339
340 nullify(this%mixfield)
341
342 safe_deallocate_a(this%fxc)
343
344 pop_sub(sternheimer_end)
345 end subroutine sternheimer_end
346
347
348 !-----------------------------------------------------------
350 subroutine sternheimer_build_fxc(this, namespace, gr, st, xc)
351 type(sternheimer_t), intent(inout) :: this
352 type(namespace_t), intent(in) :: namespace
353 type(grid_t), intent(in) :: gr
354 type(states_elec_t), intent(in) :: st
355 type(xc_t), intent(in) :: xc
356
357 real(real64), allocatable :: rho(:, :)
358
359 push_sub(sternheimer_build_fxc)
360
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))
363 call states_elec_total_density(st, gr, rho)
364 call xc_get_fxc(xc, gr, namespace, rho, st%d%ispin, this%fxc)
365 safe_deallocate_a(rho)
366
367 pop_sub(sternheimer_build_fxc)
368 end subroutine sternheimer_build_fxc
369
370
371 !-----------------------------------------------------------
372 subroutine sternheimer_build_kxc(this, namespace, mesh, st, xc)
373 type(sternheimer_t), intent(inout) :: this
374 type(namespace_t), intent(in) :: namespace
375 class(mesh_t), intent(in) :: mesh
376 type(states_elec_t), intent(in) :: st
377 type(xc_t), intent(in) :: xc
378
379 real(real64), allocatable :: rho(:, :)
380
381 push_sub(sternheimer_build_kxc)
382
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))
385 this%kxc = m_zero
386
387 safe_allocate(rho(1:mesh%np, 1:st%d%nspin))
388 call states_elec_total_density(st, mesh, rho)
389 call xc_get_kxc(xc, mesh, namespace, rho, st%d%ispin, this%kxc)
390 safe_deallocate_a(rho)
391 end if
392
393 pop_sub(sternheimer_build_kxc)
394
395 end subroutine sternheimer_build_kxc
396
397 !---------------------------------------------
398 subroutine sternheimer_unset_kxc(this)
399 type(sternheimer_t), intent(inout) :: this
400
401 push_sub(sternheimer_unset_kxc)
402
403 safe_deallocate_a(this%kxc)
404
405 pop_sub(sternheimer_unset_kxc)
406 end subroutine sternheimer_unset_kxc
407
408 !-----------------------------------------------------------
409 pure logical function sternheimer_add_fxc(this) result(rr)
410 class(sternheimer_t), intent(in) :: this
411 rr = this%add_fxc_kernel
412 end function sternheimer_add_fxc
413
414
415 !-----------------------------------------------------------
416 pure logical function sternheimer_add_hartree(this) result(rr)
417 class(sternheimer_t), intent(in) :: this
418 rr = abs(this%coeff_hartree) > m_epsilon
419 end function sternheimer_add_hartree
420
421
422 !-----------------------------------------------------------
423 logical function sternheimer_has_converged(this) result(rr)
424 type(sternheimer_t), intent(in) :: this
425 rr = this%ok
426 end function sternheimer_has_converged
427
428 !-----------------------------------------------------------
429 logical pure function sternheimer_have_rhs(this) result(have)
430 type(sternheimer_t), intent(in) :: this
431 have = associated(this%drhs) .or. associated(this%zrhs)
432 end function sternheimer_have_rhs
433
434 !-----------------------------------------------------------
435 subroutine sternheimer_unset_rhs(this)
436 type(sternheimer_t), intent(inout) :: this
437
438 push_sub(sternheimer_unset_rhs)
439
440 nullify(this%drhs)
441 nullify(this%zrhs)
442
443 pop_sub(sternheimer_unset_rhs)
444 end subroutine sternheimer_unset_rhs
446 !-----------------------------------------------------------
447 logical pure function sternheimer_have_inhomog(this) result(have)
448 type(sternheimer_t), intent(in) :: this
449 have = associated(this%dinhomog) .or. associated(this%zinhomog)
450 end function sternheimer_have_inhomog
451
452 !-----------------------------------------------------------
453 subroutine sternheimer_unset_inhomog(this)
454 type(sternheimer_t), intent(inout) :: this
455
457
458 nullify(this%dinhomog)
459 nullify(this%zinhomog)
460
462 end subroutine sternheimer_unset_inhomog
463
464 !-----------------------------------------------------------
465 integer pure function swap_sigma(sigma)
466 integer, intent(in) :: sigma
468 if (sigma == 1) then
469 swap_sigma = 2
470 else
471 swap_sigma = 1
472 end if
473
474 end function swap_sigma
475
476! ---------------------------------------------------------
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
481
482 character :: sigma_char
483
484 push_sub(wfs_tag_sigma)
485
486 select case (isigma)
487 case (1)
488 sigma_char = '+'
489 case (2)
490 sigma_char = '-'
491 case default
492 write(message(1),'(a,i2)') "Illegal integer isigma passed to wfs_tag_sigma: ", isigma
493 call messages_fatal(1, namespace=namespace)
494 end select
495
496 str = trim(base_name) // sigma_char
497
498 pop_sub(wfs_tag_sigma)
499
500 end function wfs_tag_sigma
501
502 ! --------------------------------------------------------
503
504 subroutine sternheimer_obsolete_variables(namespace, old_prefix, new_prefix)
505 type(namespace_t), intent(in) :: namespace
506 character(len=*), intent(in) :: old_prefix
507 character(len=*), intent(in) :: new_prefix
508
510
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')
513
514 call linear_solver_obsolete_variables(namespace, old_prefix, new_prefix)
515 call scf_tol_obsolete_variables(namespace, old_prefix, new_prefix)
516
519
520 !--------------------------------------------------------------
521 subroutine calc_hvar_photons(this, mesh, nspin, lr_rho, nsigma, hvar, idir)
522 type(sternheimer_t), intent(inout) :: this
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
529
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
534
535 push_sub(calc_hvar_photons)
536 call profiling_in('CALC_HVAR_PHOTONS')
537
538 nm = this%pt_modes%nmodes
539
540 ! photonic terms
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))
545
546 ! spin summed density
547 s_lr_rho = m_zero
548 do is = 1, nspin
549 s_lr_rho = s_lr_rho + lr_rho(:, is)
550 end do
551
552 ! Compute electron-photon potentials
553 vp_dip_self_ener = m_zero
554 vp_bilinear_el_pt = m_zero
555 do ii = 1, nm
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))
558
559 ! Compute photon displacement coordinate q_{\alpha}s
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
564
565 ! Compute potential for bilinear el-pt interaction
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)
568
569 ! Compute potential with dipole-self energy term
570 vp_dip_self_ener = vp_dip_self_ener + first_moments*lambda_dot_r(1:mesh%np)
571 end do
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)
574
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)
579
580 if (nsigma == 2) hvar(1:mesh%np, 1:nspin, 2) = conjg(hvar(1:mesh%np, 1:nspin, 1))
581
582 call profiling_out('CALC_HVAR_PHOTONS')
583 pop_sub(calc_hvar_photons)
584 end subroutine calc_hvar_photons
585
586
587#include "complex.F90"
588#include "sternheimer_inc.F90"
589
590#include "undef.F90"
591
592#include "real.F90"
593#include "sternheimer_inc.F90"
594
595end module sternheimer_oct_m
596
597!! Local Variables:
598!! mode: f90
599!! coding: utf-8
600!! End:
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:181
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:857
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_four
Definition: global.F90:195
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:189
integer, parameter, public independent_particles
Theory level.
Definition: global.F90:241
real(real64), parameter, public m_epsilon
Definition: global.F90:207
real(real64), parameter, public m_one
Definition: global.F90:192
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:361
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:402
A module to handle KS potential, without the external potential.
integer, parameter, public dft_u_none
Definition: lda_u.F90:205
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_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:897
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1067
character(len=512), private msg
Definition: messages.F90:166
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:409
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1039
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:593
subroutine, public mix_get_field(this, mixfield)
Definition: mix.F90:820
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 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
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
This module provices a simple timer class which can be used to trigger the writing of a restart file ...
Definition: walltimer.F90:123
subroutine, public xc_get_fxc(xcs, gr, namespace, rho, ispin, fxc, fxc_grad, fxc_grad_spin)
Returns the exchange-correlation kernel.
Definition: xc_kernel.F90:172
subroutine, public xc_get_kxc(xcs, mesh, namespace, rho, ispin, kxc)
Definition: xc_kernel.F90:514
Definition: xc.F90:116
pure logical function, public family_is_mgga(family, only_collinear)
Is the xc function part of the mGGA family.
Definition: xc.F90:580
logical pure function, public family_is_hybrid(xcs)
Returns true if the functional is an hybrid functional.
Definition: xc.F90:614
subroutine, public xc_write_fxc_info(xcs, iunit, namespace)
Definition: xc.F90:255
integer, parameter, public sic_none
no self-interaction correction
Definition: xc_sic.F90:153
subroutine, public xc_sic_write_info(sic, iunit, namespace)
Definition: xc_sic.F90:259
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:171
Describes mesh distribution to nodes.
Definition: mesh.F90:187
The states_elec_t class contains all electronic wave functions.